MUQ  0.4.3
muq::Modeling::Gaussian Class Reference

#include <Gaussian.h>

Inheritance diagram for muq::Modeling::Gaussian:

Detailed Description

@seealso GaussianBase

Definition at line 14 of file Gaussian.h.

Public Types

enum  Mode { Covariance , Precision }
 Are we specifying the mean, covariance matrix, or precision matrix. More...
 
enum  ExtraInputs {
  None = 1 << 0 , Mean = 1 << 1 , DiagCovariance = 1 << 2 , DiagPrecision = 1 << 3 ,
  FullCovariance = 1 << 4 , FullPrecision = 1 << 5
}
 
typedef uint8_t InputMask
 

Public Member Functions

 Gaussian (unsigned int dim, InputMask extraInputs=ExtraInputs::None)
 
 Gaussian (Eigen::VectorXd const &mu, InputMask extraInputs=ExtraInputs::None)
 Construct a Gaussian with scaled identity covariance/precision. More...
 
 Gaussian (Eigen::VectorXd const &mu, Eigen::MatrixXd const &obj, Gaussian::Mode mode=Gaussian::Mode::Covariance, InputMask extraInputs=ExtraInputs::None)
 Construct a Gaussian by specifying both the mean and covariance or precision matrix. More...
 
virtual ~Gaussian ()=default
 
Mode GetMode () const
 
virtual Eigen::MatrixXd ApplyCovariance (Eigen::Ref< const Eigen::MatrixXd > const &x) const override
 
virtual Eigen::MatrixXd ApplyPrecision (Eigen::Ref< const Eigen::MatrixXd > const &x) const override
 
virtual Eigen::MatrixXd ApplyCovSqrt (Eigen::Ref< const Eigen::MatrixXd > const &x) const override
 
virtual Eigen::MatrixXd ApplyPrecSqrt (Eigen::Ref< const Eigen::MatrixXd > const &x) const override
 
Eigen::MatrixXd GetCovariance () const
 Get the covariance. More...
 
Eigen::MatrixXd GetPrecision () const
 
void SetCovariance (Eigen::MatrixXd const &newCov)
 Set the covariance matrix. More...
 
void SetPrecision (Eigen::MatrixXd const &newPrec)
 Set the precision matrix. More...
 
Gaussian::InputMask GetInputTypes () const
 
std::shared_ptr< GaussianCondition (Eigen::MatrixXd const &obsMat, Eigen::VectorXd const &data, Eigen::MatrixXd const &obsCov) const
 Returns a new Gaussian distribution conditioned on a linear observation. More...
 
void ResetHyperparameters (ref_vector< Eigen::VectorXd > const &params) override
 
virtual double LogDeterminant () const override
 
- Public Member Functions inherited from muq::Modeling::GaussianBase
 GaussianBase (unsigned int dim)
 
 GaussianBase (unsigned int dim, Eigen::VectorXi const &hyperSizesIn)
 
 GaussianBase (Eigen::VectorXd const &meanIn)
 
 GaussianBase (Eigen::VectorXd const &meanIn, Eigen::VectorXi const &hyperSizesIn)
 
virtual ~GaussianBase ()=default
 
virtual unsigned int Dimension () const
 
virtual Eigen::VectorXd const & GetMean () const
 
virtual void SetMean (Eigen::VectorXd const &newMu)
 
virtual Eigen::VectorXd GradLogDensityImpl (unsigned int wrt, ref_vector< Eigen::VectorXd > const &inputs) override
 
virtual Eigen::VectorXd ApplyLogDensityHessianImpl (unsigned int wrt1, unsigned int wrt2, ref_vector< Eigen::VectorXd > const &inputs, Eigen::VectorXd const &vec) override
 
- Public Member Functions inherited from muq::Modeling::Distribution
 Distribution (unsigned int varSizeIn, Eigen::VectorXi const &hyperSizesIn=Eigen::VectorXi())
 
