gtsam  4.0.0
gtsam
RegularImplicitSchurFactor.h
Go to the documentation of this file.
1 
8 #pragma once
9 
13 #include <iosfwd>
14 
15 namespace gtsam {
16 
20 template<class CAMERA>
22 
23 public:
25  typedef boost::shared_ptr<This> shared_ptr;
26 
27 protected:
28 
29  // This factor is closely related to a CameraSet
30  typedef CameraSet<CAMERA> Set;
31 
32  typedef typename CAMERA::Measurement Z;
33  static const int D = traits<CAMERA>::dimension;
34  static const int ZDim = traits<Z>::dimension;
35 
36  typedef Eigen::Matrix<double, ZDim, D> MatrixZD;
37  typedef Eigen::Matrix<double, D, D> MatrixDD;
38 
39  const std::vector<MatrixZD, Eigen::aligned_allocator<MatrixZD> > FBlocks_;
40  const Matrix PointCovariance_;
41  const Matrix E_;
42  const Vector b_;
43 
44 public:
45 
48  }
49 
52  const std::vector<MatrixZD, Eigen::aligned_allocator<MatrixZD> >& FBlocks, const Matrix& E, const Matrix& P,
53  const Vector& b) :
54  GaussianFactor(keys), FBlocks_(FBlocks), PointCovariance_(P), E_(E), b_(b) {
55  }
56 
59  }
60 
61  std::vector<MatrixZD, Eigen::aligned_allocator<MatrixZD> >& FBlocks() const {
62  return FBlocks_;
63  }
64 
65  const Matrix& E() const {
66  return E_;
67  }
68 
69  const Vector& b() const {
70  return b_;
71  }
72 
73  const Matrix& getPointCovariance() const {
74  return PointCovariance_;
75  }
76 
78  void print(const std::string& s = "", const KeyFormatter& keyFormatter =
79  DefaultKeyFormatter) const {
80  std::cout << " RegularImplicitSchurFactor " << std::endl;
81  Factor::print(s);
82  for (size_t pos = 0; pos < size(); ++pos) {
83  std::cout << "Fblock:\n" << FBlocks_[pos] << std::endl;
84  }
85  std::cout << "PointCovariance:\n" << PointCovariance_ << std::endl;
86  std::cout << "E:\n" << E_ << std::endl;
87  std::cout << "b:\n" << b_.transpose() << std::endl;
88  }
89 
91  bool equals(const GaussianFactor& lf, double tol) const {
92  const This* f = dynamic_cast<const This*>(&lf);
93  if (!f)
94  return false;
95  for (size_t k = 0; k < FBlocks_.size(); ++k) {
96  if (keys_[k] != f->keys_[k])
97  return false;
98  if (!equal_with_abs_tol(FBlocks_[k], f->FBlocks_[k], tol))
99  return false;
100  }
101  return equal_with_abs_tol(PointCovariance_, f->PointCovariance_, tol)
102  && equal_with_abs_tol(E_, f->E_, tol)
103  && equal_with_abs_tol(b_, f->b_, tol);
104  }
105 
107  virtual DenseIndex getDim(const_iterator variable) const {
108  return D;
109  }
110 
111  virtual void updateHessian(const KeyVector& keys,
112  SymmetricBlockMatrix* info) const {
113  throw std::runtime_error(
114  "RegularImplicitSchurFactor::updateHessian non implemented");
115  }
116  virtual Matrix augmentedJacobian() const {
117  throw std::runtime_error(
118  "RegularImplicitSchurFactor::augmentedJacobian non implemented");
119  return Matrix();
120  }
121  virtual std::pair<Matrix, Vector> jacobian() const {
122  throw std::runtime_error(
123  "RegularImplicitSchurFactor::jacobian non implemented");
124  return std::make_pair(Matrix(), Vector());
125  }
126 
128  virtual Matrix augmentedInformation() const {
129 
130  // Do the Schur complement
131  SymmetricBlockMatrix augmentedHessian = //
133  return augmentedHessian.selfadjointView();
134  }
135 
137  virtual Matrix information() const {
138  Matrix augmented = augmentedInformation();
139  int m = this->keys_.size();
140  size_t M = D * m;
141  return augmented.block(0, 0, M, M);
142  }
143 
145  virtual VectorValues hessianDiagonal() const {
146  // diag(Hessian) = diag(F' * (I - E * PointCov * E') * F);
147  VectorValues d;
148 
149  for (size_t k = 0; k < size(); ++k) { // for each camera
150  Key j = keys_[k];
151 
152  // Calculate Fj'*Ej for the current camera (observing a single point)
153  // D x 3 = (D x ZDim) * (ZDim x 3)
154  const MatrixZD& Fj = FBlocks_[k];
155  Eigen::Matrix<double, D, 3> FtE = Fj.transpose()
156  * E_.block<ZDim, 3>(ZDim * k, 0);
157 
158  Eigen::Matrix<double, D, 1> dj;
159  for (int k = 0; k < D; ++k) { // for each diagonal element of the camera hessian
160  // Vector column_k_Fj = Fj.col(k);
161  dj(k) = Fj.col(k).squaredNorm(); // dot(column_k_Fj, column_k_Fj);
162  // Vector column_k_FtE = FtE.row(k);
163  // (1 x 1) = (1 x 3) * (3 * 3) * (3 x 1)
164  dj(k) -= FtE.row(k) * PointCovariance_ * FtE.row(k).transpose();
165  }
166  d.insert(j, dj);
167  }
168  return d;
169  }
170 
175  virtual void hessianDiagonal(double* d) const {
176  // diag(Hessian) = diag(F' * (I - E * PointCov * E') * F);
177  // Use eigen magic to access raw memory
178  typedef Eigen::Matrix<double, D, 1> DVector;
179  typedef Eigen::Map<DVector> DMap;
180 
181  for (size_t pos = 0; pos < size(); ++pos) { // for each camera in the factor
182  Key j = keys_[pos];
183 
184  // Calculate Fj'*Ej for the current camera (observing a single point)
185  // D x 3 = (D x ZDim) * (ZDim x 3)
186  const MatrixZD& Fj = FBlocks_[pos];
187  Eigen::Matrix<double, D, 3> FtE = Fj.transpose()
188  * E_.block<ZDim, 3>(ZDim * pos, 0);
189 
190  DVector dj;
191  for (int k = 0; k < D; ++k) { // for each diagonal element of the camera hessian
192  dj(k) = Fj.col(k).squaredNorm();
193  // (1 x 1) = (1 x 3) * (3 * 3) * (3 x 1)
194  dj(k) -= FtE.row(k) * PointCovariance_ * FtE.row(k).transpose();
195  }
196  DMap(d + D * j) += dj;
197  }
198  }
199 
201  virtual std::map<Key, Matrix> hessianBlockDiagonal() const {
202  std::map<Key, Matrix> blocks;
203  // F'*(I - E*P*E')*F
204  for (size_t pos = 0; pos < size(); ++pos) {
205  Key j = keys_[pos];
206  // F'*F - F'*E*P*E'*F e.g. (9*2)*(2*9) - (9*2)*(2*3)*(3*3)*(3*2)*(2*9)
207  const MatrixZD& Fj = FBlocks_[pos];
208  // Eigen::Matrix<double, D, 3> FtE = Fj.transpose()
209  // * E_.block<ZDim, 3>(ZDim * pos, 0);
210  // blocks[j] = Fj.transpose() * Fj
211  // - FtE * PointCovariance_ * FtE.transpose();
212 
213  const Matrix23& Ej = E_.block<ZDim, 3>(ZDim * pos, 0);
214  blocks[j] = Fj.transpose()
215  * (Fj - Ej * PointCovariance_ * Ej.transpose() * Fj);
216 
217  // F'*(I - E*P*E')*F, TODO: this should work, but it does not :-(
218  // static const Eigen::Matrix<double, ZDim, ZDim> I2 = eye(ZDim);
219  // Matrix2 Q = //
220  // I2 - E_.block<ZDim, 3>(ZDim * pos, 0) * PointCovariance_ * E_.block<ZDim, 3>(ZDim * pos, 0).transpose();
221  // blocks[j] = Fj.transpose() * Q * Fj;
222  }
223  return blocks;
224  }
225 
227  return boost::make_shared<RegularImplicitSchurFactor<CAMERA> >(keys_,
229  throw std::runtime_error(
230  "RegularImplicitSchurFactor::clone non implemented");
231  }
232  virtual bool empty() const {
233  return false;
234  }
235 
237  return boost::make_shared<RegularImplicitSchurFactor<CAMERA> >(keys_,
239  throw std::runtime_error(
240  "RegularImplicitSchurFactor::negate non implemented");
241  }
242 
243  // Raw Vector version of y += F'*alpha*(I - E*P*E')*F*x, for testing
244  static
245  void multiplyHessianAdd(const Matrix& F, const Matrix& E,
246  const Matrix& PointCovariance, double alpha, const Vector& x, Vector& y) {
247  Vector e1 = F * x;
248  Vector d1 = E.transpose() * e1;
249  Vector d2 = PointCovariance * d1;
250  Vector e2 = E * d2;
251  Vector e3 = alpha * (e1 - e2);
252  y += F.transpose() * e3;
253  }
254 
255  typedef std::vector<Vector2, Eigen::aligned_allocator<Vector2>> Error2s;
256 
260  void projectError2(const Error2s& e1, Error2s& e2) const {
261 
262  // d1 = E.transpose() * (e1-ZDim*b) = (3*2m)*2m
263  Vector3 d1;
264  d1.setZero();
265  for (size_t k = 0; k < size(); k++)
266  d1 += E_.block<ZDim, 3>(ZDim * k, 0).transpose()
267  * (e1[k] - ZDim * b_.segment<ZDim>(k * ZDim));
268 
269  // d2 = E.transpose() * e1 = (3*2m)*2m
270  Vector3 d2 = PointCovariance_ * d1;
271 
272  // e3 = alpha*(e1 - E*d2) = 1*[2m-(2m*3)*3]
273  for (size_t k = 0; k < size(); k++)
274  e2[k] = e1[k] - ZDim * b_.segment<ZDim>(k * ZDim)
275  - E_.block<ZDim, 3>(ZDim * k, 0) * d2;
276  }
277 
278  /*
279  * This definition matches the linearized error in the Hessian Factor:
280  * LinError(x) = x'*H*x - 2*x'*eta + f
281  * with:
282  * H = F' * (I-E'*P*E) * F = F' * Q * F
283  * eta = F' * (I-E'*P*E) * b = F' * Q * b
284  * f = nonlinear error
285  * (x'*H*x - 2*x'*eta + f) = x'*F'*Q*F*x - 2*x'*F'*Q *b + f = x'*F'*Q*(F*x - 2*b) + f
286  */
287  virtual double error(const VectorValues& x) const {
288 
289  // resize does not do malloc if correct size
290  e1.resize(size());
291  e2.resize(size());
292 
293  // e1 = F * x - b = (2m*dm)*dm
294  for (size_t k = 0; k < size(); ++k)
295  e1[k] = FBlocks_[k] * x.at(keys_[k]);
296  projectError2(e1, e2);
297 
298  double result = 0;
299  for (size_t k = 0; k < size(); ++k)
300  result += dot(e1[k], e2[k]);
301 
302  double f = b_.squaredNorm();
303  return 0.5 * (result + f);
304  }
305 
306  // needed to be GaussianFactor - (I - E*P*E')*(F*x - b)
307  // This is wrong and does not match the definition in Hessian,
308  // but it matches the definition of the Jacobian factor (JF)
309  double errorJF(const VectorValues& x) const {
310 
311  // resize does not do malloc if correct size
312  e1.resize(size());
313  e2.resize(size());
314 
315  // e1 = F * x - b = (2m*dm)*dm
316  for (size_t k = 0; k < size(); ++k)
317  e1[k] = FBlocks_[k] * x.at(keys_[k]) - b_.segment<ZDim>(k * ZDim);
318  projectError(e1, e2);
319 
320  double result = 0;
321  for (size_t k = 0; k < size(); ++k)
322  result += dot(e2[k], e2[k]);
323 
324  // std::cout << "implicitFactor::error result " << result << std::endl;
325  return 0.5 * result;
326  }
330  void projectError(const Error2s& e1, Error2s& e2) const {
331 
332  // d1 = E.transpose() * e1 = (3*2m)*2m
333  Vector3 d1;
334  d1.setZero();
335  for (size_t k = 0; k < size(); k++)
336  d1 += E_.block<ZDim, 3>(ZDim * k, 0).transpose() * e1[k];
337 
338  // d2 = E.transpose() * e1 = (3*2m)*2m
339  Vector3 d2 = PointCovariance_ * d1;
340 
341  // e3 = alpha*(e1 - E*d2) = 1*[2m-(2m*3)*3]
342  for (size_t k = 0; k < size(); k++)
343  e2[k] = e1[k] - E_.block<ZDim, 3>(ZDim * k, 0) * d2;
344  }
345 
347  mutable Error2s e1, e2;
348 
353  void multiplyHessianAdd(double alpha, const double* x, double* y) const {
354 
355  // Use eigen magic to access raw memory
356  typedef Eigen::Matrix<double, D, 1> DVector;
357  typedef Eigen::Map<DVector> DMap;
358  typedef Eigen::Map<const DVector> ConstDMap;
359 
360  // resize does not do malloc if correct size
361  e1.resize(size());
362  e2.resize(size());
363 
364  // e1 = F * x = (2m*dm)*dm
365  for (size_t k = 0; k < size(); ++k) {
366  Key key = keys_[k];
367  e1[k] = FBlocks_[k] * ConstDMap(x + D * key);
368  }
369 
370  projectError(e1, e2);
371 
372  // y += F.transpose()*e2 = (2d*2m)*2m
373  for (size_t k = 0; k < size(); ++k) {
374  Key key = keys_[k];
375  DMap(y + D * key) += FBlocks_[k].transpose() * alpha * e2[k];
376  }
377  }
378 
379  void multiplyHessianAdd(double alpha, const double* x, double* y,
380  std::vector<size_t> keys) const {
381  }
382  ;
383 
387  void multiplyHessianAdd(double alpha, const VectorValues& x,
388  VectorValues& y) const {
389 
390  // resize does not do malloc if correct size
391  e1.resize(size());
392  e2.resize(size());
393 
394  // e1 = F * x = (2m*dm)*dm
395  for (size_t k = 0; k < size(); ++k)
396  e1[k] = FBlocks_[k] * x.at(keys_[k]);
397 
398  projectError(e1, e2);
399 
400  // y += F.transpose()*e2 = (2d*2m)*2m
401  for (size_t k = 0; k < size(); ++k) {
402  Key key = keys_[k];
403  static const Vector empty;
404  std::pair<VectorValues::iterator, bool> it = y.tryInsert(key, empty);
405  Vector& yi = it.first->second;
406  // Create the value as a zero vector if it does not exist.
407  if (it.second)
408  yi = Vector::Zero(FBlocks_[k].cols());
409  yi += FBlocks_[k].transpose() * alpha * e2[k];
410  }
411  }
412 
416  void multiplyHessianDummy(double alpha, const VectorValues& x,
417  VectorValues& y) const {
418 
419  for (size_t k = 0; k < size(); ++k) {
420  static const Vector empty;
421  Key key = keys_[k];
422  std::pair<VectorValues::iterator, bool> it = y.tryInsert(key, empty);
423  Vector& yi = it.first->second;
424  yi = x.at(key);
425  }
426  }
427 
432  // calculate Q*b
433  e1.resize(size());
434  e2.resize(size());
435  for (size_t k = 0; k < size(); k++)
436  e1[k] = b_.segment<ZDim>(ZDim * k);
437  projectError(e1, e2);
438 
439  // g = F.transpose()*e2
440  VectorValues g;
441  for (size_t k = 0; k < size(); ++k) {
442  Key key = keys_[k];
443  g.insert(key, -FBlocks_[k].transpose() * e2[k]);
444  }
445 
446  // return it
447  return g;
448  }
449 
453  virtual void gradientAtZero(double* d) const {
454 
455  // Use eigen magic to access raw memory
456  typedef Eigen::Matrix<double, D, 1> DVector;
457  typedef Eigen::Map<DVector> DMap;
458 
459  // calculate Q*b
460  e1.resize(size());
461  e2.resize(size());
462  for (size_t k = 0; k < size(); k++)
463  e1[k] = b_.segment<ZDim>(ZDim * k);
464  projectError(e1, e2);
465 
466  for (size_t k = 0; k < size(); ++k) { // for each camera in the factor
467  Key j = keys_[k];
468  DMap(d + D * j) += -FBlocks_[k].transpose() * e2[k];
469  }
470  }
471 
473  Vector gradient(Key key, const VectorValues& x) const {
474  throw std::runtime_error(
475  "gradient for RegularImplicitSchurFactor is not implemented yet");
476  }
477 
478 };
479 // end class RegularImplicitSchurFactor
480 
481 template<class CAMERA>
483 
484 template<class CAMERA>
486 
487 // traits
488 template<class CAMERA> struct traits<RegularImplicitSchurFactor<CAMERA> > : public Testable<
489  RegularImplicitSchurFactor<CAMERA> > {
490 };
491 
492 }
493 
virtual GaussianFactor::shared_ptr negate() const
Construct the corresponding anti-factor to negate information stored stored in this factor.
Definition: RegularImplicitSchurFactor.h:236
RegularImplicitSchurFactor()
Constructor.
Definition: RegularImplicitSchurFactor.h:47
void multiplyHessianAdd(double alpha, const VectorValues &x, VectorValues &y) const
Hessian-vector multiply, i.e.
Definition: RegularImplicitSchurFactor.h:387
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
ptrdiff_t DenseIndex
The index type for Eigen objects.
Definition: types.h:63
virtual std::pair< Matrix, Vector > jacobian() const
Return the dense Jacobian and right-hand-side , with the noise models baked into A and b.
Definition: RegularImplicitSchurFactor.h:121
virtual Matrix augmentedJacobian() const
Return a dense Jacobian matrix, augmented with b with the noise models baked into A and b.
Definition: RegularImplicitSchurFactor.h:116
virtual void hessianDiagonal(double *d) const
add the contribution of this factor to the diagonal of the hessian d(output) = d(input) + deltaHessia...
Definition: RegularImplicitSchurFactor.h:175
static const int ZDim
Measurement dimension.
Definition: RegularImplicitSchurFactor.h:34
void multiplyHessianAdd(double alpha, const double *x, double *y) const
double* Hessian-vector multiply, i.e.
Definition: RegularImplicitSchurFactor.h:353
void projectError2(const Error2s &e1, Error2s &e2) const
Calculate corrected error Q*(e-ZDim*b) = (I - E*P*E')*(e-ZDim*b)
Definition: RegularImplicitSchurFactor.h:260
A set of cameras, all with their own calibration.
Definition: CameraSet.h:35
Error2s e1
Scratch space for multiplyHessianAdd.
Definition: RegularImplicitSchurFactor.h:347
bool equal_with_abs_tol(const Eigen::DenseBase< MATRIX > &A, const Eigen::DenseBase< MATRIX > &B, double tol=1e-9)
equals with a tolerance
Definition: Matrix.h:82
const std::vector< MatrixZD, Eigen::aligned_allocator< MatrixZD > > FBlocks_
All ZDim*D F blocks (one for each camera)
Definition: RegularImplicitSchurFactor.h:39
std::uint64_t Key
Integer nonlinear key type.
Definition: types.h:57
Eigen::Matrix< double, ZDim, D > MatrixZD
type of an F block
Definition: RegularImplicitSchurFactor.h:36
A helper that implements the traits interface for GTSAM types.
Definition: Testable.h:150
An abstract virtual base class for JacobianFactor and HessianFactor.
Definition: GaussianFactor.h:38
virtual bool empty() const
Test whether the factor is empty.
Definition: RegularImplicitSchurFactor.h:232
boost::shared_ptr< This > shared_ptr
shared_ptr to this class
Definition: RegularImplicitSchurFactor.h:25
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
virtual DenseIndex getDim(const_iterator variable) const
Degrees of freedom of camera.
Definition: RegularImplicitSchurFactor.h:107
Eigen::SelfAdjointView< constBlock, Eigen::Upper > selfadjointView(DenseIndex I, DenseIndex J) const
Return the square sub-matrix that contains blocks(i:j, i:j).
Definition: SymmetricBlockMatrix.h:161
RegularImplicitSchurFactor.
Definition: RegularImplicitSchurFactor.h:21
Eigen::Matrix< double, D, D > MatrixDD
camera hessian
Definition: RegularImplicitSchurFactor.h:37
virtual void gradientAtZero(double *d) const
Calculate gradient, which is -F'Q*b, see paper - RAW MEMORY ACCESS.
Definition: RegularImplicitSchurFactor.h:453
size_t size() const
Definition: Factor.h:129
virtual Matrix information() const
Compute full information matrix
Definition: RegularImplicitSchurFactor.h:137
This class represents a collection of vector-valued variables associated each with a unique integer i...
Definition: VectorValues.h:73
const Matrix PointCovariance_
the 3*3 matrix P = inv(E'E) (2*2 if degenerate)
Definition: RegularImplicitSchurFactor.h:40
void print(const std::string &s="Factor", const KeyFormatter &formatter=DefaultKeyFormatter) const
print
Definition: Factor.cpp:29
virtual Matrix augmentedInformation() const
Compute full augmented information matrix
Definition: RegularImplicitSchurFactor.h:128
FastVector< Key > KeyVector
Define collection type once and for all - also used in wrappers.
Definition: Key.h:56
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 GaussianFactor::shared_ptr clone() const
Clone a factor (make a deep copy)
Definition: RegularImplicitSchurFactor.h:226
RegularImplicitSchurFactor(const KeyVector &keys, const std::vector< MatrixZD, Eigen::aligned_allocator< MatrixZD > > &FBlocks, const Matrix &E, const Matrix &P, const Vector &b)
Construct from blocks of F, E, inv(E'*E), and RHS vector b.
Definition: RegularImplicitSchurFactor.h:51
void print(const std::string &s="", const KeyFormatter &keyFormatter=DefaultKeyFormatter) const
print
Definition: RegularImplicitSchurFactor.h:78
virtual VectorValues hessianDiagonal() const
Return the diagonal of the Hessian for this factor.
Definition: RegularImplicitSchurFactor.h:145
virtual ~RegularImplicitSchurFactor()
Destructor.
Definition: RegularImplicitSchurFactor.h:58
static SymmetricBlockMatrix SchurComplement(const FBlocks &Fs, const Matrix &E, const Eigen::Matrix< double, N, N > &P, const Vector &b)
Do Schur complement, given Jacobian as Fs,E,P, return SymmetricBlockMatrix G = F' * F - F' * E * P * ...
Definition: CameraSet.h:149
Base class to create smart factors on poses or cameras.
Definition: SymmetricBlockMatrix.h:51
const KeyVector & keys() const
Access the factor's involved variable keys.
Definition: Factor.h:118
Global functions in a separate testing namespace.
Definition: chartTesting.h:28
KeyVector keys_
The keys involved in this factor.
Definition: Factor.h:72
virtual void updateHessian(const KeyVector &keys, SymmetricBlockMatrix *info) const
Update an information matrix by adding the information corresponding to this factor (used internally ...
Definition: RegularImplicitSchurFactor.h:111
void projectError(const Error2s &e1, Error2s &e2) const
Calculate corrected error Q*e = (I - E*P*E')*e.
Definition: RegularImplicitSchurFactor.h:330
iterator insert(const std::pair< Key, Vector > &key_value)
Insert a vector value with key j.
Definition: VectorValues.cpp:78
static const int D
Camera dimension.
Definition: RegularImplicitSchurFactor.h:33
virtual std::map< Key, Matrix > hessianBlockDiagonal() const
Return the block diagonal of the Hessian for this factor.
Definition: RegularImplicitSchurFactor.h:201
void multiplyHessianDummy(double alpha, const VectorValues &x, VectorValues &y) const
Dummy version to measure overhead of key access.
Definition: RegularImplicitSchurFactor.h:416
KeyVector::const_iterator const_iterator
Const iterator over keys.
Definition: Factor.h:67
const Vector b_
2m-dimensional RHS vector
Definition: RegularImplicitSchurFactor.h:42
bool equals(const GaussianFactor &lf, double tol) const
equals
Definition: RegularImplicitSchurFactor.h:91
std::pair< iterator, bool > tryInsert(Key j, const Vector &value)
insert that mimics the STL map insert - if the value already exists, the map is not modified and an i...
Definition: VectorValues.h:199
VectorValues gradientAtZero() const
Calculate gradient, which is -F'Q*b, see paper.
Definition: RegularImplicitSchurFactor.h:431
Factor Graph Values.
Vector gradient(Key key, const VectorValues &x) const
Gradient wrt a key at any values.
Definition: RegularImplicitSchurFactor.h:473
RegularImplicitSchurFactor This
Typedef to this class.
Definition: RegularImplicitSchurFactor.h:24
virtual double error(const VectorValues &x) const
Print for testable.
Definition: RegularImplicitSchurFactor.h:287
const Matrix E_
The 2m*3 E Jacobian with respect to the point.
Definition: RegularImplicitSchurFactor.h:41
Vector & at(Key j)
Read/write access to the vector value with key j, throws std::out_of_range if j does not exist,...
Definition: VectorValues.h:136
boost::shared_ptr< This > shared_ptr
shared_ptr to this class
Definition: GaussianFactor.h:42