hops
Public Member Functions | List of all members
hops::GaussianProcess< MatrixType, VectorType, Kernel > Class Template Reference

#include <GaussianProcess.hpp>

Collaboration diagram for hops::GaussianProcess< MatrixType, VectorType, Kernel >:
Collaboration graph

Public Member Functions

 GaussianProcess (Kernel kernel, double constantPriorMean=0)
 
 GaussianProcess (Kernel kernel, std::function< double(VectorType)> priorMeanFunction)
 
GaussianProcess getPriorCopy ()
 
GaussianProcess getPosteriorCopy ()
 
void computePosterior (const MatrixType &input)
 
Eigen::VectorXd sample (const MatrixType &input, hops::RandomNumberGenerator &randomNumberGenerator, size_t &maxElement)
 
Eigen::VectorXd sample (hops::RandomNumberGenerator &randomNumberGenerator, size_t &maxElement)
 
Eigen::VectorXd sample (const MatrixType &x, hops::RandomNumberGenerator &randomNumberGenerator)
 
Eigen::VectorXd sample (hops::RandomNumberGenerator &randomNumberGenerator)
 
std::tuple< MatrixType, MatrixType, VectorTypeupdateObservations (const MatrixType &x, const VectorType &y, const VectorType &error, bool isUnique=false)
 
void updateObservations (const MatrixType &x, const VectorType &y)
 
void addObservations (const MatrixType &x, VectorType &y, VectorType &error)
 
void addObservations (const MatrixType &x, const VectorType &y)
 
const VectorTypegetPosteriorMean ()
 
const MatrixTypegetPosteriorCovariance ()
 
const MatrixTypegetSqrtPosteriorCovariance ()
 
const MatrixTypegetObservedInputs () const
 
const VectorTypegetObservedValues () const
 
const VectorTypegetObservedValueErrors () const
 
const MatrixTypegetObservedCovariance () const
 
std::function< double(VectorType)> & getPriorMeanFunction ()
 
void setKernelSigma (double m_sigma)
 
Kernel getKernel ()
 

Constructor & Destructor Documentation

◆ GaussianProcess() [1/2]

template<typename MatrixType , typename VectorType , typename Kernel >
hops::GaussianProcess< MatrixType, VectorType, Kernel >::GaussianProcess ( Kernel  kernel,
double  constantPriorMean = 0 
)
inline

◆ GaussianProcess() [2/2]

template<typename MatrixType , typename VectorType , typename Kernel >
hops::GaussianProcess< MatrixType, VectorType, Kernel >::GaussianProcess ( Kernel  kernel,
std::function< double(VectorType)>  priorMeanFunction 
)
inline

Member Function Documentation

◆ addObservations() [1/2]

template<typename MatrixType , typename VectorType , typename Kernel >
void hops::GaussianProcess< MatrixType, VectorType, Kernel >::addObservations ( const MatrixType x,
const VectorType y 
)
inline

◆ addObservations() [2/2]

template<typename MatrixType , typename VectorType , typename Kernel >
void hops::GaussianProcess< MatrixType, VectorType, Kernel >::addObservations ( const MatrixType x,
VectorType y,
VectorType error 
)
inline

◆ computePosterior()

template<typename MatrixType , typename VectorType , typename Kernel >
void hops::GaussianProcess< MatrixType, VectorType, Kernel >::computePosterior ( const MatrixType input)
inline

Compute the posterior mean and covariance for the given input. This method checks, if any data has changed before recomputing depending quantities, which makes it rather safe to call it before sampling or querying any posterior data

◆ getKernel()

template<typename MatrixType , typename VectorType , typename Kernel >
Kernel hops::GaussianProcess< MatrixType, VectorType, Kernel >::getKernel ( )
inline

◆ getObservedCovariance()

template<typename MatrixType , typename VectorType , typename Kernel >
const MatrixType& hops::GaussianProcess< MatrixType, VectorType, Kernel >::getObservedCovariance ( ) const
inline

◆ getObservedInputs()

template<typename MatrixType , typename VectorType , typename Kernel >
const MatrixType& hops::GaussianProcess< MatrixType, VectorType, Kernel >::getObservedInputs ( ) const
inline

◆ getObservedValueErrors()