virtual ~Distribution ()=default
 
virtual double LogDensity (ref_vector< Eigen::VectorXd > const &inputs)
 Evaluate the log-density. More...
 
virtual double LogDensity (std::vector< Eigen::VectorXd > const &inputs)
 
 VARIADIC_TO_REFVECTOR (LogDensity, Eigen::VectorXd, double)
 
virtual Eigen::VectorXd GradLogDensity (unsigned int wrt, std::vector< Eigen::VectorXd > const &inputs)
 
virtual Eigen::VectorXd GradLogDensity (unsigned int wrt, ref_vector< Eigen::VectorXd > const &inputs)
 
template<typename... Args>
Eigen::VectorXd GradLogDensity (unsigned int wrt, Args... args)
 
virtual Eigen::VectorXd ApplyLogDensityHessian (unsigned int const inWrt1, unsigned int const inWrt2, ref_vector< Eigen::VectorXd > const &input, Eigen::VectorXd const &vec)
 
virtual Eigen::VectorXd ApplyLogDensityHessian (unsigned int const inWrt1, unsigned int const inWrt2, std::vector< Eigen::VectorXd > const &input, Eigen::VectorXd const &vec)
 
Eigen::VectorXd Sample (ref_vector< Eigen::VectorXd > const &inputs)
 Sample the distribution. More...
 
Eigen::VectorXd Sample (std::vector< Eigen::VectorXd > const &inputs)
 
Eigen::VectorXd Sample ()
 Sample the distribution with no inputs. More...
 
 VARIADIC_TO_REFVECTOR (Sample, Eigen::VectorXd, Eigen::VectorXd)
 
std::shared_ptr< DensityAsDensity ()
 Returns a density built from this distribution. More...
 
std::shared_ptr< RandomVariableAsVariable ()
 Returns a random variable built from this distribution. More...
 

Additional Inherited Members

- Public Attributes inherited from muq::Modeling::Distribution
const unsigned int varSize
 
const Eigen::VectorXi hyperSizes
 

Member Typedef Documentation

◆ InputMask

Definition at line 35 of file Gaussian.h.

Member Enumeration Documentation

◆ ExtraInputs

Enumerator
None 
Mean 
DiagCovariance 
DiagPrecision 
FullCovariance 
FullPrecision 

Definition at line 27 of file Gaussian.h.

◆ Mode

Are we specifying the mean, covariance matrix, or precision matrix.

Enumerator
Covariance 

We are specifying the covariance.

Precision 

We are specifying the precision.

Definition at line 18 of file Gaussian.h.

Constructor & Destructor Documentation

◆ Gaussian() [1/3]

Gaussian::Gaussian ( unsigned int  dim,
InputMask  extraInputs = ExtraInputs::None 
)

Definition at line 8 of file Gaussian.cpp.

References ComputeNormalization().

◆ Gaussian() [2/3]

Gaussian::Gaussian ( Eigen::VectorXd const &  mu,
InputMask  extraInputs = ExtraInputs::None 
)

Construct a Gaussian with scaled identity covariance/precision.

Parameters
[in]muThe mean

Definition at line 17 of file Gaussian.cpp.

References ComputeNormalization().

◆ Gaussian() [3/3]

Gaussian::Gaussian ( Eigen::VectorXd const &  mu,
Eigen::MatrixXd const &  obj,
Gaussian::Mode  mode = Gaussian::Mode::Covariance,
InputMask  extraInputs = ExtraInputs::None 
)

Construct a Gaussian by specifying both the mean and covariance or precision matrix.

Parameters
[in]muThe mean vector
[in]objA matrix holding either the covariance or precision
[in]modeA flag indicating whether the matrix should be treated as the covariance or precision

Definition at line 27 of file Gaussian.cpp.

References CheckInputTypes(), ComputeNormalization(), covPrec, muq::Modeling::GaussianBase::mean, mode, and sqrtCovPrec.

◆ ~Gaussian()

