gtsam  4.0.0
gtsam
NonlinearEquality.h
1 /* ----------------------------------------------------------------------------
2 
3  * GTSAM Copyright 2010, Georgia Tech Research Corporation,
4  * Atlanta, Georgia 30332-0415
5  * All Rights Reserved
6  * Authors: Frank Dellaert, et al. (see THANKS for the full author list)
7 
8  * See LICENSE for the license information
9 
10  * -------------------------------------------------------------------------- */
11 
12 /*
13  * @file NonlinearEquality.h
14  * @brief Factor to handle enforced equality between factors
15  * @author Alex Cunningham
16  */
17 
18 #pragma once
19 
21 #include <gtsam/base/Testable.h>
22 #include <gtsam/base/Manifold.h>
23 
24 #include <boost/bind.hpp>
25 
26 #include <limits>
27 #include <iostream>
28 #include <cmath>
29 
30 namespace gtsam {
31 
44 template<class VALUE>
45 class NonlinearEquality: public NoiseModelFactor1<VALUE> {
46 
47 public:
48  typedef VALUE T;
49 
50 private:
51 
52  // feasible value
53  T feasible_;
54 
55  // error handling flag
56  bool allow_error_;
57 
58  // error gain in allow error case
59  double error_gain_;
60 
61  // typedef to this class
63 
64  // typedef to base class
66 
67 public:
68 
72  typedef boost::function<bool(const T&, const T&)> CompareFunction;
73  CompareFunction compare_;
74 // bool (*compare_)(const T& a, const T& b);
75 
78  }
79 
80  virtual ~NonlinearEquality() {
81  }
82 
85 
89  NonlinearEquality(Key j, const T& feasible,
90  const CompareFunction &_compare = boost::bind(traits<T>::Equals,_1,_2,1e-9)) :
91  Base(noiseModel::Constrained::All(traits<T>::GetDimension(feasible)),
92  j), feasible_(feasible), allow_error_(false), error_gain_(0.0), //
93  compare_(_compare) {
94  }
95 
99  NonlinearEquality(Key j, const T& feasible, double error_gain,
100  const CompareFunction &_compare = boost::bind(traits<T>::Equals,_1,_2,1e-9)) :
101  Base(noiseModel::Constrained::All(traits<T>::GetDimension(feasible)),
102  j), feasible_(feasible), allow_error_(true), error_gain_(error_gain), //
103  compare_(_compare) {
104  }
105 
109 
110  virtual void print(const std::string& s = "",
111  const KeyFormatter& keyFormatter = DefaultKeyFormatter) const {
112  std::cout << s << "Constraint: on [" << keyFormatter(this->key()) << "]\n";
113  traits<VALUE>::Print(feasible_, "Feasible Point:\n");
114  std::cout << "Variable Dimension: " << traits<T>::GetDimension(feasible_) << std::endl;
115  }
116 
118  virtual bool equals(const NonlinearFactor& f, double tol = 1e-9) const {
119  const This* e = dynamic_cast<const This*>(&f);
120  return e && Base::equals(f) && traits<T>::Equals(feasible_,e->feasible_, tol)
121  && std::abs(error_gain_ - e->error_gain_) < tol;
122  }
123 
127 
129  virtual double error(const Values& c) const {
130  const T& xj = c.at<T>(this->key());
131  Vector e = this->unwhitenedError(c);
132  if (allow_error_ || !compare_(xj, feasible_)) {
133  return error_gain_ * dot(e, e);
134  } else {
135  return 0.0;
136  }
137  }
138 
140  Vector evaluateError(const T& xj,
141  boost::optional<Matrix&> H = boost::none) const {
142  const size_t nj = traits<T>::GetDimension(feasible_);
143  if (allow_error_) {
144  if (H)
145  *H = Matrix::Identity(nj,nj); // FIXME: this is not the right linearization for nonlinear compare
146  return traits<T>::Local(xj,feasible_);
147  } else if (compare_(feasible_, xj)) {
148  if (H)
149  *H = Matrix::Identity(nj,nj);
150  return Vector::Zero(nj); // set error to zero if equal
151  } else {
152  if (H)
153  throw std::invalid_argument(
154  "Linearization point not feasible for "
155  + DefaultKeyFormatter(this->key()) + "!");
156  return Vector::Constant(nj, std::numeric_limits<double>::infinity()); // set error to infinity if not equal
157  }
158  }
159 
160  // Linearize is over-written, because base linearization tries to whiten
161  virtual GaussianFactor::shared_ptr linearize(const Values& x) const {
162  const T& xj = x.at<T>(this->key());
163  Matrix A;
164  Vector b = evaluateError(xj, A);
165  SharedDiagonal model = noiseModel::Constrained::All(b.size());
167  new JacobianFactor(this->key(), A, b, model));
168  }
169 
171  virtual gtsam::NonlinearFactor::shared_ptr clone() const {
172  return boost::static_pointer_cast<gtsam::NonlinearFactor>(
173  gtsam::NonlinearFactor::shared_ptr(new This(*this)));
174  }
175 
177 
178 private:
179 
182  template<class ARCHIVE>
183  void serialize(ARCHIVE & ar, const unsigned int /*version*/) {
184  ar
185  & boost::serialization::make_nvp("NoiseModelFactor1",
186  boost::serialization::base_object<Base>(*this));
187  ar & BOOST_SERIALIZATION_NVP(feasible_);
188  ar & BOOST_SERIALIZATION_NVP(allow_error_);
189  ar & BOOST_SERIALIZATION_NVP(error_gain_);
190  }
191 
192 };
193 // \class NonlinearEquality
194 
195 template<typename VALUE>
196 struct traits<NonlinearEquality<VALUE> > : Testable<NonlinearEquality<VALUE> > {
197 };
198 
199 /* ************************************************************************* */
203 template<class VALUE>
204 class NonlinearEquality1: public NoiseModelFactor1<VALUE> {
205 
206 public:
207  typedef VALUE X;
208 
209 protected:
212 
215  }
216 
217  X value_;
218 
219  GTSAM_CONCEPT_MANIFOLD_TYPE(X)
220 
221  GTSAM_CONCEPT_TESTABLE_TYPE(X)
222 
223 public:
224 
225  typedef boost::shared_ptr<NonlinearEquality1<VALUE> > shared_ptr;
226 
234  NonlinearEquality1(const X& value, Key key, double mu = 1000.0) :
235  Base( noiseModel::Constrained::All(traits<X>::GetDimension(value),
236  std::abs(mu)), key), value_(value) {
237  }
238 
239  virtual ~NonlinearEquality1() {
240  }
241 
243  virtual gtsam::NonlinearFactor::shared_ptr clone() const {
244  return boost::static_pointer_cast<gtsam::NonlinearFactor>(
245  gtsam::NonlinearFactor::shared_ptr(new This(*this)));
246  }
247 
249  Vector evaluateError(const X& x1,
250  boost::optional<Matrix&> H = boost::none) const {
251  if (H)
252  (*H) = Matrix::Identity(traits<X>::GetDimension(x1),traits<X>::GetDimension(x1));
253  // manifold equivalent of h(x)-z -> log(z,h(x))
254  return traits<X>::Local(value_,x1);
255  }
256 
258  virtual void print(const std::string& s = "",
259  const KeyFormatter& keyFormatter = DefaultKeyFormatter) const {
260  std::cout << s << ": NonlinearEquality1(" << keyFormatter(this->key())
261  << ")," << "\n";
262  this->noiseModel_->print();
263  traits<X>::Print(value_, "Value");
264  }
265 
266 private:
267 
270  template<class ARCHIVE>
271  void serialize(ARCHIVE & ar, const unsigned int /*version*/) {
272  ar
273  & boost::serialization::make_nvp("NoiseModelFactor1",
274  boost::serialization::base_object<Base>(*this));
275  ar & BOOST_SERIALIZATION_NVP(value_);
276  }
277 };
278 // \NonlinearEquality1
279 
280 template<typename VALUE>
281 struct traits<NonlinearEquality1<VALUE> > : Testable<NonlinearEquality1<VALUE> > {
282 };
283 
284 /* ************************************************************************* */
289 template<class VALUE>
290 class NonlinearEquality2: public NoiseModelFactor2<VALUE, VALUE> {
291 public:
292  typedef VALUE X;
293 
294 protected:
297 
298  GTSAM_CONCEPT_MANIFOLD_TYPE(X)
299 
300 
302  }
303 
304 public:
305 
306  typedef boost::shared_ptr<NonlinearEquality2<VALUE> > shared_ptr;
307 
309  NonlinearEquality2(Key key1, Key key2, double mu = 1000.0) :
310  Base(noiseModel::Constrained::All(traits<X>::dimension, std::abs(mu)), key1, key2) {
311  }
312  virtual ~NonlinearEquality2() {
313  }
314 
316  virtual gtsam::NonlinearFactor::shared_ptr clone() const {
317  return boost::static_pointer_cast<gtsam::NonlinearFactor>(
318  gtsam::NonlinearFactor::shared_ptr(new This(*this)));
319  }
320 
322  Vector evaluateError(const X& x1, const X& x2, boost::optional<Matrix&> H1 =
323  boost::none, boost::optional<Matrix&> H2 = boost::none) const {
324  static const size_t p = traits<X>::dimension;
325  if (H1) *H1 = -Matrix::Identity(p,p);
326  if (H2) *H2 = Matrix::Identity(p,p);
327  return traits<X>::Local(x1,x2);
328  }
329 
330 private:
331 
334  template<class ARCHIVE>
335  void serialize(ARCHIVE & ar, const unsigned int /*version*/) {
336  ar
337  & boost::serialization::make_nvp("NoiseModelFactor2",
338  boost::serialization::base_object<Base>(*this));
339  }
340 };
341 // \NonlinearEquality2
342 
343 template<typename VALUE>
344 struct traits<NonlinearEquality2<VALUE> > : Testable<NonlinearEquality2<VALUE> > {
345 };
346 
347 
348 }// namespace gtsam
349 
Vector evaluateError(const X &x1, const X &x2, boost::optional< Matrix & > H1=boost::none, boost::optional< Matrix & > H2=boost::none) const
g(x) with optional derivative2
Definition: NonlinearEquality.h:322
Base class and basic functions for Manifold types.
double dot(const V1 &a, const V2 &b)
Dot product.
Definition: Vector.h:162
This is the base class for all factor types.
Definition: Factor.h:54
Vector evaluateError(const T &xj, boost::optional< Matrix & > H=boost::none) const
error function
Definition: NonlinearEquality.h:140
A non-templated config holding any types of Manifold-group elements.
Definition: Values.h:70
A convenient base class for creating your own NoiseModelFactor with 2 variables.
Definition: NonlinearFactor.h:345
virtual Vector unwhitenedError(const Values &x, boost::optional< std::vector< Matrix > & > H=boost::none) const =0
Error function without the NoiseModel, .
ValueType at(Key j) const
Retrieve a variable by key j.
Definition: Values-inl.h:342
virtual gtsam::NonlinearFactor::shared_ptr clone() const
Definition: NonlinearEquality.h:316
friend class boost::serialization::access
Serialization function.
Definition: NonlinearEquality.h:269
static shared_ptr All(size_t dim)
Fully constrained variations.
Definition: NoiseModel.h:460
An equality factor that forces either one variable to a constant, or a set of variables to be equal t...
Definition: NonlinearEquality.h:45
virtual GaussianFactor::shared_ptr linearize(const Values &x) const
linearize to a GaussianFactor
Definition: NonlinearEquality.h:161
std::uint64_t Key
Integer nonlinear key type.
Definition: types.h:57
NonlinearEquality1()
default constructor to allow for serialization
Definition: NonlinearEquality.h:214
friend class boost::serialization::access
Serialization function.
Definition: NonlinearEquality.h:181
A helper that implements the traits interface for GTSAM types.
Definition: Testable.h:150
Simple unary equality constraint - fixes a value for a variable.
Definition: NonlinearEquality.h:204
NonlinearEquality()
default constructor - only for serialization
Definition: NonlinearEquality.h:77
NonlinearEquality(Key j, const T &feasible, const CompareFunction &_compare=boost::bind(traits< T >::Equals, _1, _2, 1e-9))
Constructor - forces exact evaluation.
Definition: NonlinearEquality.h:89
boost::function< bool(const T &, const T &)> CompareFunction
Function that compares two values.
Definition: NonlinearEquality.h:72
boost::function< std::string(Key)> KeyFormatter
Typedef for a function to format a key, i.e. to convert it to a string.
Definition: Key.h:33
A convenient base class for creating your own NoiseModelFactor with 1 variable.
Definition: NonlinearFactor.h:276
NonlinearEquality2(Key key1, Key key2, double mu=1000.0)
TODO: comment.
Definition: NonlinearEquality.h:309
Nonlinear factor base class.
Definition: NonlinearFactor.h:50
virtual void print(const std::string &s="", const KeyFormatter &keyFormatter=DefaultKeyFormatter) const
Print.
Definition: NonlinearEquality.h:110
Simple binary equality constraint - this constraint forces two factors to be the same.
Definition: NonlinearEquality.h:290
const SharedNoiseModel & noiseModel() const
access to the noise model
Definition: NonlinearFactor.h:210
Key key1() const
methods to retrieve both keys
Definition: NonlinearFactor.h:377
Vector evaluateError(const X &x1, boost::optional< Matrix & > H=boost::none) const
g(x) with optional derivative
Definition: NonlinearEquality.h:249
A manifold defines a space in which there is a notion of a linear tangent space that can be centered ...
Definition: concepts.h:30
virtual gtsam::NonlinearFactor::shared_ptr clone() const
Definition: NonlinearEquality.h:243
Non-linear factor base classes.
A Gaussian factor in the squared-error form.
Definition: JacobianFactor.h:87
virtual bool equals(const NonlinearFactor &f, double tol=1e-9) const
Check if two factors are equal.
Definition: NonlinearEquality.h:118
virtual double error(const Values &c) const
actual error function calculation
Definition: NonlinearEquality.h:129
NonlinearEquality1(const X &value, Key key, double mu=1000.0)
Constructor.
Definition: NonlinearEquality.h:234
Global functions in a separate testing namespace.
Definition: chartTesting.h:28
Concept check for values that can be used in unit tests.
boost::shared_ptr< NonlinearEquality1< VALUE > > shared_ptr
fixed value for variable
Definition: NonlinearEquality.h:225
virtual void print(const std::string &s="", const KeyFormatter &keyFormatter=DefaultKeyFormatter) const
Print.
Definition: NonlinearEquality.h:258
virtual gtsam::NonlinearFactor::shared_ptr clone() const
Definition: NonlinearEquality.h:171
friend class boost::serialization::access
Serialization function.
Definition: NonlinearEquality.h:333
virtual bool equals(const NonlinearFactor &f, double tol=1e-9) const
Check if two factors are equal.
Definition: NonlinearFactor.cpp:71
NonlinearEquality(Key j, const T &feasible, double error_gain, const CompareFunction &_compare=boost::bind(traits< T >::Equals, _1, _2, 1e-9))
Constructor - allows inexact evaluation.
Definition: NonlinearEquality.h:99
boost::shared_ptr< This > shared_ptr
shared_ptr to this class
Definition: GaussianFactor.h:42