Marginals.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 
19 #pragma once
20 
23 #include <gtsam/nonlinear/Values.h>
24 
25 namespace gtsam {
26 
27 class JointMarginal;
28 
32 class GTSAM_EXPORT Marginals {
33 
34 public:
35 
39  QR
40  };
41 
42 protected:
43 
48 
49 public:
50 
53 
59  Marginals(const NonlinearFactorGraph& graph, const Values& solution, Factorization factorization = CHOLESKY);
60 
67  Marginals(const NonlinearFactorGraph& graph, const Values& solution, const Ordering& ordering,
68  Factorization factorization = CHOLESKY);
69 
75  Marginals(const GaussianFactorGraph& graph, const Values& solution, Factorization factorization = CHOLESKY);
76 
83  Marginals(const GaussianFactorGraph& graph, const Values& solution, const Ordering& ordering,
84  Factorization factorization = CHOLESKY);
85 
92  Marginals(const GaussianFactorGraph& graph, const VectorValues& solution, Factorization factorization = CHOLESKY);
93 
100  Marginals(const GaussianFactorGraph& graph, const VectorValues& solution, const Ordering& ordering,
101  Factorization factorization = CHOLESKY);
102 
104  void print(const std::string& str = "Marginals: ", const KeyFormatter& keyFormatter = DefaultKeyFormatter) const;
105 
107  GaussianFactor::shared_ptr marginalFactor(Key variable) const;
108 
111  Matrix marginalInformation(Key variable) const;
112 
114  Matrix marginalCovariance(Key variable) const;
115 
117  JointMarginal jointMarginalCovariance(const KeyVector& variables) const;
118 
120  JointMarginal jointMarginalInformation(const KeyVector& variables) const;
121 
123  VectorValues optimize() const;
124 
125 protected:
126 
128  void computeBayesTree();
129 
131  void computeBayesTree(const Ordering& ordering);
132 };
133 
137 class GTSAM_EXPORT JointMarginal {
138 
139 protected:
143 
144 public:
147 
161  Matrix operator()(Key iVariable, Key jVariable) const {
162  const auto indexI = indices_.at(iVariable);
163  const auto indexJ = indices_.at(jVariable);
164  return blockMatrix_.block(indexI, indexJ);
165  }
166 
168  Matrix at(Key iVariable, Key jVariable) const {
169  return (*this)(iVariable, jVariable);
170  }
171 
173  Matrix fullMatrix() const {
174  return blockMatrix_.selfadjointView();
175  }
176 
178  void print(const std::string& s = "", const KeyFormatter& formatter = DefaultKeyFormatter) const;
179 
180 protected:
181  JointMarginal(const Matrix& fullMatrix, const std::vector<size_t>& dims, const KeyVector& keys) :
182  blockMatrix_(dims, fullMatrix), keys_(keys), indices_(Ordering(keys).invert()) {}
183 
184  friend class Marginals;
185 
186 };
187 
188 } /* namespace gtsam */
gtsam::JointMarginal::JointMarginal
JointMarginal(const Matrix &fullMatrix, const std::vector< size_t > &dims, const KeyVector &keys)
Definition: Marginals.h:181
gtsam::Marginals::Marginals
Marginals()
Default constructor only for wrappers.
Definition: Marginals.h:52
gtsam::JointMarginal::blockMatrix_
SymmetricBlockMatrix blockMatrix_
Definition: Marginals.h:140
s
RealScalar s
Definition: level1_cplx_impl.h:126
gtsam::Marginals::factorization_
Factorization factorization_
Definition: Marginals.h:46
keys
const KeyVector keys
Definition: testRegularImplicitSchurFactor.cpp:40
gtsam::optimize
Point3 optimize(const NonlinearFactorGraph &graph, const Values &values, Key landmarkKey)
Definition: triangulation.cpp:177
gtsam::JointMarginal::operator()
Matrix operator()(Key iVariable, Key jVariable) const
Definition: Marginals.h:161
gtsam::FastMap< Key, size_t >
gtsam::Marginals
Definition: Marginals.h:32
formatter
const KeyFormatter & formatter
Definition: treeTraversal-inst.h:204
gtsam::SymmetricBlockMatrix::selfadjointView
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:158
gtsam::Matrix
Eigen::MatrixXd Matrix
Definition: base/Matrix.h:39
gtsam::JointMarginal::at
Matrix at(Key iVariable, Key jVariable) const
Definition: Marginals.h:168
gtsam::KeyVector
FastVector< Key > KeyVector
Define collection type once and for all - also used in wrappers.
Definition: Key.h:92
gtsam::Marginals::CHOLESKY
@ CHOLESKY
Definition: Marginals.h:38
gtsam::DefaultKeyFormatter
KeyFormatter DefaultKeyFormatter
Assign default key formatter.
Definition: Key.cpp:30
gtsam::Marginals::values_
Values values_
Definition: Marginals.h:45
gtsam::JointMarginal::fullMatrix
Matrix fullMatrix() const
Definition: Marginals.h:173
gtsam::GaussianFactorGraph
Definition: GaussianFactorGraph.h:73
gtsam::print
void print(const Matrix &A, const string &s, ostream &stream)
Definition: Matrix.cpp:155
gtsam::JointMarginal::JointMarginal
JointMarginal()
Default constructor only for wrappers.
Definition: Marginals.h:146
gtsam::VectorValues
Definition: VectorValues.h:74
gtsam::KeyFormatter
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
gtsam::NonlinearFactorGraph
Definition: NonlinearFactorGraph.h:55
gtsam::GaussianFactor::shared_ptr
std::shared_ptr< This > shared_ptr
shared_ptr to this class
Definition: GaussianFactor.h:42
ordering
static enum @1096 ordering
str
Definition: pytypes.h:1524
gtsam::Marginals::graph_
GaussianFactorGraph graph_
Definition: Marginals.h:44
gtsam::JointMarginal::keys_
KeyVector keys_
Definition: Marginals.h:141
gtsam
traits
Definition: chartTesting.h:28
gtsam::GaussianBayesTree
Definition: GaussianBayesTree.h:49
NonlinearFactorGraph.h
Factor Graph consisting of non-linear factors.
gtsam::Values
Definition: Values.h:65
gtsam::JointMarginal::indices_
FastMap< Key, size_t > indices_
Definition: Marginals.h:142
GaussianBayesTree.h
Gaussian Bayes Tree, the result of eliminating a GaussianJunctionTree.
gtsam::SymmetricBlockMatrix::block
Matrix block(DenseIndex I, DenseIndex J) const
Definition: SymmetricBlockMatrix.cpp:60
gtsam::Marginals::Factorization
Factorization
Definition: Marginals.h:37
gtsam::Marginals::bayesTree_
GaussianBayesTree bayesTree_
Definition: Marginals.h:47
graph
NonlinearFactorGraph graph
Definition: doc/Code/OdometryExample.cpp:2
gtsam::JointMarginal
Definition: Marginals.h:137
gtsam::Key
std::uint64_t Key
Integer nonlinear key type.
Definition: types.h:97
gtsam::SymmetricBlockMatrix
Definition: SymmetricBlockMatrix.h:53
gtsam::Ordering
Definition: inference/Ordering.h:33
Values.h
A non-templated config holding any types of Manifold-group elements.


gtsam
Author(s):
autogenerated on Thu Jun 13 2024 03:03:32