virtual muq::Modeling::Gaussian::~Gaussian ( )
virtualdefault

Member Function Documentation

◆ ApplyCovariance()

Eigen::MatrixXd Gaussian::ApplyCovariance ( Eigen::Ref< const Eigen::MatrixXd > const &  x) const
overridevirtual

Applies the covariance matrix to a vector \(x\)

Parameters
[in]xA reference to the input vector
Returns
\(\Sigma x\), where \(\Sigma\) is the covariance matrix.

Implements muq::Modeling::GaussianBase.

Definition at line 248 of file Gaussian.cpp.

References covPrec, mode, and sqrtCovPrec.

Referenced by Condition(), muq::SamplingAlgorithms::SMMALAProposal::LogDensity(), and muq::SamplingAlgorithms::SMMALAProposal::Sample().

◆ ApplyCovSqrt()

Eigen::MatrixXd Gaussian::ApplyCovSqrt ( Eigen::Ref< const Eigen::MatrixXd > const &  x) const
overridevirtual

Applies a matrix square root of the covariance to a vector x.

Parameters
[in]xA reference to the input vector
Returns
\(\Sigma^{1/2} x\), where \(\Sigma^{1/2}\) is a square root of the covariance matrix (e.g., Cholesky factor).

Implements muq::Modeling::GaussianBase.

Definition at line 215 of file Gaussian.cpp.

References covPrec, mode, and sqrtCovPrec.

◆ ApplyPrecision()

Eigen::MatrixXd Gaussian::ApplyPrecision ( Eigen::Ref< const Eigen::MatrixXd > const &  x) const
overridevirtual

Applies the precision matrix (inverse covariance) to a vector \(x\)

Parameters
[in]xA reference to the input vector
Returns
\(\Sigma^{-1} x\), where \(\Sigma\) is the covariance matrix.

Implements muq::Modeling::GaussianBase.

Definition at line 264 of file Gaussian.cpp.

References covPrec, mode, and sqrtCovPrec.

◆ ApplyPrecSqrt()

Eigen::MatrixXd Gaussian::ApplyPrecSqrt ( Eigen::Ref< const Eigen::MatrixXd > const &  x) const
overridevirtual

Applies a matrix square root of the precision to a vector x. If \(\Sigma= L L^T\) for a covariance matrix \(\Sigma\) and a square root \(L\), then \(\Sigma^{-1} = L^{-T}L^{-1}\) is a decomposition of the precision matrix and \(L^{-T}\) is a square root of the precision matrix. This function returns the action of \(L^{-T}\) on a vector \(x\).

Parameters
[in]xA reference to the input vector
Returns
\(L^{-T} x\)

Implements muq::Modeling::GaussianBase.

Definition at line 231 of file Gaussian.cpp.

References covPrec, mode, and sqrtCovPrec.

◆ CheckInputTypes()

void Gaussian::CheckInputTypes ( InputMask  extraInputs,
Mode  mode 
)
staticprotected

Definition at line 48 of file Gaussian.cpp.

References mode.

Referenced by Gaussian().

◆ ComputeNormalization()

void Gaussian::ComputeNormalization ( )
protected

Compute the distribution's scaling constant.

Returns
Scaling constant

Definition at line 126 of file Gaussian.cpp.

References covPrec, logDet, mode, and sqrtCovPrec.

Referenced by Gaussian(), ResetHyperparameters(), SetCovariance(), and SetPrecision().

◆ Condition()

std::shared_ptr< Gaussian > Gaussian::Condition ( Eigen::MatrixXd const &  obsMat,
Eigen::VectorXd const &  data,
Eigen::MatrixXd const &  obsCov 
) const

Returns a new Gaussian distribution conditioned on a linear observation.

Definition at line 282 of file Gaussian.cpp.

References ApplyCovariance(), muq::Modeling::GaussianBase::Dimension(), GetCovariance(), muq::Modeling::GaussianBase::mean, and nlohmann::to_string().

