MUQ  0.4.3
Gaussian.h
Go to the documentation of this file.
1 #ifndef GAUSSIAN_H_
2 #define GAUSSIAN_H_
3 
5 
6 namespace muq {
7  namespace Modeling {
8 
14  class Gaussian : public GaussianBase {
15  public:
16 
18  enum Mode {
19 
22 
24  Precision
25  };
26 
27  enum ExtraInputs {
28  None = 1 << 0,
29  Mean = 1 << 1,
30  DiagCovariance = 1 << 2,
31  DiagPrecision = 1 << 3,
32  FullCovariance = 1 << 4,
33  FullPrecision = 1 << 5,
34  };
35  typedef uint8_t InputMask;
36 
37  Gaussian(unsigned int dim,
38  InputMask extraInputs = ExtraInputs::None);
39 
41 
44  Gaussian(Eigen::VectorXd const& mu,
45  InputMask extraInputs = ExtraInputs::None);
46 
48 
53  Gaussian(Eigen::VectorXd const& mu,
54  Eigen::MatrixXd const& obj,
55  Gaussian::Mode mode = Gaussian::Mode::Covariance,
56  InputMask extraInputs = ExtraInputs::None);
57 
58 
59  virtual ~Gaussian() = default;
60 
61  Mode GetMode() const{return mode;};
62 
63  virtual Eigen::MatrixXd ApplyCovariance(Eigen::Ref<const Eigen::MatrixXd> const& x) const override;
64  virtual Eigen::MatrixXd ApplyPrecision(Eigen::Ref<const Eigen::MatrixXd> const& x) const override;
65 
66  virtual Eigen::MatrixXd ApplyCovSqrt(Eigen::Ref<const Eigen::MatrixXd> const& x) const override;
67  virtual Eigen::MatrixXd ApplyPrecSqrt(Eigen::Ref<const Eigen::MatrixXd> const& x) const override;
68 
70 
73  Eigen::MatrixXd GetCovariance() const;
74  Eigen::MatrixXd GetPrecision() const;
75 
77 
80  void SetCovariance(Eigen::MatrixXd const& newCov);
81 
83 
86  void SetPrecision(Eigen::MatrixXd const& newPrec);
87 
89 
91  std::shared_ptr<Gaussian> Condition(Eigen::MatrixXd const& obsMat,
92  Eigen::VectorXd const& data,
93  Eigen::MatrixXd const& obsCov) const;
94 
95  void ResetHyperparameters(ref_vector<Eigen::VectorXd> const& params) override;
96 
97  virtual double LogDeterminant() const override{return logDet;};
98 
99  protected:
100 
102 
105  void ComputeNormalization();
106 
107  static Eigen::VectorXi GetExtraSizes(unsigned dim, InputMask extraInputs);
108 
109  static Gaussian::Mode ModeFromExtras(InputMask extraInputs);
110  static void CheckInputTypes(InputMask extraInputs, Mode mode);
111 
114 
117 
118  // Space to store either the covariance or precision matrix (depending on mode)
119  Eigen::MatrixXd covPrec;
120 
121  // Space to tore the matrix square root of the covariance or precision
122  Eigen::LLT<Eigen::MatrixXd> sqrtCovPrec;
123 
125  double logDet;
126 
127  };
128  } // namespace Modeling
129 } // namespace muq
130 
131 #endif
Defines an abstract Gaussian class.@seealso Gaussian.
Definition: GaussianBase.h:15
Gaussian::InputMask inputTypes
What form do the extra inputs take? Just the mean, or the mean and covariance?
Definition: Gaussian.h:116
Eigen::MatrixXd GetPrecision() const
Definition: Gaussian.cpp:108
Gaussian::Mode mode
Have we specified the covariance or the precision.
Definition: Gaussian.h:113
void SetCovariance(Eigen::MatrixXd const &newCov)
Set the covariance matrix.
Definition: Gaussian.cpp:182
void SetPrecision(Eigen::MatrixXd const &newPrec)
Set the precision matrix.
Definition: Gaussian.cpp:199
Mode GetMode() const
Definition: Gaussian.h:61
static Gaussian::Mode ModeFromExtras(InputMask extraInputs)
Definition: Gaussian.cpp:56
Gaussian::InputMask GetInputTypes() const
Definition: Gaussian.h:88
static void CheckInputTypes(InputMask extraInputs, Mode mode)
Definition: Gaussian.cpp:48
virtual Eigen::MatrixXd ApplyCovSqrt(Eigen::Ref< const Eigen::MatrixXd > const &x) const override
Definition: Gaussian.cpp:215
void ComputeNormalization()
Compute the distribution's scaling constant.
Definition: Gaussian.cpp:126
Eigen::LLT< Eigen::MatrixXd > sqrtCovPrec
Definition: Gaussian.h:122
Eigen::MatrixXd covPrec
Definition: Gaussian.h:119
std::shared_ptr< Gaussian > Condition(Eigen::MatrixXd const &obsMat, Eigen::VectorXd const &data, Eigen::MatrixXd const &obsCov) const
Returns a new Gaussian distribution conditioned on a linear observation.
Definition: Gaussian.cpp:282
static Eigen::VectorXi GetExtraSizes(unsigned dim, InputMask extraInputs)
Definition: Gaussian.cpp:65
Mode
Are we specifying the mean, covariance matrix, or precision matrix.
Definition: Gaussian.h:18
@ Precision
We are specifying the precision.
Definition: Gaussian.h:24
@ Covariance
We are specifying the covariance.
Definition: Gaussian.h:21
Gaussian(unsigned int dim, InputMask extraInputs=ExtraInputs::None)
Definition: Gaussian.cpp:8
virtual double LogDeterminant() const override
Definition: Gaussian.h:97
virtual Eigen::MatrixXd ApplyPrecSqrt(Eigen::Ref< const Eigen::MatrixXd > const &x) const override
Definition: Gaussian.cpp:231
double logDet
The log determinant of the covariance matrix.
Definition: Gaussian.h:125
virtual Eigen::MatrixXd ApplyCovariance(Eigen::Ref< const Eigen::MatrixXd > const &x) const override
Definition: Gaussian.cpp:248
virtual Eigen::MatrixXd ApplyPrecision(Eigen::Ref< const Eigen::MatrixXd > const &x) const override
Definition: Gaussian.cpp:264
Eigen::MatrixXd GetCovariance() const
Get the covariance.
Definition: Gaussian.cpp:91
void ResetHyperparameters(ref_vector< Eigen::VectorXd > const &params) override
Definition: Gaussian.cpp:153
virtual ~Gaussian()=default
std::vector< std::reference_wrapper< const T > > ref_vector
A vector of references to something ...
Definition: WorkPiece.h:37