MUQ  0.4.3
muq::Approximation::GaussianProcess Class Reference

#include <GaussianProcess.h>

Inheritance diagram for muq::Approximation::GaussianProcess:

Detailed Description

Definition at line 280 of file GaussianProcess.h.

Public Types

enum  CovarianceType { DiagonalCov , BlockCov , FullCov , NoCov }
 

Public Member Functions

 GaussianProcess (MeanFunctionBase &meanIn, KernelBase &kernelIn)
 
 GaussianProcess (std::shared_ptr< MeanFunctionBase > meanIn, std::shared_ptr< KernelBase > covKernelIn)
 
virtual ~GaussianProcess ()=default
 
virtual GaussianProcessCondition (Eigen::Ref< const Eigen::MatrixXd > const &loc, Eigen::Ref< const Eigen::MatrixXd > const &vals)
 
virtual GaussianProcessCondition (Eigen::Ref< const Eigen::MatrixXd > const &loc, Eigen::Ref< const Eigen::MatrixXd > const &vals, double obsVar)
 
virtual GaussianProcessCondition (std::shared_ptr< ObservationInformation > obs)
 
std::shared_ptr< muq::Modeling::GaussianDiscretize (Eigen::MatrixXd const &pts)
 
virtual void Optimize ()
 
virtual std::pair< Eigen::MatrixXd, Eigen::MatrixXd > Predict (Eigen::MatrixXd const &newLocs, CovarianceType covType)
 
virtual Eigen::MatrixXd PredictMean (Eigen::MatrixXd const &newPts)
 
virtual Eigen::MatrixXd Sample (Eigen::MatrixXd const &newPts)
 
virtual double LogLikelihood (Eigen::MatrixXd const &xs, Eigen::MatrixXd const &vals)
 
virtual double MarginalLogLikelihood ()
 
virtual double MarginalLogLikelihood (Eigen::Ref< Eigen::VectorXd > grad)
 
virtual double MarginalLogLikelihood (Eigen::Ref< Eigen::VectorXd > grad, bool computeGrad)
 
std::shared_ptr< MeanFunctionBaseMean ()
 
std::shared_ptr< KernelBaseKernel ()
 

Member Enumeration Documentation

◆ CovarianceType

Enumerator
DiagonalCov 
BlockCov 
FullCov 
NoCov 

Definition at line 285 of file GaussianProcess.h.

Constructor & Destructor Documentation

◆ GaussianProcess() [1/2]

muq::Approximation::GaussianProcess::GaussianProcess ( MeanFunctionBase meanIn,
KernelBase kernelIn 
)
inline

Definition at line 293 of file GaussianProcess.h.

◆ GaussianProcess() [2/2]

GaussianProcess::GaussianProcess ( std::shared_ptr< MeanFunctionBase meanIn,
std::shared_ptr< KernelBase covKernelIn 
)

Definition at line 30 of file GaussianProcess.cpp.

◆ ~GaussianProcess()

virtual muq::Approximation::GaussianProcess::~GaussianProcess ( )
virtualdefault

Member Function Documentation

◆ BuildCrossCov()

Eigen::MatrixXd GaussianProcess::BuildCrossCov ( Eigen::MatrixXd const &  newLocs)
protected

Definition at line 191 of file GaussianProcess.cpp.

References coDim, covKernel, obsDim, and observations.

Referenced by Predict(), and PredictMean().

◆ Condition() [1/3]

virtual GaussianProcess& muq::Approximation::GaussianProcess::Condition ( Eigen::Ref< const Eigen::MatrixXd > const &  loc,
Eigen::Ref< const Eigen::MatrixXd > const &  vals 
)
inlinevirtual

Update this Gaussian process with with direct observations of the field at the columns of loc.

Definition at line 302 of file GaussianProcess.h.

Referenced by Condition().

◆ Condition() [2/3]

GaussianProcess & GaussianProcess::Condition ( Eigen::Ref< const Eigen::MatrixXd > const &  loc,
Eigen::Ref< const Eigen::MatrixXd > const &  vals,
double  obsVar 
)
virtual

Definition at line 39 of file GaussianProcess.cpp.

References coDim, and Condition().

◆ Condition() [3/3]

GaussianProcess & GaussianProcess::Condition ( std::shared_ptr< ObservationInformation obs)
virtual

Definition at line 62 of file GaussianProcess.cpp.

References hasNewObs, and observations.

◆ Discretize()

std::shared_ptr< muq::Modeling::Gaussian > GaussianProcess::Discretize ( Eigen::MatrixXd const &  pts)

Construct a Gaussian distribution (finite dimensional) by evaluating the mean function and kernel of this Gaussian process at the provided locations.

Definition at line 70 of file GaussianProcess.cpp.

References FullCov, mean, and Predict().

◆ Kernel()

std::shared_ptr<KernelBase> muq::Approximation::GaussianProcess::Kernel ( )
inline

Definition at line 348 of file GaussianProcess.h.

References covKernel.

Referenced by muq::Approximation::nlopt_obj().

◆ LogLikelihood()

double GaussianProcess::LogLikelihood ( Eigen::MatrixXd const &  xs,
Eigen::MatrixXd const &  vals 
)
virtual

Reimplemented in muq::Approximation::StateSpaceGP.

Definition at line 363 of file GaussianProcess.cpp.

