#include <GaussianProcess.h>
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 GaussianProcess & | Condition (Eigen::Ref< const Eigen::MatrixXd > const &loc, Eigen::Ref< const Eigen::MatrixXd > const &vals) |
virtual GaussianProcess & | Condition (Eigen::Ref< const Eigen::MatrixXd > const &loc, Eigen::Ref< const Eigen::MatrixXd > const &vals, double obsVar) |
virtual GaussianProcess & | Condition (std::shared_ptr< ObservationInformation > obs) |
std::shared_ptr< muq::Modeling::Gaussian > | Discretize (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< MeanFunctionBase > | Mean () |
std::shared_ptr< KernelBase > | Kernel () |
Enumerator | |
---|---|
DiagonalCov | |
BlockCov | |
FullCov | |
NoCov |
Definition at line 285 of file GaussianProcess.h.
|
inline |
Definition at line 293 of file GaussianProcess.h.
GaussianProcess::GaussianProcess | ( | std::shared_ptr< MeanFunctionBase > | meanIn, |
std::shared_ptr< KernelBase > | covKernelIn | ||
) |
Definition at line 30 of file GaussianProcess.cpp.
|
virtualdefault |
|
protected |
Definition at line 191 of file GaussianProcess.cpp.
References coDim, covKernel, obsDim, and observations.
Referenced by Predict(), and PredictMean().
|
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().
|
virtual |
Definition at line 39 of file GaussianProcess.cpp.
References coDim, and Condition().
|
virtual |
Definition at line 62 of file GaussianProcess.cpp.
References hasNewObs, and observations.
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.
|
inline |
Definition at line 348 of file GaussianProcess.h.
References covKernel.
Referenced by muq::Approximation::nlopt_obj().
|
virtual |
Reimplemented in muq::Approximation::StateSpaceGP.
Definition at line 363 of file GaussianProcess.cpp.
References FullCov, mean, Predict(), and ProcessObservations().
|
virtual |
Definition at line 386 of file GaussianProcess.cpp.
Referenced by muq::Approximation::nlopt_obj().
|
inlinevirtual |
Definition at line 343 of file GaussianProcess.h.
References MarginalLogLikelihood().
Referenced by MarginalLogLikelihood().
|
virtual |
Reimplemented in muq::Approximation::StateSpaceGP.
Definition at line 393 of file GaussianProcess.cpp.
|
inline |
Definition at line 347 of file GaussianProcess.h.
References mean.
|
virtual |
CURRENTLY NOT IMPLEMENTED
Definition at line 160 of file GaussianProcess.cpp.
References covKernel, muq::Approximation::OptInfo::gp, and ProcessObservations().
|
virtual |
Evaluate the mean and covariance at the locations in newLocs
Reimplemented in muq::Approximation::StateSpaceGP.
Definition at line 210 of file GaussianProcess.cpp.
References BlockCov, BuildCrossCov(), coDim, covKernel, covSolver, DiagonalCov, FullCov, mean, nugget, observations, ProcessObservations(), sigmaTrainDiff, SolveFromEig(), and useEig.
Referenced by Discretize(), LogLikelihood(), muq::Approximation::StateSpaceGP::Predict(), and Sample().
|
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.
|
protected |
Definition at line 101 of file GaussianProcess.cpp.
References covKernel, covSolver, covSolverEig, hasNewObs, mean, obsDim, observations, sigmaTrainDiff, SolveFromEig(), trainDiff, and useEig.
Referenced by LogLikelihood(), Optimize(), Predict(), PredictMean(), and Sample().
|
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().
|
protected |
Definition at line 82 of file GaussianProcess.cpp.
References covSolverEig, and nlohmann::detail::dtoa_impl::e.
Referenced by Predict(), and ProcessObservations().
|
protected |
Definition at line 375 of file GaussianProcess.h.
Referenced by BuildCrossCov(), Condition(), Predict(), muq::Approximation::StateSpaceGP::Predict(), and PredictMean().
|
protected |
Definition at line 360 of file GaussianProcess.h.
Referenced by BuildCrossCov(), Kernel(), Optimize(), Predict(), and ProcessObservations().
|
protected |
Definition at line 368 of file GaussianProcess.h.
Referenced by Predict(), and ProcessObservations().
|
protected |
Definition at line 370 of file GaussianProcess.h.
Referenced by ProcessObservations(), and SolveFromEig().
|
protected |
Definition at line 378 of file GaussianProcess.h.
Referenced by Condition(), ProcessObservations(), and muq::Approximation::StateSpaceGP::SortObservations().
|
protected |
Definition at line 374 of file GaussianProcess.h.
|
protected |
Definition at line 359 of file GaussianProcess.h.
Referenced by Discretize(), LogLikelihood(), Mean(), Predict(), PredictMean(), ProcessObservations(), and Sample().
|
protected |
Definition at line 381 of file GaussianProcess.h.
Referenced by Predict().
|
protected |
Definition at line 373 of file GaussianProcess.h.
Referenced by BuildCrossCov(), and ProcessObservations().
|
protected |
Definition at line 362 of file GaussianProcess.h.
Referenced by BuildCrossCov(), Condition(), Predict(), muq::Approximation::StateSpaceGP::Predict(), PredictMean(), ProcessObservations(), muq::Approximation::StateSpaceGP::Sample(), and muq::Approximation::StateSpaceGP::SortObservations().
|
protected |
Definition at line 380 of file GaussianProcess.h.
|
protected |
Definition at line 366 of file GaussianProcess.h.
Referenced by Predict(), PredictMean(), and ProcessObservations().
|
protected |
Definition at line 365 of file GaussianProcess.h.
Referenced by ProcessObservations().
|
protected |
Definition at line 371 of file GaussianProcess.h.
Referenced by Predict(), and ProcessObservations().