MUQ  0.4.3
muq::Approximation::KernelBase Class Referenceabstract

Base class for all covariance kernels. More...

#include <KernelBase.h>

Inheritance diagram for muq::Approximation::KernelBase:

Detailed Description

Base class for all covariance kernels.

Definition at line 35 of file KernelBase.h.

Public Member Functions

 KernelBase (unsigned int inputDimIn, unsigned int coDimIn, unsigned int numParamsIn)
 
 KernelBase (unsigned int inputDimIn, std::vector< unsigned int > dimIndsIn, unsigned int coDimIn, unsigned int numParamsIn)
 
virtual ~KernelBase ()
 
virtual std::vector< std::shared_ptr< KernelBase > > GetSeperableComponents ()
 Overridden by ProductKernel. More...
 
virtual Eigen::MatrixXd Evaluate (Eigen::VectorXd const &x1, Eigen::VectorXd const &x2) const
 
virtual Eigen::MatrixXd BuildCovariance (Eigen::MatrixXd const &x) const
 
virtual Eigen::MatrixXd BuildCovariance (Eigen::MatrixXd const &x1, Eigen::MatrixXd const &x2) const
 
virtual void FillCovariance (Eigen::MatrixXd const &xs, Eigen::MatrixXd const &ys, Eigen::Ref< Eigen::MatrixXd > cov) const
 
virtual void FillCovariance (Eigen::MatrixXd const &xs, Eigen::Ref< Eigen::MatrixXd > cov) const
 
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 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. More...
 
virtual Eigen::MatrixXd GetParamBounds () const
 
virtual Eigen::VectorXd GetParams () const
 
virtual void SetParams (Eigen::VectorXd const &params)
 
virtual std::shared_ptr< KernelBaseClone () const =0
 
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. More...
 
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 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
 

Public Attributes

const std::vector< unsigned int > dimInds
 
const unsigned int inputDim
 
const unsigned int coDim
 
const unsigned int numParams
 

Constructor & Destructor Documentation

◆ KernelBase() [1/2]

muq::Approximation::KernelBase::KernelBase ( unsigned int  inputDimIn,
unsigned int  coDimIn,
unsigned int  numParamsIn 
)
inline

Definition at line 40 of file KernelBase.h.

◆ KernelBase() [2/2]

muq::Approximation::KernelBase::KernelBase ( unsigned int  inputDimIn,
std::vector< unsigned int >  dimIndsIn,
unsigned int  coDimIn,
unsigned int  numParamsIn 
)
inline

Definition at line 45 of file KernelBase.h.

References coDim, and inputDim.

◆ ~KernelBase()

virtual muq::Approximation::KernelBase::~KernelBase ( )
inlinevirtual

Definition at line 55 of file KernelBase.h.

Member Function Documentation

◆ BuildCovariance() [1/2]

Eigen::MatrixXd KernelBase::BuildCovariance ( Eigen::MatrixXd const &  x) const
virtual

Definition at line 24 of file KernelBase.cpp.

References coDim, and FillCovariance().

◆ BuildCovariance() [2/2]

Eigen::MatrixXd KernelBase::BuildCovariance ( Eigen::MatrixXd const &  x1,
Eigen::MatrixXd const &  x2 
) const
virtual

Definition at line 31 of file KernelBase.cpp.

References coDim, and FillCovariance().

◆ BuildDimInds()

static std::vector<unsigned> muq::Approximation::KernelBase::BuildDimInds ( unsigned  dim)
inlinestaticprivate

Definition at line 162 of file KernelBase.h.

◆ Clone()

◆ Evaluate()

Eigen::MatrixXd KernelBase::Evaluate ( Eigen::VectorXd const &  x1,
Eigen::VectorXd const &  x2 
) const
virtual

Definition at line 5 of file KernelBase.cpp.

References cachedParams, coDim, dimInds, and FillBlock().

◆ FillBlock()

virtual void muq::Approximation::KernelBase::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
pure virtual

◆ FillCovariance() [1/2]

void KernelBase::FillCovariance ( Eigen::MatrixXd const &  xs,
Eigen::MatrixXd const &  ys,
Eigen::Ref< Eigen::MatrixXd >  cov 
) const
virtual

Definition at line 64 of file KernelBase.cpp.

