MUQ  0.4.3
muq::Approximation::MaternKernel Class Reference

#include <MaternKernel.h>

Inheritance diagram for muq::Approximation::MaternKernel:

Detailed Description

This class implements a kernel of the form

\[ k(x,y) = \sigma^2 \frac{2^{1-\nu}}{\Gamma(\nu)}\left(\frac{\sqrt{2\nu}}{l}\tau\right)^\nu K_\nu\left(\frac{\sqrt{2\nu}}{l}\tau\right), \]

where \(\Gamma(\cdot)\) is the gamma function, \(K_\nu(\cdot)\) is a modified Bessel function of the second kind, \(\nu\) is a smoothness parameter, \(l\) is a lengthscale, and \(\sigma^2\) is the variance. Note that we only allow values of \(\nu\) that take the form \(i-0.5\) for some positive integer \(i\). Typical choices are \(\nu=1/2\), \(\nu=3/2\), and \(\nu=5/2\). In the limit, \(\nu\rightarrow\infty\), this Matern kernel converges to the squared exponential kernel (muq::Approximation::SquaredExpKernel) and for \(\nu=1/2\), the Matern kernel is equivalent to the exponential kernel associated with the Ornstein-Uhlenbeck process.

Note: the smoothness parameter \(\nu\) is not optimized as a hyperparameter.

Definition at line 33 of file MaternKernel.h.

Public Member Functions

 MaternKernel (unsigned dimIn, std::vector< unsigned > dimInds, double sigma2In, double lengthIn, double nuIn)
 
 MaternKernel (unsigned dimIn, std::vector< unsigned > dimInds, double sigma2In, double lengthIn, double nuIn, Eigen::Vector2d sigmaBounds, Eigen::Vector2d lengthBounds)
 
 MaternKernel (unsigned dimIn, double sigma2In, double lengthIn, double nuIn)
 
 MaternKernel (unsigned dimIn, double sigma2In, double lengthIn, double nuIn, Eigen::Vector2d sigmaBounds, Eigen::Vector2d lengthBounds)
 
virtual ~MaternKernel ()
 
template<typename ScalarType1 , typename ScalarType2 , typename ScalarType3 >
void FillBlockImpl (Eigen::Ref< const Eigen::Matrix< ScalarType1, Eigen::Dynamic, 1 >> const &x1, Eigen::Ref< const Eigen::Matrix< ScalarType1, Eigen::Dynamic, 1 >> const &x2, Eigen::Ref< const Eigen::Matrix< ScalarType2, Eigen::Dynamic, 1 >> const &params, Eigen::Ref< Eigen::Matrix< ScalarType3, Eigen::Dynamic, Eigen::Dynamic >> block) 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 override
 Returns a state space representation of the covariance kernel. More...
 
- Public Member Functions inherited from muq::Approximation::KernelImpl< MaternKernel >
 KernelImpl (unsigned inputDimIn, unsigned coDimIn, unsigned numParamsIn)
 
 KernelImpl (unsigned inputDimIn, std::vector< unsigned > dimIndsIn, unsigned coDimIn, unsigned numParamsIn)
 
virtual ~KernelImpl ()
 
virtual std::shared_ptr< KernelBaseClone () const override
 
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 override
 
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 override
 
void FillPosDerivBlockImpl (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
 
Eigen::Matrix< ScalarType, Eigen::Dynamic, 1 > GetSegment (Eigen::Ref< const Eigen::Matrix< ScalarType, Eigen::Dynamic, 1 >> const &input) const
 
- Public Member Functions inherited from muq::Approximation::KernelBase
 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)
 

Additional Inherited Members

- Public Attributes inherited from muq::Approximation::KernelBase
const std::vector< unsigned int > dimInds
 
const unsigned int inputDim
 
const unsigned int coDim
 
const unsigned int numParams
 

Constructor & Destructor Documentation

◆ MaternKernel() [1/4]

