gtsam 4.1.1
gtsam
GaussianFactorGraph.h
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
22#pragma once
23
26#include <gtsam/linear/Errors.h> // Included here instead of fw-declared so we can use Errors::iterator
31
32namespace gtsam {
33
34 // Forward declarations
35 class GaussianFactorGraph;
36 class GaussianFactor;
37 class GaussianConditional;
38 class GaussianBayesNet;
39 class GaussianEliminationTree;
40 class GaussianBayesTree;
41 class GaussianJunctionTree;
42
43 /* ************************************************************************* */
45 {
54 static std::pair<boost::shared_ptr<ConditionalType>, boost::shared_ptr<FactorType> >
55 DefaultEliminate(const FactorGraphType& factors, const Ordering& keys) {
56 return EliminatePreferCholesky(factors, keys); }
57 };
58
59 /* ************************************************************************* */
66 class GTSAM_EXPORT GaussianFactorGraph :
67 public FactorGraph<GaussianFactor>,
68 public EliminateableFactorGraph<GaussianFactorGraph>
69 {
70 public:
71
75 typedef boost::shared_ptr<This> shared_ptr;
76
79
81 template<typename ITERATOR>
82 GaussianFactorGraph(ITERATOR firstFactor, ITERATOR lastFactor) : Base(firstFactor, lastFactor) {}
83
85 template<class CONTAINER>
86 explicit GaussianFactorGraph(const CONTAINER& factors) : Base(factors) {}
87
89 template<class DERIVEDFACTOR>
91
94
97
98 bool equals(const This& fg, double tol = 1e-9) const;
99
101
103 void add(const GaussianFactor& factor) { push_back(factor.clone()); }
104
106 void add(const sharedFactor& factor) { push_back(factor); }
107
109 void add(const Vector& b) {
110 add(JacobianFactor(b)); }
111
113 void add(Key key1, const Matrix& A1,
114 const Vector& b, const SharedDiagonal& model = SharedDiagonal()) {
115 add(JacobianFactor(key1,A1,b,model)); }
116
118 void add(Key key1, const Matrix& A1,
119 Key key2, const Matrix& A2,
120 const Vector& b, const SharedDiagonal& model = SharedDiagonal()) {
121 add(JacobianFactor(key1,A1,key2,A2,b,model)); }
122
124 void add(Key key1, const Matrix& A1,
125 Key key2, const Matrix& A2,
126 Key key3, const Matrix& A3,
127 const Vector& b, const SharedDiagonal& model = SharedDiagonal()) {
128 add(JacobianFactor(key1,A1,key2,A2,key3,A3,b,model)); }
129
131 template<class TERMS>
132 void add(const TERMS& terms, const Vector &b, const SharedDiagonal& model = SharedDiagonal()) {
133 add(JacobianFactor(terms,b,model)); }
134
139 typedef KeySet Keys;
140 Keys keys() const;
141
142 /* return a map of (Key, dimension) */
143 std::map<Key, size_t> getKeyDimMap() const;
144
146 double error(const VectorValues& x) const {
147 double total_error = 0.;
148 for(const sharedFactor& factor: *this){
149 if(factor)
150 total_error += factor->error(x);
151 }
152 return total_error;
153 }
154
156 double probPrime(const VectorValues& c) const {
157 return exp(-0.5 * error(c));
158 }
159
165 virtual GaussianFactorGraph clone() const;
166
171 virtual GaussianFactorGraph::shared_ptr cloneToPtr() const;
172
179 GaussianFactorGraph negate() const;
180
183
194 std::vector<std::tuple<int, int, double> > sparseJacobian(
195 const Ordering& ordering, size_t& nrows, size_t& ncols) const;
196
198 std::vector<std::tuple<int, int, double> > sparseJacobian() const;
199
206 Matrix sparseJacobian_() const;
207
215 Matrix augmentedJacobian(const Ordering& ordering) const;
216
224 Matrix augmentedJacobian() const;
225
233 std::pair<Matrix,Vector> jacobian(const Ordering& ordering) const;
234
242 std::pair<Matrix,Vector> jacobian() const;
243
255 Matrix augmentedHessian(const Ordering& ordering) const;
256
268 Matrix augmentedHessian() const;
269
276 std::pair<Matrix,Vector> hessian(const Ordering& ordering) const;
277
284 std::pair<Matrix,Vector> hessian() const;
285
287 virtual VectorValues hessianDiagonal() const;
288
290 virtual std::map<Key,Matrix> hessianBlockDiagonal() const;
291
297 const Eliminate& function = EliminationTraitsType::DefaultEliminate) const;
298
304 const Eliminate& function = EliminationTraitsType::DefaultEliminate) const;
305
309 VectorValues optimizeDensely() const;
310
320 VectorValues gradient(const VectorValues& x0) const;
321
329 virtual VectorValues gradientAtZero() const;
330
355 VectorValues optimizeGradientSearch() const;
356
358 VectorValues transposeMultiply(const Errors& e) const;
359
361 void transposeMultiplyAdd(double alpha, const Errors& e, VectorValues& x) const;
362
364 Errors gaussianErrors(const VectorValues& x) const;
365
367 Errors operator*(const VectorValues& x) const;
368
370 void multiplyHessianAdd(double alpha, const VectorValues& x,
371 VectorValues& y) const;
372
374 void multiplyInPlace(const VectorValues& x, Errors& e) const;
375
377 void multiplyInPlace(const VectorValues& x, const Errors::iterator& e) const;
378
379 void printErrors(
380 const VectorValues& x,
381 const std::string& str = "GaussianFactorGraph: ",
382 const KeyFormatter& keyFormatter = DefaultKeyFormatter,
383 const std::function<bool(const Factor* /*factor*/,
384 double /*whitenedError*/, size_t /*index*/)>&
385 printCondition =
386 [](const Factor*, double, size_t) { return true; }) const;
388
389 private:
391 friend class boost::serialization::access;
392 template<class ARCHIVE>
393 void serialize(ARCHIVE & ar, const unsigned int /*version*/) {
394 ar & BOOST_SERIALIZATION_BASE_OBJECT_NVP(Base);
395 }
396
397 public:
398
399#ifdef GTSAM_ALLOW_DEPRECATED_SINCE_V41
401 VectorValues optimize(boost::none_t,
402 const Eliminate& function =
403 EliminationTraitsType::DefaultEliminate) const {
404 return optimize(function);
405 }
406#endif
407
408 };
409
414 GTSAM_EXPORT bool hasConstraints(const GaussianFactorGraph& factors);
415
416 /****** Linear Algebra Opeations ******/
417
419 //GTSAM_EXPORT void residual(const GaussianFactorGraph& fg, const VectorValues &x, VectorValues &r);
420 //GTSAM_EXPORT void multiply(const GaussianFactorGraph& fg, const VectorValues &x, VectorValues &r);
421
423template<>
424struct traits<GaussianFactorGraph> : public Testable<GaussianFactorGraph> {
425};
426
427} // \ namespace gtsam
Factor Graph Base Class.
Variable elimination algorithms for factor graphs.
Contains the HessianFactor class, a general quadratic factor.
vector of errors
A factor with a quadratic error function - a Gaussian.
Factor Graph Values.
Global functions in a separate testing namespace.
Definition: chartTesting.h:28
std::string serialize(const T &input)
serializes to a string
Definition: serialization.h:112
bool hasConstraints(const GaussianFactorGraph &factors)
Evaluates whether linear factors have any constrained noise models.
Definition: GaussianFactorGraph.cpp:426
Point3 optimize(const NonlinearFactorGraph &graph, const Values &values, Key landmarkKey)
Optimize for triangulation.
Definition: triangulation.cpp:73
Point2 operator*(double s, const Point2 &p)
multiply with scalar
Definition: Point2.h:47
std::uint64_t Key
Integer nonlinear key type.
Definition: types.h:69
std::function< std::string(Key)> KeyFormatter
Typedef for a function to format a key, i.e. to convert it to a string.
Definition: Key.h:35
A manifold defines a space in which there is a notion of a linear tangent space that can be centered ...
Definition: concepts.h:30
Template to create a binary predicate.
Definition: Testable.h:111
A helper that implements the traits interface for GTSAM types.
Definition: Testable.h:151
A factor graph is a bipartite graph with factor nodes connected to variable nodes.
Definition: FactorGraph.h:93
boost::shared_ptr< GaussianFactor > sharedFactor
Shared pointer to a factor.
Definition: FactorGraph.h:97
Traits class for eliminateable factor graphs, specifies the types that result from elimination,...
Definition: EliminateableFactorGraph.h:36
EliminateableFactorGraph is a base class for factor graphs that contains elimination algorithms.
Definition: EliminateableFactorGraph.h:57
This is the base class for all factor types.
Definition: Factor.h:56
Definition: Ordering.h:34
vector of errors
Definition: Errors.h:34
A Bayes net made from linear-Gaussian densities.
Definition: GaussianBayesNet.h:31
A Bayes tree representing a Gaussian density.
Definition: GaussianBayesTree.h:52
A conditional Gaussian functions as the node in a Bayes network It has a set of parents y,...
Definition: GaussianConditional.h:39
Definition: GaussianEliminationTree.h:29
An abstract virtual base class for JacobianFactor and HessianFactor.
Definition: GaussianFactor.h:39
virtual GaussianFactor::shared_ptr clone() const =0
Clone a factor (make a deep copy)
GaussianBayesTree BayesTreeType
Type of Bayes tree.
Definition: GaussianFactorGraph.h:51
GaussianConditional ConditionalType
Type of conditionals from elimination.
Definition: GaussianFactorGraph.h:48
GaussianFactor FactorType
Type of factors in factor graph.
Definition: GaussianFactorGraph.h:46
GaussianEliminationTree EliminationTreeType
Type of elimination tree.
Definition: GaussianFactorGraph.h:50
GaussianFactorGraph FactorGraphType
Type of the factor graph (e.g. GaussianFactorGraph)
Definition: GaussianFactorGraph.h:47
GaussianBayesNet BayesNetType
Type of Bayes net from sequential elimination.
Definition: GaussianFactorGraph.h:49
GaussianJunctionTree JunctionTreeType
Type of Junction tree.
Definition: GaussianFactorGraph.h:52
static std::pair< boost::shared_ptr< ConditionalType >, boost::shared_ptr< FactorType > > DefaultEliminate(const FactorGraphType &factors, const Ordering &keys)
The default dense elimination function.
Definition: GaussianFactorGraph.h:55
A Linear Factor Graph is a factor graph where all factors are Gaussian, i.e.
Definition: GaussianFactorGraph.h:69
EliminateableFactorGraph< This > BaseEliminateable
Typedef to base elimination class.
Definition: GaussianFactorGraph.h:74
boost::shared_ptr< This > shared_ptr
shared_ptr to this class
Definition: GaussianFactorGraph.h:75
void add(const TERMS &terms, const Vector &b, const SharedDiagonal &model=SharedDiagonal())
Add an n-ary factor.
Definition: GaussianFactorGraph.h:132
GaussianFactorGraph()
Default constructor.
Definition: GaussianFactorGraph.h:78
void add(const GaussianFactor &factor)
Add a factor by value - makes a copy.
Definition: GaussianFactorGraph.h:103
double probPrime(const VectorValues &c) const
Unnormalized probability.
Definition: GaussianFactorGraph.h:156
void add(Key key1, const Matrix &A1, Key key2, const Matrix &A2, Key key3, const Matrix &A3, const Vector &b, const SharedDiagonal &model=SharedDiagonal())
Add a ternary factor.
Definition: GaussianFactorGraph.h:124
void add(const sharedFactor &factor)
Add a factor by pointer - stores pointer without copying the factor.
Definition: GaussianFactorGraph.h:106
GaussianFactorGraph This
Typedef to this class.
Definition: GaussianFactorGraph.h:72
KeySet Keys
Return the set of variables involved in the factors (computes a set union).
Definition: GaussianFactorGraph.h:139
double error(const VectorValues &x) const
unnormalized error
Definition: GaussianFactorGraph.h:146
void add(const Vector &b)
Add a null factor.
Definition: GaussianFactorGraph.h:109
void add(Key key1, const Matrix &A1, const Vector &b, const SharedDiagonal &model=SharedDiagonal())
Add a unary factor.
Definition: GaussianFactorGraph.h:113
virtual ~GaussianFactorGraph()
Virtual destructor.
Definition: GaussianFactorGraph.h:93
void add(Key key1, const Matrix &A1, Key key2, const Matrix &A2, const Vector &b, const SharedDiagonal &model=SharedDiagonal())
Add a binary factor.
Definition: GaussianFactorGraph.h:118
GaussianFactorGraph(ITERATOR firstFactor, ITERATOR lastFactor)
Construct from iterator over factors.
Definition: GaussianFactorGraph.h:82
FactorGraph< GaussianFactor > Base
Typedef to base factor graph type.
Definition: GaussianFactorGraph.h:73
GaussianFactorGraph(const FactorGraph< DERIVEDFACTOR > &graph)
Implicit copy/downcast constructor to override explicit template container constructor.
Definition: GaussianFactorGraph.h:90
GaussianFactorGraph(const CONTAINER &factors)
Construct from container of factors (shared_ptr or plain objects)
Definition: GaussianFactorGraph.h:86
Definition: GaussianJunctionTree.h:37
A Gaussian factor in the squared-error form.
Definition: JacobianFactor.h:91
This class represents a collection of vector-valued variables associated each with a unique integer i...
Definition: VectorValues.h:74