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 
24 #include <cstddef>
27 #include <gtsam/linear/Errors.h> // Included here instead of fw-declared so we can use Errors::iterator
32 
33 namespace gtsam {
34 
35  // Forward declarations
36  class GaussianFactorGraph;
37  class GaussianFactor;
38  class GaussianConditional;
39  class GaussianBayesNet;
40  class GaussianEliminationTree;
41  class GaussianBayesTree;
42  class GaussianJunctionTree;
43 
44  /* ************************************************************************* */
46  {
54  static std::pair<std::shared_ptr<ConditionalType>, std::shared_ptr<FactorType> >
56  DefaultEliminate(const FactorGraphType& factors, const Ordering& keys) {
57  return EliminatePreferCholesky(factors, keys); }
60  const FactorGraphType& graph,
61  std::optional<std::reference_wrapper<const VariableIndex>> variableIndex) {
62  return Ordering::Colamd((*variableIndex).get());
63  }
64  };
65 
66  /* ************************************************************************* */
73  class GTSAM_EXPORT GaussianFactorGraph :
74  public FactorGraph<GaussianFactor>,
75  public EliminateableFactorGraph<GaussianFactorGraph>
76  {
77  public:
78 
82  typedef std::shared_ptr<This> shared_ptr;
83 
86 
89 
95  GaussianFactorGraph(std::initializer_list<sharedFactor> factors) : Base(factors) {}
96 
97 
99  template<typename ITERATOR>
100  GaussianFactorGraph(ITERATOR firstFactor, ITERATOR lastFactor) : Base(firstFactor, lastFactor) {}
101 
103  template<class CONTAINER>
104  explicit GaussianFactorGraph(const CONTAINER& factors) : Base(factors) {}
105 
107  template<class DERIVEDFACTOR>
109 
113 
114  bool equals(const This& fg, double tol = 1e-9) const;
115 
117 
119  friend bool operator==(const GaussianFactorGraph& lhs,
120  const GaussianFactorGraph& rhs) {
121  return lhs.isEqual(rhs);
122  }
123 
125  void add(const GaussianFactor& factor) { push_back(factor.clone()); }
126 
128  void add(const sharedFactor& factor) { push_back(factor); }
129 
131  void add(const Vector& b) {
132  add(JacobianFactor(b)); }
133 
135  void add(Key key1, const Matrix& A1,
136  const Vector& b, const SharedDiagonal& model = SharedDiagonal()) {
137  add(JacobianFactor(key1,A1,b,model)); }
138 
140  void add(Key key1, const Matrix& A1,
141  Key key2, const Matrix& A2,
142  const Vector& b, const SharedDiagonal& model = SharedDiagonal()) {
143  add(JacobianFactor(key1,A1,key2,A2,b,model)); }
144 
146  void add(Key key1, const Matrix& A1,
147  Key key2, const Matrix& A2,
148  Key key3, const Matrix& A3,
149  const Vector& b, const SharedDiagonal& model = SharedDiagonal()) {
150  add(JacobianFactor(key1,A1,key2,A2,key3,A3,b,model)); }
151 
153  template<class TERMS>
154  void add(const TERMS& terms, const Vector &b, const SharedDiagonal& model = SharedDiagonal()) {
155  add(JacobianFactor(terms,b,model)); }
156 
161  typedef KeySet Keys;
162  Keys keys() const;
163 
164  /* return a map of (Key, dimension) */
165  std::map<Key, size_t> getKeyDimMap() const;
166 
168  double error(const VectorValues& x) const;
169 
171  double probPrime(const VectorValues& c) const;
172 
178  virtual GaussianFactorGraph clone() const;
179 
184  virtual GaussianFactorGraph::shared_ptr cloneToPtr() const;
185 
192  GaussianFactorGraph negate() const;
193 
196 
207  std::vector<std::tuple<int, int, double> > sparseJacobian(
208  const Ordering& ordering, size_t& nrows, size_t& ncols) const;
209 
211  std::vector<std::tuple<int, int, double> > sparseJacobian() const;
212 
219  Matrix sparseJacobian_() const;
220 
228  Matrix augmentedJacobian(const Ordering& ordering) const;
229 
237  Matrix augmentedJacobian() const;
238 
246  std::pair<Matrix,Vector> jacobian(const Ordering& ordering) const;
247 
255  std::pair<Matrix,Vector> jacobian() const;
256 
268  Matrix augmentedHessian(const Ordering& ordering) const;
269 
281  Matrix augmentedHessian() const;
282 
289  std::pair<Matrix,Vector> hessian(const Ordering& ordering) const;
290 
297  std::pair<Matrix,Vector> hessian() const;
298 
300  virtual VectorValues hessianDiagonal() const;
301 
303  virtual std::map<Key,Matrix> hessianBlockDiagonal() const;
304 
310  const Eliminate& function = EliminationTraitsType::DefaultEliminate) const;
311 
317  const Eliminate& function = EliminationTraitsType::DefaultEliminate) const;
318 
322  VectorValues optimizeDensely() const;
323 
333  VectorValues gradient(const VectorValues& x0) const;
334 
342  virtual VectorValues gradientAtZero() const;
343 
368  VectorValues optimizeGradientSearch() const;
369 
371  VectorValues transposeMultiply(const Errors& e) const;
372 
374  void transposeMultiplyAdd(double alpha, const Errors& e, VectorValues& x) const;
375 
377  Errors gaussianErrors(const VectorValues& x) const;
378 
380  Errors operator*(const VectorValues& x) const;
381 
383  void multiplyHessianAdd(double alpha, const VectorValues& x,
384  VectorValues& y) const;
385 
387  void multiplyInPlace(const VectorValues& x, Errors& e) const;
388 
390  void multiplyInPlace(const VectorValues& x, const Errors::iterator& e) const;
391 
392  void printErrors(
393  const VectorValues& x,
394  const std::string& str = "GaussianFactorGraph: ",
395  const KeyFormatter& keyFormatter = DefaultKeyFormatter,
396  const std::function<bool(const Factor* /*factor*/,
397  double /*whitenedError*/, size_t /*index*/)>&
398  printCondition =
399  [](const Factor*, double, size_t) { return true; }) const;
401 
402  private:
403 #ifdef GTSAM_ENABLE_BOOST_SERIALIZATION
404 
405  friend class boost::serialization::access;
406  template<class ARCHIVE>
407  void serialize(ARCHIVE & ar, const unsigned int /*version*/) {
408  ar & BOOST_SERIALIZATION_BASE_OBJECT_NVP(Base);
409  }
410 #endif
411  };
412 
417  GTSAM_EXPORT bool hasConstraints(const GaussianFactorGraph& factors);
418 
419  /****** Linear Algebra Operations ******/
420 
422  //GTSAM_EXPORT void residual(const GaussianFactorGraph& fg, const VectorValues &x, VectorValues &r);
423  //GTSAM_EXPORT void multiply(const GaussianFactorGraph& fg, const VectorValues &x, VectorValues &r);
424 
426 template<>
427 struct traits<GaussianFactorGraph> : public Testable<GaussianFactorGraph> {
428 };
429 
430 } // \ namespace gtsam
static Matrix A1
const Symbol key3('v', 3)
GaussianFactorGraph FactorGraphType
Type of the factor graph (e.g. GaussianFactorGraph)
vector of errors
Scalar * y
std::string serialize(const T &input)
serializes to a string
friend bool operator==(const GaussianFactorGraph &lhs, const GaussianFactorGraph &rhs)
Check exact equality.
Variable elimination algorithms for factor graphs.
noiseModel::Diagonal::shared_ptr model
GaussianEliminationTree EliminationTreeType
Type of elimination tree.
bool isEqual(const FactorGraph &other) const
Check exact equality of the factor pointers. Useful for derived ==.
Definition: FactorGraph.h:138
Eigen::MatrixXd Matrix
Definition: base/Matrix.h:39
A factor with a quadratic error function - a Gaussian.
Scalar Scalar * c
Definition: benchVecAdd.cpp:17
const GaussianFactorGraph factors
std::pair< std::shared_ptr< GaussianConditional >, std::shared_ptr< GaussianFactor > > EliminatePreferCholesky(const GaussianFactorGraph &factors, const Ordering &keys)
GaussianFactorGraph(const FactorGraph< DERIVEDFACTOR > &graph)
virtual GaussianFactor::shared_ptr clone() const =0
static std::pair< std::shared_ptr< ConditionalType >, std::shared_ptr< FactorType > > DefaultEliminate(const FactorGraphType &factors, const Ordering &keys)
The default dense elimination function.
Point2 operator*(double s, const Point2 &p)
multiply with scalar
Definition: Point2.h:52
NonlinearFactorGraph graph
static const KeyFormatter DefaultKeyFormatter
Definition: Key.h:43
static enum @1107 ordering
Continuous multi-dimensional vectors for.
static Ordering Colamd(const FACTOR_GRAPH &graph)
GaussianFactorGraph(std::initializer_list< sharedFactor > factors)
static Ordering DefaultOrderingFunc(const FactorGraphType &graph, std::optional< std::reference_wrapper< const VariableIndex >> variableIndex)
The default ordering generation function.
Factor Graph Values.
void add(const GaussianFactor &factor)
GaussianFactorGraph(ITERATOR firstFactor, ITERATOR lastFactor)
Eigen::VectorXd Vector
Definition: Vector.h:38
const Symbol key1('v', 1)
EliminateableFactorGraph< This > BaseEliminateable
Typedef to base elimination class.
void add(Key key1, const Matrix &A1, Key key2, const Matrix &A2, Key key3, const Matrix &A3, const Vector &b, const SharedDiagonal &model=SharedDiagonal())
Definition: pytypes.h:1403
Point3 optimize(const NonlinearFactorGraph &graph, const Values &values, Key landmarkKey)
RealScalar alpha
Contains the HessianFactor class, a general quadratic factor.
Array< double, 1, 3 > e(1./3., 0.5, 2.)
std::shared_ptr< This > shared_ptr
shared_ptr to this class
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
const G & b
Definition: Group.h:86
noiseModel::Diagonal::shared_ptr SharedDiagonal
Definition: NoiseModel.h:743
GaussianFactorGraph(const CONTAINER &factors)
traits
Definition: chartTesting.h:28
GaussianConditional ConditionalType
Type of conditionals from elimination.
static Symbol x0('x', 0)
GaussianBayesTree BayesTreeType
Type of Bayes tree.
void add(Key key1, const Matrix &A1, Key key2, const Matrix &A2, const Vector &b, const SharedDiagonal &model=SharedDiagonal())
Factor Graph Base Class.
void add(Key key1, const Matrix &A1, const Vector &b, const SharedDiagonal &model=SharedDiagonal())
GaussianBayesNet BayesNetType
Type of Bayes net from sequential elimination.
GaussianFactor FactorType
Type of factors in factor graph.
GaussianFactorGraph This
Typedef to this class.
graph add(PriorFactor< Pose2 >(1, priorMean, priorNoise))
static double error
Definition: testRot3.cpp:37
FactorGraph< GaussianFactor > Base
Typedef to base factor graph type.
void add(const sharedFactor &factor)
const G double tol
Definition: Group.h:86
std::shared_ptr< GaussianFactor > sharedFactor
Shared pointer to a factor.
Definition: FactorGraph.h:105
void add(const TERMS &terms, const Vector &b, const SharedDiagonal &model=SharedDiagonal())
const KeyVector keys
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
bool hasConstraints(const GaussianFactorGraph &factors)
std::uint64_t Key
Integer nonlinear key type.
Definition: types.h:102
const Symbol key2('v', 2)


gtsam
Author(s):
autogenerated on Tue Jul 4 2023 02:34:15