MUQ  0.4.3
ConstantKernel.cpp
Go to the documentation of this file.
5 
6 using namespace muq::Approximation;
7 
9  const double sigma2In,
10  const Eigen::Vector2d sigmaBounds) : ConstantKernel(dim, sigma2In*Eigen::MatrixXd::Ones(1,1), sigmaBounds){};
11 
13  std::vector<unsigned> dimInds,
14  const double sigma2In,
15  const Eigen::Vector2d sigmaBounds) : ConstantKernel(dim, dimInds, sigma2In*Eigen::MatrixXd::Ones(1,1), sigmaBounds){};
16 
17 
19  Eigen::MatrixXd const& sigma2In,
20  const Eigen::Vector2d sigmaBounds) : KernelImpl<ConstantKernel>(dim, sigma2In.rows(), GetNumParams(sigma2In))
21 {
22  paramBounds.resize(2,1);
23  paramBounds(0,0) = sigmaBounds(0);
24  paramBounds(1,0) = sigmaBounds(1);
25 
26  cachedParams.resize(numParams);
27  int ind = 0;
28  for(int i=0; i<sigma2In.rows(); ++i){
29  for(int j=0; j<=i; ++j){
30  cachedParams(ind) = sigma2In(i,j);
31  ind++;
32  }
33  }
34 };
35 
37  std::vector<unsigned> dimInds,
38  Eigen::MatrixXd const& sigma2In,
39  const Eigen::Vector2d sigmaBounds) : KernelImpl<ConstantKernel>(dim, dimInds, sigma2In.rows(), GetNumParams(sigma2In))
40 {
41  paramBounds.resize(2,1);
42  paramBounds(0,0) = sigmaBounds(0);
43  paramBounds(1,0) = sigmaBounds(1);
44 
45  cachedParams.resize(numParams);
46  int ind = 0;
47  for(int i=0; i<sigma2In.rows(); ++i){
48  for(int j=0; j<=i; ++j){
49  cachedParams(ind) = sigma2In(i,j);
50  ind++;
51  }
52  }
53 };
54 
55 
56 std::tuple<std::shared_ptr<muq::Modeling::LinearSDE>, std::shared_ptr<muq::Modeling::LinearOperator>, Eigen::MatrixXd> ConstantKernel::GetStateSpace(boost::property_tree::ptree sdeOptions) const
57 {
58  std::shared_ptr<muq::Modeling::LinearOperator> id = std::make_shared<muq::Modeling::IdentityOperator>(coDim);
59 
60  Eigen::MatrixXd marginalVar(coDim, coDim);
61  Eigen::VectorXd x(inputDim);
62  FillBlockImpl<double,double,double>(x,x, cachedParams, marginalVar);
63 
64  boost::property_tree::ptree options;
65  options.put("SDE.dt", 1e6); // large step size to ensure that we only ever take one step
66 
67  std::shared_ptr<muq::Modeling::LinearOperator> zo = std::make_shared<muq::Modeling::ZeroOperator>(coDim,coDim);
68  auto sde = std::make_shared<muq::Modeling::LinearSDE>(zo, options);
69 
70  return std::make_tuple(sde, id, marginalVar);
71 }
72 
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 override
Returns a state space representation of the covariance kernel.
ConstantKernel(unsigned dim, const double sigma2In, const Eigen::Vector2d sigmaBounds={0.0, std::numeric_limits< double >::infinity()})
const unsigned int inputDim
Definition: KernelBase.h:132
const unsigned int coDim
Definition: KernelBase.h:133
const unsigned int numParams
Definition: KernelBase.h:134
Eigen::VectorXd cachedParams
Definition: KernelBase.h:156
Base class in CRTP pattern for covariance kernels.
Definition: KernelImpl.h:26