MUQ  0.4.3
StateSpaceGP.h
Go to the documentation of this file.
1 #ifndef STATESPACEGP_H
2 #define STATESPACEGP_H
3 
4 #include <memory>
5 
8 
10 
12 
13 #include <Eigen/Core>
14 
15 namespace muq
16 {
17 namespace Approximation
18 {
19 
20 
22 {
23 public:
24 
26  KernelBase& kernelIn,
27  boost::property_tree::ptree options = boost::property_tree::ptree()) : StateSpaceGP(meanIn.Clone(), kernelIn.Clone(), options){};
28 
29  StateSpaceGP(std::shared_ptr<MeanFunctionBase> meanIn,
30  std::shared_ptr<KernelBase> covKernelIn,
31  boost::property_tree::ptree options = boost::property_tree::ptree());
32 
33  virtual ~StateSpaceGP() = default;
34 
35  virtual Eigen::MatrixXd Sample(Eigen::MatrixXd const& times) override;
36 
37 
38  virtual std::pair<Eigen::MatrixXd, Eigen::MatrixXd> Predict(Eigen::MatrixXd const& newLocs,
39  CovarianceType covType) override;
40 
41  virtual Eigen::MatrixXd PredictMean(Eigen::MatrixXd const& newPts) override;
42 
43 
44  virtual double LogLikelihood(Eigen::MatrixXd const& xs,
45  Eigen::MatrixXd const& vals) override;
46 
47  virtual double MarginalLogLikelihood(Eigen::Ref<Eigen::VectorXd> grad,
48  bool computeGrad = true) override;
49 
50 
51  std::shared_ptr<muq::Modeling::LinearSDE> GetSDE(){return sde;};
52 
53  std::shared_ptr<muq::Modeling::LinearOperator> GetObs(){return obsOp;};
54 
55  void SetObs(std::shared_ptr<muq::Modeling::LinearOperator> newObs);
56 
57  Eigen::MatrixXd GetCov(){return L.triangularView<Eigen::Lower>()*L.transpose();};
58 
59 
60 
61  const int stateDim;
62 
63  //static std::shared_ptr<StateSpaceGP> Concatenate(std::vector<std::shared_ptr<StateSpaceGP>> const& gps);
64 
65 private:
66 
67  StateSpaceGP(std::tuple<std::shared_ptr<muq::Modeling::LinearSDE>, std::shared_ptr<muq::Modeling::LinearOperator>, Eigen::MatrixXd> ssInfo,
68  std::shared_ptr<MeanFunctionBase> meanIn,
69  std::shared_ptr<KernelBase> covKernelIn);
70 
71 
72  void SortObservations();
73 
74  bool ComputeAQ(double dt);
75  Eigen::MatrixXd sdeA, sdeQ;
76  double dtAQ; // the last deltat passed to the ComputeAQ function
77 
78  // Stocastic Differential equation describing correlations
79  std::shared_ptr<muq::Modeling::LinearSDE> sde;
80 
81  // Observation operator
82  std::shared_ptr<muq::Modeling::LinearOperator> obsOp;
83 
84  // Cholesky factor of the stationary covariance matrix (used to initialize SDE integration)
85  Eigen::MatrixXd L;
86 
87  std::shared_ptr<MeanFunctionBase> mean;
88  std::shared_ptr<KernelBase> covKernel;
89 
90 
91 };
92 
93 
94 } // namespace Approximation
95 } // namespace GP
96 
97 
98 
99 
100 
101 #endif
Base class for all covariance kernels.
Definition: KernelBase.h:36
virtual Eigen::MatrixXd PredictMean(Eigen::MatrixXd const &newPts) override
virtual std::pair< Eigen::MatrixXd, Eigen::MatrixXd > Predict(Eigen::MatrixXd const &newLocs, CovarianceType covType) override
std::shared_ptr< muq::Modeling::LinearOperator > GetObs()
Definition: StateSpaceGP.h:53
virtual Eigen::MatrixXd Sample(Eigen::MatrixXd const &times) override
void SetObs(std::shared_ptr< muq::Modeling::LinearOperator > newObs)
std::shared_ptr< MeanFunctionBase > mean
Definition: StateSpaceGP.h:87
std::shared_ptr< muq::Modeling::LinearSDE > sde
Definition: StateSpaceGP.h:79
std::shared_ptr< muq::Modeling::LinearOperator > obsOp
Definition: StateSpaceGP.h:82
std::shared_ptr< muq::Modeling::LinearSDE > GetSDE()
Definition: StateSpaceGP.h:51
virtual double LogLikelihood(Eigen::MatrixXd const &xs, Eigen::MatrixXd const &vals) override
std::shared_ptr< KernelBase > covKernel
Definition: StateSpaceGP.h:88
StateSpaceGP(MeanFunctionBase &meanIn, KernelBase &kernelIn, boost::property_tree::ptree options=boost::property_tree::ptree())
Definition: StateSpaceGP.h:25