◆ GetCovariance()

Eigen::MatrixXd Gaussian::GetCovariance ( ) const

Get the covariance.

Returns
The covariance

Definition at line 91 of file Gaussian.cpp.

References covPrec, muq::Modeling::GaussianBase::mean, and mode.

Referenced by Condition().

◆ GetExtraSizes()

Eigen::VectorXi Gaussian::GetExtraSizes ( unsigned  dim,
InputMask  extraInputs 
)
staticprotected

Definition at line 65 of file Gaussian.cpp.

◆ GetInputTypes()

Gaussian::InputMask muq::Modeling::Gaussian::GetInputTypes ( ) const
inline

Definition at line 88 of file Gaussian.h.

References inputTypes.

◆ GetMode()

Mode muq::Modeling::Gaussian::GetMode ( ) const
inline

Definition at line 61 of file Gaussian.h.

References mode.

◆ GetPrecision()

Eigen::MatrixXd Gaussian::GetPrecision ( ) const

Definition at line 108 of file Gaussian.cpp.

References covPrec, muq::Modeling::GaussianBase::mean, and mode.

◆ LogDeterminant()

virtual double muq::Modeling::Gaussian::LogDeterminant ( ) const
inlineoverridevirtual

Return the log determinant of the covariance matrix. Defaults to 0.0, so this should be overridden by child classes that are used in applications where the density needs to known completely, not only up to a normalizing constant.

Reimplemented from muq::Modeling::GaussianBase.

Definition at line 97 of file Gaussian.h.

References logDet.

◆ ModeFromExtras()

Gaussian::Mode Gaussian::ModeFromExtras ( InputMask  extraInputs)
staticprotected

Definition at line 56 of file Gaussian.cpp.

◆ ResetHyperparameters()

void Gaussian::ResetHyperparameters ( ref_vector< Eigen::VectorXd > const &  params)
overridevirtual

Process the hyperparameter inputs. This should be overridden by children that have extra inputs.

Parameters
[in]paramsA vector of hyperparameter vectors.

Reimplemented from muq::Modeling::GaussianBase.

Definition at line 153 of file Gaussian.cpp.

References ComputeNormalization(), covPrec, inputTypes, muq::Modeling::GaussianBase::mean, sqrtCovPrec, and nlohmann::to_string().

◆ SetCovariance()

void Gaussian::SetCovariance ( Eigen::MatrixXd const &  newCov)

Set the covariance matrix.

Parameters
[in]newcovThe new covariance

Definition at line 182 of file Gaussian.cpp.

References ComputeNormalization(), covPrec, muq::Modeling::GaussianBase::mean, mode, and sqrtCovPrec.

◆ SetPrecision()

void Gaussian::SetPrecision ( Eigen::MatrixXd const &  newPrec)

Set the precision matrix.

Parameters
[in]newprecThe new precision

Definition at line 199 of file Gaussian.cpp.

References ComputeNormalization(), covPrec, muq::Modeling::GaussianBase::mean, mode, and sqrtCovPrec.

Member Data Documentation

◆ covPrec

Eigen::MatrixXd muq::Modeling::Gaussian::covPrec
protected

◆ inputTypes

Gaussian::InputMask muq::Modeling::Gaussian::inputTypes
protected

What form do the extra inputs take? Just the mean, or the mean and covariance?

Definition at line 116 of file Gaussian.h.

Referenced by GetInputTypes(), and ResetHyperparameters().

◆ logDet

double muq::Modeling::Gaussian::logDet
protected

The log determinant of the covariance matrix.

Definition at line 125 of file Gaussian.h.

Referenced by ComputeNormalization(), and LogDeterminant().

◆ mode

Gaussian::Mode muq::Modeling::Gaussian::mode
protected

◆ sqrtCovPrec

Eigen::LLT<Eigen::MatrixXd> muq::Modeling::Gaussian::sqrtCovPrec
protected

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