muq::Approximation::MaternKernel::MaternKernel ( unsigned  dimIn,
std::vector< unsigned >  dimInds,
double  sigma2In,
double  lengthIn,
double  nuIn 
)
inline

Definition at line 38 of file MaternKernel.h.

◆ MaternKernel() [2/4]

MaternKernel::MaternKernel ( unsigned  dimIn,
std::vector< unsigned >  dimInds,
double  sigma2In,
double  lengthIn,
double  nuIn,
Eigen::Vector2d  sigmaBounds,
Eigen::Vector2d  lengthBounds 
)

◆ MaternKernel() [3/4]

muq::Approximation::MaternKernel::MaternKernel ( unsigned  dimIn,
double  sigma2In,
double  lengthIn,
double  nuIn 
)
inline

Definition at line 60 of file MaternKernel.h.

◆ MaternKernel() [4/4]

MaternKernel::MaternKernel ( unsigned  dimIn,
double  sigma2In,
double  lengthIn,
double  nuIn,
Eigen::Vector2d  sigmaBounds,
Eigen::Vector2d  lengthBounds 
)

◆ ~MaternKernel()

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

Definition at line 79 of file MaternKernel.h.

Member Function Documentation

◆ BuildWeights()

Eigen::VectorXd MaternKernel::BuildWeights ( int  p)
staticprivate

Definition at line 136 of file MaternKernel.cpp.

References Factorial().

◆ CheckNu()

void MaternKernel::CheckNu ( ) const
private

Definition at line 73 of file MaternKernel.cpp.

References nu.

Referenced by MaternKernel().

◆ Factorial()

static int muq::Approximation::MaternKernel::Factorial ( int  n)
inlinestaticprivate

Definition at line 112 of file MaternKernel.h.

References nlohmann::detail::dtoa_impl::n.

Referenced by BuildWeights().

◆ FillBlockImpl()

template<typename ScalarType1 , typename ScalarType2 , typename ScalarType3 >
void muq::Approximation::MaternKernel::FillBlockImpl ( Eigen::Ref< const Eigen::Matrix< ScalarType1, Eigen::Dynamic, 1 >> const &  x1,
Eigen::Ref< const Eigen::Matrix< ScalarType1, Eigen::Dynamic, 1 >> const &  x2,
Eigen::Ref< const Eigen::Matrix< ScalarType2, Eigen::Dynamic, 1 >> const &  params,
Eigen::Ref< Eigen::Matrix< ScalarType3, Eigen::Dynamic, Eigen::Dynamic >>  block 
) const
inline

Definition at line 83 of file MaternKernel.h.

References nu, scale, and weights.

◆ GetStateSpace()

std::tuple< std::shared_ptr< muq::Modeling::LinearSDE >, std::shared_ptr< muq::Modeling::LinearOperator >, Eigen::MatrixXd > MaternKernel::GetStateSpace ( boost::property_tree::ptree  sdeOptions = boost::property_tree::ptree()) const
overridevirtual

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 from muq::Approximation::KernelBase.

Definition at line 84 of file MaternKernel.cpp.

References muq::Approximation::KernelBase::cachedParams, muq::Modeling::LyapunovSolver< ScalarType, FixedRows, FixedCols >::compute(), muq::Modeling::LyapunovSolver< ScalarType, FixedRows, FixedCols >::matrixX(), and nu.

Member Data Documentation

◆ nu

const double muq::Approximation::MaternKernel::nu
private

Definition at line 103 of file MaternKernel.h.

Referenced by CheckNu(), FillBlockImpl(), and GetStateSpace().

◆ scale

const double muq::Approximation::MaternKernel::scale
private

Definition at line 104 of file MaternKernel.h.

Referenced by FillBlockImpl().

◆ weights

const Eigen::VectorXd muq::Approximation::MaternKernel::weights
private

Definition at line 106 of file MaternKernel.h.

Referenced by FillBlockImpl().


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