References cachedParams, coDim, dimInds, and FillBlock().

Referenced by BuildCovariance().

◆ FillCovariance() [2/2]

void KernelBase::FillCovariance ( Eigen::MatrixXd const &  xs,
Eigen::Ref< Eigen::MatrixXd >  cov 
) const
virtual

Definition at line 39 of file KernelBase.cpp.

References cachedParams, coDim, dimInds, and FillBlock().

◆ FillDerivCovariance()

void KernelBase::FillDerivCovariance ( Eigen::MatrixXd const &  xs,
Eigen::MatrixXd const &  ys,
std::vector< int > const &  wrts,
Eigen::Ref< Eigen::MatrixXd >  cov 
) const
virtual

Definition at line 87 of file KernelBase.cpp.

References cachedParams, coDim, dimInds, and FillPosDerivBlock().

Referenced by GetPosDerivative().

◆ FillPosDerivBlock()

virtual void muq::Approximation::KernelBase::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
pure virtual

◆ GetParamBounds()

virtual Eigen::MatrixXd muq::Approximation::KernelBase::GetParamBounds ( ) const
inlinevirtual

Definition at line 103 of file KernelBase.h.

References paramBounds.

◆ GetParams()

virtual Eigen::VectorXd muq::Approximation::KernelBase::GetParams ( ) const
inlinevirtual

Definition at line 109 of file KernelBase.h.

References cachedParams.

◆ GetPosDerivative()

Eigen::MatrixXd KernelBase::GetPosDerivative ( Eigen::VectorXd const &  x1,
Eigen::VectorXd const &  x2,
std::vector< int > const &  wrts 
) const
virtual

Returns derivatives of the kernel with respect to the first input, x1.

Parameters
[in]x1The first position passed to the kernel.
[in]x2The second position passed to the kernel.
[in]wrtsA vector defining the order and directions of the spatial derivatives. wrts.size() is the derivative order.
Returns
A matrix containing the derivatives of the kernel output with respect to the dimensions defined by wrts.

Definition at line 112 of file KernelBase.cpp.

References coDim, and FillDerivCovariance().

◆ GetSeperableComponents()

virtual std::vector<std::shared_ptr<KernelBase> > muq::Approximation::KernelBase::GetSeperableComponents ( )
inlinevirtual

Overridden by ProductKernel.

Reimplemented in muq::Approximation::ProductKernel.

Definition at line 63 of file KernelBase.h.

References Clone().

◆ GetStateSpace()

virtual std::tuple<std::shared_ptr<muq::Modeling::LinearSDE>, std::shared_ptr<muq::Modeling::LinearOperator>, Eigen::MatrixXd> muq::Approximation::KernelBase::GetStateSpace ( boost::property_tree::ptree  sdeOptions = boost::property_tree::ptree()) const
inlinevirtual

Returns a state space representation of the covariance kernel.

If this is a one dimensional kernel (i.e., inputDim=1 and coDim=1), this function returns a state space representation of the covariance kernel. In particular, it returns a linear time invariant stochastic differential equation, whose solution, when started with the returned stationary covariance, provides the same information as this Gaussian process. The first component of the vector-valued stochastic differential equation is related to the Gaussian process. See "Kalman filtering and smoothing solutions to temporal Gaussian process regression models," by Jouni Hartikainen and Simo Sarkka, for more information.

Reimplemented in muq::Approximation::SumKernel, muq::Approximation::ProductKernel, muq::Approximation::MaternKernel, muq::Approximation::ConstantKernel, and muq::Approximation::ConcatenateKernel.

Definition at line 124 of file KernelBase.h.

◆ SetParams()

virtual void muq::Approximation::KernelBase::SetParams ( Eigen::VectorXd const &  params)
inlinevirtual

Definition at line 111 of file KernelBase.h.

References cachedParams, and numParams.

Member Data Documentation

◆ cachedParams

◆ coDim

◆ dimInds

const std::vector<unsigned int> muq::Approximation::KernelBase::dimInds

◆ inputDim

const unsigned int muq::Approximation::KernelBase::inputDim

Definition at line 132 of file KernelBase.h.

Referenced by muq::Approximation::ConstantKernel::GetStateSpace(), and KernelBase().

◆ numParams

◆ paramBounds


The documentation for this class was generated from the following files: