MUQ  0.4.3
GaussianBase.h
Go to the documentation of this file.
1 #ifndef GAUSSIANBASE_H_
2 #define GAUSSIANBASE_H_
3 
5 
6 namespace muq {
7  namespace Modeling {
8 
15  class GaussianBase : public Distribution {
16  public:
17 
21  GaussianBase(unsigned int dim);
22 
28  GaussianBase(unsigned int dim,
29  Eigen::VectorXi const& hyperSizesIn);
30 
31 
36  GaussianBase(Eigen::VectorXd const& meanIn);
37 
44  GaussianBase(Eigen::VectorXd const& meanIn,
45  Eigen::VectorXi const& hyperSizesIn);
46 
47 
48  virtual ~GaussianBase() = default;
49 
53  virtual unsigned int Dimension() const;
54 
60  virtual Eigen::MatrixXd ApplyCovariance(Eigen::Ref<const Eigen::MatrixXd> const& x) const = 0;
61 
67  virtual Eigen::MatrixXd ApplyPrecision(Eigen::Ref<const Eigen::MatrixXd> const& x) const = 0;
68 
74  virtual Eigen::MatrixXd ApplyCovSqrt(Eigen::Ref<const Eigen::MatrixXd> const& x) const = 0;
75 
85  virtual Eigen::MatrixXd ApplyPrecSqrt(Eigen::Ref<const Eigen::MatrixXd> const& x) const = 0;
86 
91  virtual Eigen::VectorXd const& GetMean() const{return mean;};
92 
97  virtual void SetMean(Eigen::VectorXd const& newMu);
98 
105  virtual double LogDeterminant() const{return 0.0;};
106 
111  virtual void ResetHyperparameters(ref_vector<Eigen::VectorXd> const& params){};
112 
120  virtual Eigen::VectorXd GradLogDensityImpl(unsigned int wrt,
121  ref_vector<Eigen::VectorXd> const& inputs) override;
122 
123  virtual Eigen::VectorXd ApplyLogDensityHessianImpl(unsigned int wrt1,
124  unsigned int wrt2,
125  ref_vector<Eigen::VectorXd> const& inputs,
126  Eigen::VectorXd const& vec) override;
127  protected:
128 
134  virtual double LogDensityImpl(ref_vector<Eigen::VectorXd> const& inputs) override;
135 
137  virtual Eigen::VectorXd SampleImpl(ref_vector<Eigen::VectorXd> const& inputs) override;
138 
139  // Space to store the mean of the distribution
140  Eigen::VectorXd mean;
141 
142  };
143  } // namespace Modeling
144 } // namespace muq
145 
146 #endif
Defines an abstract Gaussian class.@seealso Gaussian.
Definition: GaussianBase.h:15
virtual Eigen::VectorXd SampleImpl(ref_vector< Eigen::VectorXd > const &inputs) override
Sample the distribution.
virtual ~GaussianBase()=default
virtual Eigen::VectorXd const & GetMean() const
Definition: GaussianBase.h:91
virtual unsigned int Dimension() const
virtual double LogDeterminant() const
Definition: GaussianBase.h:105
virtual Eigen::MatrixXd ApplyPrecSqrt(Eigen::Ref< const Eigen::MatrixXd > const &x) const =0
virtual Eigen::VectorXd GradLogDensityImpl(unsigned int wrt, ref_vector< Eigen::VectorXd > const &inputs) override
virtual Eigen::MatrixXd ApplyCovariance(Eigen::Ref< const Eigen::MatrixXd > const &x) const =0
virtual Eigen::MatrixXd ApplyPrecision(Eigen::Ref< const Eigen::MatrixXd > const &x) const =0
GaussianBase(unsigned int dim)
Definition: GaussianBase.cpp:8
virtual void SetMean(Eigen::VectorXd const &newMu)
virtual double LogDensityImpl(ref_vector< Eigen::VectorXd > const &inputs) override
virtual Eigen::VectorXd ApplyLogDensityHessianImpl(unsigned int wrt1, unsigned int wrt2, ref_vector< Eigen::VectorXd > const &inputs, Eigen::VectorXd const &vec) override
virtual Eigen::MatrixXd ApplyCovSqrt(Eigen::Ref< const Eigen::MatrixXd > const &x) const =0
virtual void ResetHyperparameters(ref_vector< Eigen::VectorXd > const &params)
Definition: GaussianBase.h:111
std::vector< std::reference_wrapper< const T > > ref_vector
A vector of references to something ...
Definition: WorkPiece.h:37