GaussianFactorGraph.cpp
Go to the documentation of this file.
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 
29 #include <gtsam/base/debug.h>
30 #include <gtsam/base/timing.h>
31 #include <gtsam/base/cholesky.h>
32 
33 using namespace std;
34 using namespace gtsam;
35 
36 namespace gtsam {
37 
38  // Instantiate base classes
39  template class FactorGraph<GaussianFactor>;
41 
42  /* ************************************************************************* */
43  bool GaussianFactorGraph::equals(const This& fg, double tol) const
44  {
45  return Base::equals(fg, tol);
46  }
47 
48  /* ************************************************************************* */
50  KeySet keys;
51  for (const sharedFactor& factor: *this)
52  if (factor)
53  keys.insert(factor->begin(), factor->end());
54  return keys;
55  }
56 
57  /* ************************************************************************* */
58  std::map<Key, size_t> GaussianFactorGraph::getKeyDimMap() const {
59  map<Key, size_t> spec;
60  for (const GaussianFactor::shared_ptr& gf : *this) {
61  for (GaussianFactor::const_iterator it = gf->begin(); it != gf->end(); it++) {
62  map<Key,size_t>::iterator it2 = spec.find(*it);
63  if ( it2 == spec.end() ) {
64  spec.emplace(*it, gf->getDim(it));
65  }
66  }
67  }
68  return spec;
69  }
70 
71  /* ************************************************************************* */
72  GaussianFactorGraph::shared_ptr GaussianFactorGraph::cloneToPtr() const {
74  *result = *this;
75  return result;
76  }
77 
78  /* ************************************************************************* */
79  GaussianFactorGraph GaussianFactorGraph::clone() const {
81  for (const sharedFactor& f : *this) {
82  if (f)
83  result.push_back(f->clone());
84  else
85  result.push_back(sharedFactor()); // Passes on null factors so indices remain valid
86  }
87  return result;
88  }
89 
90  /* ************************************************************************* */
93  for (const sharedFactor& factor: *this) {
94  if (factor)
95  result.push_back(factor->negate());
96  else
97  result.push_back(sharedFactor()); // Passes on null factors so indices remain valid
98  }
99  return result;
100  }
101 
102  /* ************************************************************************* */
103  using SparseTriplets = std::vector<std::tuple<int, int, double> >;
104  SparseTriplets GaussianFactorGraph::sparseJacobian(const Ordering& ordering,
105  size_t& nrows,
106  size_t& ncols) const {
107  gttic_(GaussianFactorGraph_sparseJacobian);
108  // First find dimensions of each variable
109  typedef std::map<Key, size_t> KeySizeMap;
110  KeySizeMap dims;
111  for (const auto& factor : *this) {
112  if (!static_cast<bool>(factor)) continue;
113 
114  for (auto it = factor->begin(); it != factor->end(); ++it) {
115  dims[*it] = factor->getDim(it);
116  }
117  }
118 
119  // Compute first scalar column of each variable
120  ncols = 0;
121  KeySizeMap columnIndices = dims;
122  for (const auto key : ordering) {
123  columnIndices[key] = ncols;
124  ncols += dims[key];
125  }
126 
127  // Iterate over all factors, adding sparse scalar entries
128  SparseTriplets entries;
129  entries.reserve(30 * size());
130  nrows = 0;
131  for (const auto& factor : *this) {
132  if (!static_cast<bool>(factor)) continue;
133 
134  // Convert to JacobianFactor if necessary
135  JacobianFactor::shared_ptr jacobianFactor(
136  boost::dynamic_pointer_cast<JacobianFactor>(factor));
137  if (!jacobianFactor) {
139  boost::dynamic_pointer_cast<HessianFactor>(factor));
140  if (hessian)
141  jacobianFactor.reset(new JacobianFactor(*hessian));
142  else
143  throw std::invalid_argument(
144  "GaussianFactorGraph contains a factor that is neither a "
145  "JacobianFactor nor a HessianFactor.");
146  }
147 
148  // Whiten the factor and add entries for it
149  // iterate over all variables in the factor
150  const JacobianFactor whitened(jacobianFactor->whiten());
151  for (auto key = whitened.begin(); key < whitened.end(); ++key) {
152  JacobianFactor::constABlock whitenedA = whitened.getA(key);
153  // find first column index for this key
154  size_t column_start = columnIndices[*key];
155  for (size_t i = 0; i < (size_t)whitenedA.rows(); i++)
156  for (size_t j = 0; j < (size_t)whitenedA.cols(); j++) {
157  double s = whitenedA(i, j);
158  if (std::abs(s) > 1e-12)
159  entries.emplace_back(nrows + i, column_start + j, s);
160  }
161  }
162 
163  JacobianFactor::constBVector whitenedb(whitened.getb());
164  for (size_t i = 0; i < (size_t)whitenedb.size(); i++) {
165  double s = whitenedb(i);
166  if (std::abs(s) > 1e-12) entries.emplace_back(nrows + i, ncols, s);
167  }
168 
169  // Increment row index
170  nrows += jacobianFactor->rows();
171  }
172 
173  ncols++; // +1 for b-column
174  return entries;
175  }
176 
177  /* ************************************************************************* */
178  SparseTriplets GaussianFactorGraph::sparseJacobian() const {
179  size_t nrows, ncols;
180  return sparseJacobian(Ordering(this->keys()), nrows, ncols);
181  }
182 
183  /* ************************************************************************* */
184  Matrix GaussianFactorGraph::sparseJacobian_() const {
185  gttic_(GaussianFactorGraph_sparseJacobian_matrix);
186  // call sparseJacobian
187  auto result = sparseJacobian();
188 
189  // translate to base 1 matrix
190  size_t nzmax = result.size();
191  Matrix IJS(3, nzmax);
192  for (size_t k = 0; k < result.size(); k++) {
193  const auto& entry = result[k];
194  IJS(0, k) = double(std::get<0>(entry) + 1);
195  IJS(1, k) = double(std::get<1>(entry) + 1);
196  IJS(2, k) = std::get<2>(entry);
197  }
198  return IJS;
199  }
200 
201  /* ************************************************************************* */
202  Matrix GaussianFactorGraph::augmentedJacobian(
203  const Ordering& ordering) const {
204  // combine all factors
205  JacobianFactor combined(*this, ordering);
206  return combined.augmentedJacobian();
207  }
208 
209  /* ************************************************************************* */
210  Matrix GaussianFactorGraph::augmentedJacobian() const {
211  // combine all factors
212  JacobianFactor combined(*this);
213  return combined.augmentedJacobian();
214  }
215 
216  /* ************************************************************************* */
217  pair<Matrix, Vector> GaussianFactorGraph::jacobian(
218  const Ordering& ordering) const {
219  Matrix augmented = augmentedJacobian(ordering);
220  return make_pair(augmented.leftCols(augmented.cols() - 1),
221  augmented.col(augmented.cols() - 1));
222  }
223 
224  /* ************************************************************************* */
225  pair<Matrix, Vector> GaussianFactorGraph::jacobian() const {
226  Matrix augmented = augmentedJacobian();
227  return make_pair(augmented.leftCols(augmented.cols() - 1),
228  augmented.col(augmented.cols() - 1));
229  }
230 
231  /* ************************************************************************* */
232  Matrix GaussianFactorGraph::augmentedHessian(
233  const Ordering& ordering) const {
234  // combine all factors and get upper-triangular part of Hessian
235  Scatter scatter(*this, ordering);
236  HessianFactor combined(*this, scatter);
237  return combined.info().selfadjointView();;
238  }
239 
240  /* ************************************************************************* */
241  Matrix GaussianFactorGraph::augmentedHessian() const {
242  // combine all factors and get upper-triangular part of Hessian
243  Scatter scatter(*this);
244  HessianFactor combined(*this, scatter);
245  return combined.info().selfadjointView();;
246  }
247 
248  /* ************************************************************************* */
249  pair<Matrix, Vector> GaussianFactorGraph::hessian(
250  const Ordering& ordering) const {
251  Matrix augmented = augmentedHessian(ordering);
252  size_t n = augmented.rows() - 1;
253  return make_pair(augmented.topLeftCorner(n, n), augmented.topRightCorner(n, 1));
254  }
255 
256  /* ************************************************************************* */
257  pair<Matrix, Vector> GaussianFactorGraph::hessian() const {
258  Matrix augmented = augmentedHessian();
259  size_t n = augmented.rows() - 1;
260  return make_pair(augmented.topLeftCorner(n, n), augmented.topRightCorner(n, 1));
261  }
262 
263  /* ************************************************************************* */
264  VectorValues GaussianFactorGraph::hessianDiagonal() const {
265  VectorValues d;
266  for (const sharedFactor& factor : *this) {
267  if(factor){
268  factor->hessianDiagonalAdd(d);
269  }
270  }
271  return d;
272  }
273 
274  /* ************************************************************************* */
275  map<Key,Matrix> GaussianFactorGraph::hessianBlockDiagonal() const {
276  map<Key,Matrix> blocks;
277  for (const sharedFactor& factor : *this) {
278  if (!factor) continue;
279  map<Key,Matrix> BD = factor->hessianBlockDiagonal();
280  map<Key,Matrix>::const_iterator it = BD.begin();
281  for (;it!=BD.end();++it) {
282  Key j = it->first; // variable key for this block
283  const Matrix& Bj = it->second;
284  if (blocks.count(j))
285  blocks[j] += Bj;
286  else
287  blocks.emplace(j,Bj);
288  }
289  }
290  return blocks;
291  }
292 
293  /* ************************************************************************* */
295  gttic(GaussianFactorGraph_optimize);
296  return BaseEliminateable::eliminateMultifrontal(function)->optimize();
297  }
298 
299  /* ************************************************************************* */
301  gttic(GaussianFactorGraph_optimize);
302  return BaseEliminateable::eliminateMultifrontal(ordering, function)->optimize();
303  }
304 
305  /* ************************************************************************* */
306  // TODO(frank): can we cache memory across invocations
307  VectorValues GaussianFactorGraph::optimizeDensely() const {
308  gttic(GaussianFactorGraph_optimizeDensely);
309 
310  // Combine all factors in a single HessianFactor (as done in augmentedHessian)
311  Scatter scatter(*this);
312  HessianFactor combined(*this, scatter);
313 
314  // TODO(frank): cast to large dynamic matrix :-(
315  // NOTE(frank): info only valid (I think) in upper triangle. No problem for LLT...
316  Matrix augmented = combined.info().selfadjointView();
317 
318  // Do Cholesky Factorization
319  size_t n = augmented.rows() - 1;
320  auto AtA = augmented.topLeftCorner(n, n);
321  auto eta = augmented.topRightCorner(n, 1);
323 
324  // Solve and convert, re-using scatter data structure
325  Vector solution = llt.solve(eta);
326  return VectorValues(solution, scatter);
327  }
328 
329  /* ************************************************************************* */
330  namespace {
331  JacobianFactor::shared_ptr convertToJacobianFactorPtr(const GaussianFactor::shared_ptr &gf) {
332  JacobianFactor::shared_ptr result = boost::dynamic_pointer_cast<JacobianFactor>(gf);
333  if( !result )
334  // Convert any non-Jacobian factors to Jacobians (e.g. Hessian -> Jacobian with Cholesky)
335  result = boost::make_shared<JacobianFactor>(*gf);
336  return result;
337  }
338  }
339 
340  /* ************************************************************************* */
341  VectorValues GaussianFactorGraph::gradient(const VectorValues& x0) const
342  {
343  VectorValues g = VectorValues::Zero(x0);
344  for (const sharedFactor& factor: *this) {
345  JacobianFactor::shared_ptr Ai = convertToJacobianFactorPtr(factor);
346  Vector e = Ai->error_vector(x0);
347  Ai->transposeMultiplyAdd(1.0, e, g);
348  }
349  return g;
350  }
351 
352  /* ************************************************************************* */
353  VectorValues GaussianFactorGraph::gradientAtZero() const {
354  // Zero-out the gradient
355  VectorValues g;
356  for (const sharedFactor& factor: *this) {
357  if (!factor) continue;
358  VectorValues gi = factor->gradientAtZero();
359  g.addInPlace_(gi);
360  }
361  return g;
362  }
363 
364  /* ************************************************************************* */
365  VectorValues GaussianFactorGraph::optimizeGradientSearch() const
366  {
367  gttic(GaussianFactorGraph_optimizeGradientSearch);
368 
369  gttic(Compute_Gradient);
370  // Compute gradient (call gradientAtZero function, which is defined for various linear systems)
371  VectorValues grad = gradientAtZero();
372  double gradientSqNorm = grad.dot(grad);
373  gttoc(Compute_Gradient);
374 
375  gttic(Compute_Rg);
376  // Compute R * g
377  Errors Rg = *this * grad;
378  gttoc(Compute_Rg);
379 
380  gttic(Compute_minimizing_step_size);
381  // Compute minimizing step size
382  double step = -gradientSqNorm / dot(Rg, Rg);
383  gttoc(Compute_minimizing_step_size);
384 
385  gttic(Compute_point);
386  // Compute steepest descent point
387  grad *= step;
388  gttoc(Compute_point);
389 
390  return grad;
391  }
392 
393  /* ************************************************************************* */
395  Errors e;
396  for (const GaussianFactor::shared_ptr& factor: *this) {
397  JacobianFactor::shared_ptr Ai = convertToJacobianFactorPtr(factor);
398  e.push_back((*Ai) * x);
399  }
400  return e;
401  }
402 
403  /* ************************************************************************* */
404  void GaussianFactorGraph::multiplyHessianAdd(double alpha,
405  const VectorValues& x, VectorValues& y) const {
406  for (const GaussianFactor::shared_ptr& f: *this)
407  f->multiplyHessianAdd(alpha, x, y);
408  }
409 
410  /* ************************************************************************* */
411  void GaussianFactorGraph::multiplyInPlace(const VectorValues& x, Errors& e) const {
412  multiplyInPlace(x, e.begin());
413  }
414 
415  /* ************************************************************************* */
416  void GaussianFactorGraph::multiplyInPlace(const VectorValues& x, const Errors::iterator& e) const {
417  Errors::iterator ei = e;
418  for (const GaussianFactor::shared_ptr& factor: *this) {
419  JacobianFactor::shared_ptr Ai = convertToJacobianFactorPtr(factor);
420  *ei = (*Ai)*x;
421  ei++;
422  }
423  }
424 
425  /* ************************************************************************* */
427  typedef JacobianFactor J;
428  for (const GaussianFactor::shared_ptr& factor: factors) {
429  J::shared_ptr jacobian(boost::dynamic_pointer_cast<J>(factor));
430  if (jacobian && jacobian->get_model() && jacobian->get_model()->isConstrained()) {
431  return true;
432  }
433  }
434  return false;
435  }
436 
437  /* ************************************************************************* */
438  // x += alpha*A'*e
439  void GaussianFactorGraph::transposeMultiplyAdd(double alpha, const Errors& e,
440  VectorValues& x) const {
441  // For each factor add the gradient contribution
442  Errors::const_iterator ei = e.begin();
443  for (const sharedFactor& factor: *this) {
444  JacobianFactor::shared_ptr Ai = convertToJacobianFactorPtr(factor);
445  Ai->transposeMultiplyAdd(alpha, *(ei++), x);
446  }
447  }
448 
450  //void residual(const GaussianFactorGraph& fg, const VectorValues &x, VectorValues &r) {
451  // Key i = 0 ;
452  // for (const GaussianFactor::shared_ptr& factor, fg) {
453  // JacobianFactor::shared_ptr Ai = convertToJacobianFactorPtr(factor);
454  // r[i] = Ai->getb();
455  // i++;
456  // }
457  // VectorValues Ax = VectorValues::SameStructure(r);
458  // multiply(fg,x,Ax);
459  // axpy(-1.0,Ax,r);
460  //}
461 
463  //void multiply(const GaussianFactorGraph& fg, const VectorValues &x, VectorValues &r) {
464  // r.setZero();
465  // Key i = 0;
466  // for (const GaussianFactor::shared_ptr& factor, fg) {
467  // JacobianFactor::shared_ptr Ai = convertToJacobianFactorPtr(factor);
468  // Vector &y = r[i];
469  // for (JacobianFactor::const_iterator j = Ai->begin(); j != Ai->end(); ++j) {
470  // y += Ai->getA(j) * x[*j];
471  // }
472  // ++i;
473  // }
474  //}
475 
476  /* ************************************************************************* */
477  VectorValues GaussianFactorGraph::transposeMultiply(const Errors& e) const
478  {
479  VectorValues x;
480  Errors::const_iterator ei = e.begin();
481  for (const sharedFactor& factor: *this) {
482  JacobianFactor::shared_ptr Ai = convertToJacobianFactorPtr(factor);
483  for (JacobianFactor::const_iterator j = Ai->begin(); j != Ai->end(); ++j) {
484  // Create the value as a zero vector if it does not exist.
485  pair<VectorValues::iterator, bool> xi = x.tryInsert(*j, Vector());
486  if(xi.second)
487  xi.first->second = Vector::Zero(Ai->getDim(j));
488  xi.first->second += Ai->getA(j).transpose() * *ei;
489  }
490  ++ ei;
491  }
492  return x;
493  }
494 
495  /* ************************************************************************* */
496  Errors GaussianFactorGraph::gaussianErrors(const VectorValues& x) const
497  {
498  Errors e;
499  for (const sharedFactor& factor: *this) {
500  JacobianFactor::shared_ptr Ai = convertToJacobianFactorPtr(factor);
501  e.push_back(Ai->error_vector(x));
502  }
503  return e;
504  }
505 
506  /* ************************************************************************* */
509  const Eliminate& function) const {
510  return optimize(function);
511  }
512 
513 } // namespace gtsam
Vector3_ operator*(const Double_ &s, const Vector3_ &v)
def step(data, isam, result, truth, currPoseIndex)
Definition: visual_isam.py:82
Scalar * y
#define gttic_(label)
Definition: timing.h:230
Global debugging flags.
static enum @843 ordering
boost::shared_ptr< This > shared_ptr
shared_ptr to this class
VectorValues & addInPlace_(const VectorValues &c)
int n
Eigen::MatrixXd Matrix
Definition: base/Matrix.h:43
Definition: Half.h:150
GaussianFactorGraph factors(list_of(factor1)(factor2)(factor3))
const Solve< LLT, Rhs > solve(const MatrixBase< Rhs > &b) const
Definition: LLT.h:144
IsDerived< DERIVEDFACTOR > push_back(boost::shared_ptr< DERIVEDFACTOR > factor)
Add a factor directly using a shared_ptr.
Definition: FactorGraph.h:166
void g(const string &key, int i)
Definition: testBTree.cpp:43
std::pair< iterator, bool > tryInsert(Key j, const Vector &value)
Definition: VectorValues.h:207
const SymmetricBlockMatrix & info() const
Return underlying information matrix.
Scalar Scalar int size
Definition: benchVecAdd.cpp:17
T negate(const T &x)
Definition: packetmath.cpp:27
#define gttic(label)
Definition: timing.h:280
Factor Graph Values.
Scalar EIGEN_BLAS_FUNC() dot(int *n, RealScalar *px, int *incx, RealScalar *py, int *incy)
Eigen::VectorXd Vector
Definition: Vector.h:38
size_t size() const
Definition: Values.h:236
Values result
constABlock::ConstColXpr constBVector
Standard Cholesky decomposition (LL^T) of a matrix and associated features.
Definition: LLT.h:56
boost::shared_ptr< GaussianFactor > sharedFactor
Shared pointer to a factor.
Definition: FactorGraph.h:98
boost::shared_ptr< This > shared_ptr
A shared_ptr to this class.
Point2(* f)(const Point3 &, OptionalJacobian< 2, 3 >)
RealScalar alpha
Contains the HessianFactor class, a general quadratic factor.
Array< double, 1, 3 > e(1./3., 0.5, 2.)
RealScalar s
Matrix augmentedJacobian() const override
Linear Factor Graph where all factors are Gaussians.
double dot(const VectorValues &v) const
JacobiRotation< float > J
boost::shared_ptr< This > shared_ptr
shared_ptr to this class
Vector xi
Definition: testPose2.cpp:150
static Symbol x0('x', 0)
traits
Definition: chartTesting.h:28
bool hasConstraints(const GaussianFactorGraph &factors)
std::vector< std::tuple< int, int, double > > SparseTriplets
Expression of a fixed-size or dynamic-size block.
Definition: Block.h:103
#define gttoc(label)
Definition: timing.h:281
int optimize(const SfmData &db, const NonlinearFactorGraph &graph, const Values &initial, bool separateCalibration=false)
Definition: timeSFMBAL.h:63
Eigen::SelfAdjointView< constBlock, Eigen::Upper > selfadjointView(DenseIndex I, DenseIndex J) const
Return the square sub-matrix that contains blocks(i:j, i:j).
boost::function< EliminationResult(const FactorGraphType &, const Ordering &)> Eliminate
The function type that does a single dense elimination step on a subgraph.
A Gaussian factor using the canonical parameters (information form)
const G double tol
Definition: Group.h:83
KeyVector::const_iterator const_iterator
Const iterator over keys.
Definition: Factor.h:67
const KeyVector keys
boost::shared_ptr< This > shared_ptr
shared_ptr to this class
set noclip points set clip one set noclip two set bar set border lt lw set xdata set ydata set zdata set x2data set y2data set boxwidth set dummy x
#define abs(x)
Definition: datatypes.h:17
Eigen::Matrix< double, Eigen::Dynamic, 1 > Vector
Gaussian Bayes Tree, the result of eliminating a GaussianJunctionTree.
std::uint64_t Key
Integer nonlinear key type.
Definition: types.h:61
std::ptrdiff_t j
Timing utilities.
Efficient incomplete Cholesky on rank-deficient matrices, todo: constrained Cholesky.


gtsam
Author(s):
autogenerated on Sat May 8 2021 02:42:05