References FullCov, mean, Predict(), and ProcessObservations().

◆ MarginalLogLikelihood() [1/3]

double GaussianProcess::MarginalLogLikelihood ( )
virtual

Definition at line 386 of file GaussianProcess.cpp.

Referenced by muq::Approximation::nlopt_obj().

◆ MarginalLogLikelihood() [2/3]

virtual double muq::Approximation::GaussianProcess::MarginalLogLikelihood ( Eigen::Ref< Eigen::VectorXd >  grad)
inlinevirtual

Definition at line 343 of file GaussianProcess.h.

References MarginalLogLikelihood().

Referenced by MarginalLogLikelihood().

◆ MarginalLogLikelihood() [3/3]

double GaussianProcess::MarginalLogLikelihood ( Eigen::Ref< Eigen::VectorXd >  grad,
bool  computeGrad 
)
virtual

Reimplemented in muq::Approximation::StateSpaceGP.

Definition at line 393 of file GaussianProcess.cpp.

◆ Mean()

std::shared_ptr<MeanFunctionBase> muq::Approximation::GaussianProcess::Mean ( )
inline

Definition at line 347 of file GaussianProcess.h.

References mean.

◆ Optimize()

void GaussianProcess::Optimize ( )
virtual

CURRENTLY NOT IMPLEMENTED

Definition at line 160 of file GaussianProcess.cpp.

References covKernel, muq::Approximation::OptInfo::gp, and ProcessObservations().

◆ Predict()

std::pair< Eigen::MatrixXd, Eigen::MatrixXd > GaussianProcess::Predict ( Eigen::MatrixXd const &  newLocs,
CovarianceType  covType 
)
virtual

◆ PredictMean()

Eigen::MatrixXd GaussianProcess::PredictMean ( Eigen::MatrixXd const &  newPts)
virtual

Evaluate the GP mean at the locations in newPts

Reimplemented in muq::Approximation::StateSpaceGP.

Definition at line 304 of file GaussianProcess.cpp.

References BuildCrossCov(), coDim, mean, observations, ProcessObservations(), and sigmaTrainDiff.

◆ ProcessObservations()

void GaussianProcess::ProcessObservations ( )
protected

◆ Sample()

Eigen::MatrixXd GaussianProcess::Sample ( Eigen::MatrixXd const &  newPts)
virtual

Draw a random sample from the GP at the specified points.

Reimplemented in muq::Approximation::StateSpaceGP.

Definition at line 324 of file GaussianProcess.cpp.

References FullCov, mean, Predict(), and ProcessObservations().

◆ SolveFromEig()

Eigen::MatrixXd GaussianProcess::SolveFromEig ( Eigen::MatrixXd const &  b) const
protected

Definition at line 82 of file GaussianProcess.cpp.

References covSolverEig, and nlohmann::detail::dtoa_impl::e.

Referenced by Predict(), and ProcessObservations().

Member Data Documentation

◆ coDim

const int muq::Approximation::GaussianProcess::coDim
protected

◆ covKernel

std::shared_ptr<KernelBase> muq::Approximation::GaussianProcess::covKernel
protected

Definition at line 360 of file GaussianProcess.h.

Referenced by BuildCrossCov(), Kernel(), Optimize(), Predict(), and ProcessObservations().

◆ covSolver

Eigen::LDLT<Eigen::MatrixXd> muq::Approximation::GaussianProcess::covSolver
protected

Definition at line 368 of file GaussianProcess.h.

Referenced by Predict(), and ProcessObservations().

◆ covSolverEig

Eigen::SelfAdjointEigenSolver<Eigen::MatrixXd> muq::Approximation::GaussianProcess::covSolverEig
protected

Definition at line 370 of file GaussianProcess.h.

Referenced by ProcessObservations(), and SolveFromEig().

◆ hasNewObs

bool muq::Approximation::GaussianProcess::hasNewObs
protected

◆ inputDim

const int muq::Approximation::GaussianProcess::inputDim
protected

Definition at line 374 of file GaussianProcess.h.

◆ mean

std::shared_ptr<MeanFunctionBase> muq::Approximation::GaussianProcess::mean
protected

◆ nugget

const double muq::Approximation::GaussianProcess::nugget = 1e-14
protected

Definition at line 381 of file GaussianProcess.h.

Referenced by Predict().

◆ obsDim

int muq::Approximation::GaussianProcess::obsDim
protected

Definition at line 373 of file GaussianProcess.h.

Referenced by BuildCrossCov(), and ProcessObservations().

◆ observations

◆ pi

const double muq::Approximation::GaussianProcess::pi = 4.0 * atan(1.0)
protected

Definition at line 380 of file GaussianProcess.h.

◆ sigmaTrainDiff

Eigen::VectorXd muq::Approximation::GaussianProcess::sigmaTrainDiff
protected

Definition at line 366 of file GaussianProcess.h.

Referenced by Predict(), PredictMean(), and ProcessObservations().

◆ trainDiff

Eigen::VectorXd muq::Approximation::GaussianProcess::trainDiff
protected

Definition at line 365 of file GaussianProcess.h.

Referenced by ProcessObservations().

◆ useEig

bool muq::Approximation::GaussianProcess::useEig = false
protected

Definition at line 371 of file GaussianProcess.h.

Referenced by Predict(), and ProcessObservations().


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