MUQ  0.4.3
GaussianBase.cpp
Go to the documentation of this file.
2 
4 
5 using namespace muq::Utilities;
6 using namespace muq::Modeling;
7 
8 GaussianBase::GaussianBase(unsigned int dim) : Distribution(dim, Eigen::VectorXi()),
9  mean(Eigen::VectorXd::Zero(dim))
10 {}
11 
12 GaussianBase::GaussianBase(unsigned int dim,
13  Eigen::VectorXi const& hyperSizesIn) : Distribution(dim, hyperSizesIn),
14  mean(Eigen::VectorXd::Zero(dim))
15 {}
16 
17 
18 
19 GaussianBase::GaussianBase(Eigen::VectorXd const& muIn) : Distribution(muIn.size(), Eigen::VectorXi()),
20  mean(muIn)
21 {}
22 
23 
24 GaussianBase::GaussianBase(Eigen::VectorXd const& muIn,
25  Eigen::VectorXi const& hyperSizesIn) : Distribution(muIn.size(), hyperSizesIn),
26  mean(muIn)
27 {}
28 
29 
31 
32  ResetHyperparameters(ref_vector<Eigen::VectorXd>(inputs.begin()+1, inputs.end()));
33 
34  Eigen::VectorXd delta = inputs.at(0).get() - mean;
35 
36  return -0.5*varSize*std::log(2.0*M_PI) - 0.5*LogDeterminant() - 0.5 * delta.dot( ApplyPrecision(delta).col(0) );
37 }
38 
39 Eigen::VectorXd GaussianBase::SampleImpl(ref_vector<Eigen::VectorXd> const& inputs) {
40 
41  ResetHyperparameters(ref_vector<Eigen::VectorXd>(inputs.begin(), inputs.end()));
42 
43  Eigen::VectorXd z = RandomGenerator::GetNormal(mean.rows());
44 
45  return mean + ApplyCovSqrt(z);
46 }
47 
48 unsigned int GaussianBase::Dimension() const {
49  return mean.rows();
50 }
51 
52 void GaussianBase::SetMean(Eigen::VectorXd const& newMu)
53 {
54  assert(newMu.rows() == mean.rows());
55  mean = newMu;
56 }
57 
58 Eigen::VectorXd GaussianBase::GradLogDensityImpl(unsigned int wrt, ref_vector<Eigen::VectorXd> const& inputs)
59 {
60  Eigen::VectorXd delta = inputs.at(0).get() - mean;
61  if(wrt==0){
62  return -1.0 * ApplyPrecision(delta);
63  }else{
64  std::cerr << "ERROR: Gradient wrt mean and covariance has not been implemented." << std::endl;
65  assert(false);
66  }
67 }
68 
69 Eigen::VectorXd GaussianBase::ApplyLogDensityHessianImpl(unsigned int wrt1,
70  unsigned int wrt2,
71  ref_vector<Eigen::VectorXd> const& inputs,
72  Eigen::VectorXd const& vec)
73 {
74  if((wrt1==0)&&(wrt2==0)){
75  return -1.0 * ApplyPrecision(vec);
76  }else{
77  std::cerr << "ERROR: ApplyHessian wrt mean and covariance has not been implemented." << std::endl;
78  assert(false);
79  }
80 }
const unsigned int varSize
Definition: Distribution.h:145
virtual Eigen::VectorXd SampleImpl(ref_vector< Eigen::VectorXd > const &inputs) override
Sample the distribution.
virtual unsigned int Dimension() const
virtual double LogDeterminant() const
Definition: GaussianBase.h:105
virtual Eigen::VectorXd GradLogDensityImpl(unsigned int wrt, ref_vector< Eigen::VectorXd > const &inputs) override
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