MUQ  0.4.3
KernelBase.h
Go to the documentation of this file.
1 #ifndef KERNELBASE_H
2 #define KERNELBASE_H
3 
4 
5 #include <assert.h>
6 #include <memory>
7 #include <iostream>
8 #include <fstream>
9 
11 
13 
14 //#include "MUQ/Approximation/GaussianProcesses/StateSpaceGP.h"
15 
17 
18 #include <boost/property_tree/ptree.hpp>
19 
20 
21 namespace muq
22 {
23 namespace Modeling{
24  class LinearSDE;
25  class LinearOperator;
26 }
27 
28 namespace Approximation
29 {
30 
35 class KernelBase : public std::enable_shared_from_this<muq::Approximation::KernelBase>
36 {
37 
38 public:
39 
40  KernelBase(unsigned int inputDimIn,
41  unsigned int coDimIn,
42  unsigned int numParamsIn) : KernelBase(inputDimIn, BuildDimInds(inputDimIn), coDimIn, numParamsIn)
43  {};
44 
45  KernelBase(unsigned int inputDimIn,
46  std::vector<unsigned int> dimIndsIn,
47  unsigned int coDimIn,
48  unsigned int numParamsIn) : dimInds(dimIndsIn), inputDim(inputDimIn), coDim(coDimIn), numParams(numParamsIn)
49  {
50  assert(inputDim>0);
51  assert(coDim>0);
52  };
53 
54 
55  virtual ~KernelBase(){};
56 
57 
58  //virtual std::shared_ptr<muq::Approximation::KernelBase> GetPtr() {
59  // return shared_from_this();
60  //}
61 
63  virtual std::vector<std::shared_ptr<KernelBase>> GetSeperableComponents() {return std::vector<std::shared_ptr<KernelBase>>(1,Clone()); };
64 
65 
66  virtual Eigen::MatrixXd Evaluate(Eigen::VectorXd const& x1,
67  Eigen::VectorXd const& x2) const;
68 
69  virtual Eigen::MatrixXd BuildCovariance(Eigen::MatrixXd const& x) const;
70 
71  virtual Eigen::MatrixXd BuildCovariance(Eigen::MatrixXd const& x1,
72  Eigen::MatrixXd const& x2) const;
73 
74  virtual void FillCovariance(Eigen::MatrixXd const& xs,
75  Eigen::MatrixXd const& ys,
76  Eigen::Ref<Eigen::MatrixXd> cov) const;
77 
78  virtual void FillCovariance(Eigen::MatrixXd const& xs,
79  Eigen::Ref<Eigen::MatrixXd> cov) const;
80 
81 
82  virtual void FillDerivCovariance(Eigen::MatrixXd const& xs,
83  Eigen::MatrixXd const& ys,
84  std::vector<int> const& wrts,
85  Eigen::Ref<Eigen::MatrixXd> cov) const;
86 
95  virtual Eigen::MatrixXd GetPosDerivative(Eigen::VectorXd const& x1,
96  Eigen::VectorXd const& x2,
97  std::vector<int> const& wrts) const;
98 
99  // virtual Eigen::MatrixXd GetParamDerivative(Eigen::VectorXd const& x1,
100  // Eigen::VectorXd const& x2,
101  // std::vector<int> const& dimWrts) const = 0;
102 
103  virtual Eigen::MatrixXd GetParamBounds() const
104  {
105  return paramBounds;
106  };
107 
108 
109  virtual Eigen::VectorXd GetParams() const{return cachedParams;};
110 
111  virtual void SetParams(Eigen::VectorXd const& params){assert(params.size()==numParams); cachedParams = params;};
112 
113  virtual std::shared_ptr<KernelBase> Clone() const = 0;
114 
115 
124  virtual std::tuple<std::shared_ptr<muq::Modeling::LinearSDE>, std::shared_ptr<muq::Modeling::LinearOperator>, Eigen::MatrixXd> GetStateSpace(boost::property_tree::ptree sdeOptions=boost::property_tree::ptree()) const{
125  throw muq::NotImplementedError("ERROR. The GetStateSpace() function has not been implemented in this child of muq::Approximation::KernelBase.");
126  };
127 
128 
129  const std::vector<unsigned int> dimInds;
130 
131 
132  const unsigned int inputDim;
133  const unsigned int coDim;
134  const unsigned int numParams;
135 
139  virtual void FillPosDerivBlock(Eigen::Ref<const Eigen::VectorXd> const& x1,
140  Eigen::Ref<const Eigen::VectorXd> const& x2,
141  Eigen::Ref<const Eigen::VectorXd> const& params,
142  std::vector<int> const& wrts,
143  Eigen::Ref<Eigen::MatrixXd> block) const = 0;
144 
145 
149  virtual void FillBlock(Eigen::Ref<const Eigen::VectorXd> const& x1,
150  Eigen::Ref<const Eigen::VectorXd> const& x2,
151  Eigen::Ref<const Eigen::VectorXd> const& params,
152  Eigen::Ref<Eigen::MatrixXd> block) const = 0;
153 
154 protected:
155 
156  Eigen::VectorXd cachedParams;
157 
158  Eigen::MatrixXd paramBounds;
159 
160 private:
161 
162  static std::vector<unsigned> BuildDimInds(unsigned dim)
163  {
164  std::vector<unsigned int> output(dim);
165  for(unsigned int i=0; i<dim; ++i)
166  output[i] = i;
167  return output;
168  }
169 };
170 
171 }
172 }
173 
174 
175 
176 #endif
Base class for all covariance kernels.
Definition: KernelBase.h:36
const unsigned int inputDim
Definition: KernelBase.h:132
virtual Eigen::MatrixXd BuildCovariance(Eigen::MatrixXd const &x) const
Definition: KernelBase.cpp:24
virtual std::tuple< std::shared_ptr< muq::Modeling::LinearSDE >, std::shared_ptr< muq::Modeling::LinearOperator >, Eigen::MatrixXd > GetStateSpace(boost::property_tree::ptree sdeOptions=boost::property_tree::ptree()) const
Returns a state space representation of the covariance kernel.
Definition: KernelBase.h:124
virtual Eigen::MatrixXd GetPosDerivative(Eigen::VectorXd const &x1, Eigen::VectorXd const &x2, std::vector< int > const &wrts) const
Returns derivatives of the kernel with respect to the first input, x1.
Definition: KernelBase.cpp:112
const unsigned int coDim
Definition: KernelBase.h:133
virtual void SetParams(Eigen::VectorXd const &params)
Definition: KernelBase.h:111
KernelBase(unsigned int inputDimIn, std::vector< unsigned int > dimIndsIn, unsigned int coDimIn, unsigned int numParamsIn)
Definition: KernelBase.h:45
const unsigned int numParams
Definition: KernelBase.h:134
virtual void FillPosDerivBlock(Eigen::Ref< const Eigen::VectorXd > const &x1, Eigen::Ref< const Eigen::VectorXd > const &x2, Eigen::Ref< const Eigen::VectorXd > const &params, std::vector< int > const &wrts, Eigen::Ref< Eigen::MatrixXd > block) const =0
virtual Eigen::MatrixXd Evaluate(Eigen::VectorXd const &x1, Eigen::VectorXd const &x2) const
Definition: KernelBase.cpp:5
KernelBase(unsigned int inputDimIn, unsigned int coDimIn, unsigned int numParamsIn)
Definition: KernelBase.h:40
Eigen::VectorXd cachedParams
Definition: KernelBase.h:156
virtual Eigen::VectorXd GetParams() const
Definition: KernelBase.h:109
const std::vector< unsigned int > dimInds
Definition: KernelBase.h:126
static std::vector< unsigned > BuildDimInds(unsigned dim)
Definition: KernelBase.h:162
virtual void FillDerivCovariance(Eigen::MatrixXd const &xs, Eigen::MatrixXd const &ys, std::vector< int > const &wrts, Eigen::Ref< Eigen::MatrixXd > cov) const
Definition: KernelBase.cpp:87
virtual Eigen::MatrixXd GetParamBounds() const
Definition: KernelBase.h:103
virtual std::vector< std::shared_ptr< KernelBase > > GetSeperableComponents()
Overridden by ProductKernel.
Definition: KernelBase.h:63
virtual void FillCovariance(Eigen::MatrixXd const &xs, Eigen::MatrixXd const &ys, Eigen::Ref< Eigen::MatrixXd > cov) const
Definition: KernelBase.cpp:64
virtual std::shared_ptr< KernelBase > Clone() const =0
virtual void FillBlock(Eigen::Ref< const Eigen::VectorXd > const &x1, Eigen::Ref< const Eigen::VectorXd > const &x2, Eigen::Ref< const Eigen::VectorXd > const &params, Eigen::Ref< Eigen::MatrixXd > block) const =0
Class for virtual base functions that are not implemented.
Definition: Exceptions.h:22