18 #include <boost/property_tree/ptree.hpp>
28 namespace Approximation
35 class KernelBase :
public std::enable_shared_from_this<muq::Approximation::KernelBase>
46 std::vector<unsigned int> dimIndsIn,
66 virtual Eigen::MatrixXd
Evaluate(Eigen::VectorXd
const& x1,
67 Eigen::VectorXd
const& x2)
const;
69 virtual Eigen::MatrixXd
BuildCovariance(Eigen::MatrixXd
const& x)
const;
72 Eigen::MatrixXd
const& x2)
const;
75 Eigen::MatrixXd
const& ys,
76 Eigen::Ref<Eigen::MatrixXd> cov)
const;
79 Eigen::Ref<Eigen::MatrixXd> cov)
const;
83 Eigen::MatrixXd
const& ys,
84 std::vector<int>
const& wrts,
85 Eigen::Ref<Eigen::MatrixXd> cov)
const;
96 Eigen::VectorXd
const& x2,
97 std::vector<int>
const& wrts)
const;
113 virtual std::shared_ptr<KernelBase>
Clone()
const = 0;
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.");
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;
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;
164 std::vector<unsigned int> output(dim);
165 for(
unsigned int i=0; i<dim; ++i)
Base class for all covariance kernels.
const unsigned int inputDim
virtual Eigen::MatrixXd BuildCovariance(Eigen::MatrixXd const &x) const
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.
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.
virtual void SetParams(Eigen::VectorXd const ¶ms)
KernelBase(unsigned int inputDimIn, std::vector< unsigned int > dimIndsIn, unsigned int coDimIn, unsigned int numParamsIn)
const unsigned int numParams
virtual void FillPosDerivBlock(Eigen::Ref< const Eigen::VectorXd > const &x1, Eigen::Ref< const Eigen::VectorXd > const &x2, Eigen::Ref< const Eigen::VectorXd > const ¶ms, 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
KernelBase(unsigned int inputDimIn, unsigned int coDimIn, unsigned int numParamsIn)
Eigen::VectorXd cachedParams
virtual Eigen::VectorXd GetParams() const
const std::vector< unsigned int > dimInds
static std::vector< unsigned > BuildDimInds(unsigned dim)
virtual void FillDerivCovariance(Eigen::MatrixXd const &xs, Eigen::MatrixXd const &ys, std::vector< int > const &wrts, Eigen::Ref< Eigen::MatrixXd > cov) const
virtual Eigen::MatrixXd GetParamBounds() const
virtual std::vector< std::shared_ptr< KernelBase > > GetSeperableComponents()
Overridden by ProductKernel.
virtual void FillCovariance(Eigen::MatrixXd const &xs, Eigen::MatrixXd const &ys, Eigen::Ref< Eigen::MatrixXd > cov) const
Eigen::MatrixXd paramBounds
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 ¶ms, Eigen::Ref< Eigen::MatrixXd > block) const =0
Class for virtual base functions that are not implemented.