template<typename MatrixType , typename VectorType , typename Kernel >
const VectorType& hops::GaussianProcess< MatrixType, VectorType, Kernel >::getObservedValueErrors ( ) const
inline

◆ getObservedValues()

template<typename MatrixType , typename VectorType , typename Kernel >
const VectorType& hops::GaussianProcess< MatrixType, VectorType, Kernel >::getObservedValues ( ) const
inline

◆ getPosteriorCopy()

template<typename MatrixType , typename VectorType , typename Kernel >
GaussianProcess hops::GaussianProcess< MatrixType, VectorType, Kernel >::getPosteriorCopy ( )
inline

◆ getPosteriorCovariance()

template<typename MatrixType , typename VectorType , typename Kernel >
const MatrixType& hops::GaussianProcess< MatrixType, VectorType, Kernel >::getPosteriorCovariance ( )
inline

◆ getPosteriorMean()

template<typename MatrixType , typename VectorType , typename Kernel >
const VectorType& hops::GaussianProcess< MatrixType, VectorType, Kernel >::getPosteriorMean ( )
inline

◆ getPriorCopy()

template<typename MatrixType , typename VectorType , typename Kernel >
GaussianProcess hops::GaussianProcess< MatrixType, VectorType, Kernel >::getPriorCopy ( )
inline

◆ getPriorMeanFunction()

template<typename MatrixType , typename VectorType , typename Kernel >
std::function<double (VectorType)>& hops::GaussianProcess< MatrixType, VectorType, Kernel >::getPriorMeanFunction ( )
inline

◆ getSqrtPosteriorCovariance()

template<typename MatrixType , typename VectorType , typename Kernel >
const MatrixType& hops::GaussianProcess< MatrixType, VectorType, Kernel >::getSqrtPosteriorCovariance ( )
inline

◆ sample() [1/4]

template<typename MatrixType , typename VectorType , typename Kernel >
Eigen::VectorXd hops::GaussianProcess< MatrixType, VectorType, Kernel >::sample ( const MatrixType input,
hops::RandomNumberGenerator randomNumberGenerator,
size_t &  maxElement 
)
inline

Sample from posterior at x

◆ sample() [2/4]

template<typename MatrixType , typename VectorType , typename Kernel >
Eigen::VectorXd hops::GaussianProcess< MatrixType, VectorType, Kernel >::sample ( const MatrixType x,
hops::RandomNumberGenerator randomNumberGenerator 
)
inline

◆ sample() [3/4]

template<typename MatrixType , typename VectorType , typename Kernel >
Eigen::VectorXd hops::GaussianProcess< MatrixType, VectorType, Kernel >::sample ( hops::RandomNumberGenerator randomNumberGenerator)
inline

◆ sample() [4/4]

template<typename MatrixType , typename VectorType , typename Kernel >
Eigen::VectorXd hops::GaussianProcess< MatrixType, VectorType, Kernel >::sample ( hops::RandomNumberGenerator randomNumberGenerator,
size_t &  maxElement 
)
inline

◆ setKernelSigma()

template<typename MatrixType , typename VectorType , typename Kernel >
void hops::GaussianProcess< MatrixType, VectorType, Kernel >::setKernelSigma ( double  m_sigma)
inline

◆ updateObservations() [1/2]

template<typename MatrixType , typename VectorType , typename Kernel >
void hops::GaussianProcess< MatrixType, VectorType, Kernel >::updateObservations ( const MatrixType x,
const VectorType y 
)
inline

◆ updateObservations() [2/2]

template<typename MatrixType , typename VectorType , typename Kernel >
std::tuple<MatrixType, MatrixType, VectorType> hops::GaussianProcess< MatrixType, VectorType, Kernel >::updateObservations ( const MatrixType x,
const VectorType y,
const VectorType error,
bool  isUnique = false 
)
inline

Given an observation (x_i, y_i, eps_i) stored in *this, if there is an x_j in the argument x passed, s.t. x_i = x_j, then y_i := x_j and eps_i := y_j

All x_j and the respective y_j and eps_j from the passed arguments, which were not found in *this, will be stored in the reference arguments x, y and error.

The isUnique argument controls, whether x may be assumed to be unique in *this. If isUnique == true, then only the data at the first occurence of x_i in *this will be updated.


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