GTSAM  4.0.2
C++ library for smoothing and mapping (SAM)
GaussianConditional.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 
18 // \callgraph
19 
20 #pragma once
21 
22 #include <gtsam/global_includes.h>
25 #include <gtsam/inference/Conditional-inst.h>
27 
28 #include <random> // for std::mt19937_64
29 
30 namespace gtsam {
31 
38  class GTSAM_EXPORT GaussianConditional :
39  public JacobianFactor,
40  public Conditional<JacobianFactor, GaussianConditional>
41  {
42  public:
44  typedef std::shared_ptr<This> shared_ptr;
47 
50 
53 
55  GaussianConditional(Key key, const Vector& d, const Matrix& R,
56  const SharedDiagonal& sigmas = SharedDiagonal());
57 
59  GaussianConditional(Key key, const Vector& d, const Matrix& R, Key parent1,
60  const Matrix& S,
61  const SharedDiagonal& sigmas = SharedDiagonal());
62 
64  GaussianConditional(Key key, const Vector& d, const Matrix& R, Key parent1,
65  const Matrix& S, Key parent2, const Matrix& T,
66  const SharedDiagonal& sigmas = SharedDiagonal());
67 
71  template<typename TERMS>
72  GaussianConditional(const TERMS& terms,
73  size_t nrFrontals, const Vector& d,
74  const SharedDiagonal& sigmas = SharedDiagonal());
75 
80  template<typename KEYS>
82  const KEYS& keys, size_t nrFrontals, const VerticalBlockMatrix& augmentedMatrix,
83  const SharedDiagonal& sigmas = SharedDiagonal());
84 
86  static GaussianConditional FromMeanAndStddev(Key key, const Vector& mu,
87  double sigma);
88 
90  static GaussianConditional FromMeanAndStddev(Key key, const Matrix& A,
91  Key parent, const Vector& b,
92  double sigma);
93 
96  static GaussianConditional FromMeanAndStddev(Key key, //
97  const Matrix& A1, Key parent1,
98  const Matrix& A2, Key parent2,
99  const Vector& b, double sigma);
100 
102  template<typename... Args>
103  static shared_ptr sharedMeanAndStddev(Args&&... args) {
104  return std::make_shared<This>(FromMeanAndStddev(std::forward<Args>(args)...));
105  }
106 
114  template<typename ITERATOR>
115  static shared_ptr Combine(ITERATOR firstConditional, ITERATOR lastConditional);
116 
120 
122  void print(
123  const std::string& = "GaussianConditional",
124  const KeyFormatter& formatter = DefaultKeyFormatter) const override;
125 
127  bool equals(const GaussianFactor&cg, double tol = 1e-9) const override;
128 
132 
137  double logNormalizationConstant() const override;
138 
146  double logProbability(const VectorValues& x) const;
147 
153  double evaluate(const VectorValues& x) const;
154 
156  double operator()(const VectorValues& x) const {
157  return evaluate(x);
158  }
159 
173  VectorValues solve(const VectorValues& parents) const;
174 
175  VectorValues solveOtherRHS(const VectorValues& parents, const VectorValues& rhs) const;
176 
178  void solveTransposeInPlace(VectorValues& gy) const;
179 
181  JacobianFactor::shared_ptr likelihood(
182  const VectorValues& frontalValues) const;
183 
185  JacobianFactor::shared_ptr likelihood(const Vector& frontal) const;
186 
193  VectorValues sample(std::mt19937_64* rng) const;
194 
202  VectorValues sample(const VectorValues& parentsValues,
203  std::mt19937_64* rng) const;
204 
206  VectorValues sample() const;
207 
209  VectorValues sample(const VectorValues& parentsValues) const;
210 
214 
216  constABlock R() const { return Ab_.range(0, nrFrontals()); }
217 
219  constABlock S() const { return Ab_.range(nrFrontals(), size()); }
220 
222  constABlock S(const_iterator it) const { return BaseFactor::getA(it); }
223 
225  const constBVector d() const { return BaseFactor::getb(); }
226 
238  inline double determinant() const { return exp(logDeterminant()); }
239 
251  double logDeterminant() const;
252 
256 
261  double logProbability(const HybridValues& x) const override;
262 
267  double evaluate(const HybridValues& x) const override;
268 
269  using Conditional::operator(); // Expose evaluate(const HybridValues&) method..
270  using JacobianFactor::error; // Expose error(const HybridValues&) method..
271 
273 
274  private:
275 #ifdef GTSAM_ENABLE_BOOST_SERIALIZATION
276 
277  friend class boost::serialization::access;
278  template<class Archive>
279  void serialize(Archive & ar, const unsigned int /*version*/) {
280  ar & BOOST_SERIALIZATION_BASE_OBJECT_NVP(BaseFactor);
281  ar & BOOST_SERIALIZATION_BASE_OBJECT_NVP(BaseConditional);
282  }
283 #endif
284  }; // GaussianConditional
285 
287 template<>
288 struct traits<GaussianConditional> : public Testable<GaussianConditional> {};
289 
290 } // \ namespace gtsam
291 
293 
Definition: HybridValues.h:38
std::shared_ptr< This > shared_ptr
shared_ptr to this class
Definition: GaussianConditional.h:44
Base class for conditional densities.
Conditional< BaseFactor, This > BaseConditional
Typedef to our conditional base class.
Definition: GaussianConditional.h:46
std::string serialize(const T &input)
serializes to a string
Definition: serialization.h:113
Definition: VerticalBlockMatrix.h:42
Definition: Testable.h:152
Definition: Group.h:43
JacobianFactor BaseFactor
Typedef to our factor base class.
Definition: GaussianConditional.h:45
double determinant() const
Compute the determinant of the R matrix.
Definition: GaussianConditional.h:238
Conditional Gaussian Base class.
const constBVector d() const
Definition: GaussianConditional.h:225
GaussianConditional()
Definition: GaussianConditional.h:52
Definition: VectorValues.h:74
constABlock S() const
Definition: GaussianConditional.h:219
constABlock R() const
Definition: GaussianConditional.h:216
Included from all GTSAM files.
Definition: GaussianFactor.h:38
Definition: GaussianConditional.h:38
Factor Graph Values.
Definition: Testable.h:112
GaussianConditional This
Typedef to this class.
Definition: GaussianConditional.h:43
GTSAM_EXPORT void print(const Matrix &A, const std::string &s, std::ostream &stream)
Definition: Conditional.h:61
std::shared_ptr< This > shared_ptr
shared_ptr to this class
Definition: JacobianFactor.h:97
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
static shared_ptr sharedMeanAndStddev(Args &&... args)
Create shared pointer by forwarding arguments to fromMeanAndStddev.
Definition: GaussianConditional.h:103
Definition: JacobianFactor.h:91
Definition: chartTesting.h:28
constABlock S(const_iterator it) const
Definition: GaussianConditional.h:222
double operator()(const VectorValues &x) const
Evaluate probability density, sugar.
Definition: GaussianConditional.h:156
KeyVector::const_iterator const_iterator
Const iterator over keys.
Definition: Factor.h:82
std::uint64_t Key
Integer nonlinear key type.
Definition: types.h:102