#include <Gaussian.h>
@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< 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. More... | |
void | ResetHyperparameters (ref_vector< Eigen::VectorXd > const ¶ms) 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< Density > | AsDensity () |
Returns a density built from this distribution. More... | |
std::shared_ptr< RandomVariable > | AsVariable () |
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 |
typedef uint8_t muq::Modeling::Gaussian::InputMask |
Definition at line 35 of file Gaussian.h.
Enumerator | |
---|---|
None | |
Mean | |
DiagCovariance | |
DiagPrecision | |
FullCovariance | |
FullPrecision |
Definition at line 27 of file Gaussian.h.
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.
Gaussian::Gaussian | ( | unsigned int | dim, |
InputMask | extraInputs = ExtraInputs::None |
||
) |
Definition at line 8 of file Gaussian.cpp.
References ComputeNormalization().
Gaussian::Gaussian | ( | Eigen::VectorXd const & | mu, |
InputMask | extraInputs = ExtraInputs::None |
||
) |
Construct a Gaussian with scaled identity covariance/precision.
[in] | mu | The mean |
Definition at line 17 of file Gaussian.cpp.
References ComputeNormalization().
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.
[in] | mu | The mean vector |
[in] | obj | A matrix holding either the covariance or precision |
[in] | mode | A 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.
|
virtualdefault |
|
overridevirtual |
Applies the covariance matrix to a vector \(x\)
[in] | x | A reference to the input vector |
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().
|
overridevirtual |
Applies a matrix square root of the covariance to a vector x.
[in] | x | A reference to the input vector |
Implements muq::Modeling::GaussianBase.
Definition at line 215 of file Gaussian.cpp.
References covPrec, mode, and sqrtCovPrec.
|
overridevirtual |
Applies the precision matrix (inverse covariance) to a vector \(x\)
[in] | x | A reference to the input vector |
Implements muq::Modeling::GaussianBase.
Definition at line 264 of file Gaussian.cpp.
References covPrec, mode, and sqrtCovPrec.
|
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\).
[in] | x | A reference to the input vector |
Implements muq::Modeling::GaussianBase.
Definition at line 231 of file Gaussian.cpp.
References covPrec, mode, and sqrtCovPrec.
|
protected |
Compute the distribution's scaling constant.
Definition at line 126 of file Gaussian.cpp.
References covPrec, logDet, mode, and sqrtCovPrec.
Referenced by Gaussian(), ResetHyperparameters(), SetCovariance(), and SetPrecision().
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().
Eigen::MatrixXd Gaussian::GetCovariance | ( | ) | const |
Get the covariance.
Definition at line 91 of file Gaussian.cpp.
References covPrec, muq::Modeling::GaussianBase::mean, and mode.
Referenced by Condition().
|
staticprotected |
Definition at line 65 of file Gaussian.cpp.
|
inline |
Definition at line 88 of file Gaussian.h.
References inputTypes.
|
inline |
Definition at line 61 of file Gaussian.h.
References mode.
Eigen::MatrixXd Gaussian::GetPrecision | ( | ) | const |
Definition at line 108 of file Gaussian.cpp.
References covPrec, muq::Modeling::GaussianBase::mean, and mode.
|
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.
|
staticprotected |
Definition at line 56 of file Gaussian.cpp.
|
overridevirtual |
Process the hyperparameter inputs. This should be overridden by children that have extra inputs.
[in] | params | A 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().
void Gaussian::SetCovariance | ( | Eigen::MatrixXd const & | newCov | ) |
Set the covariance matrix.
[in] | newcov | The new covariance |
Definition at line 182 of file Gaussian.cpp.
References ComputeNormalization(), covPrec, muq::Modeling::GaussianBase::mean, mode, and sqrtCovPrec.
void Gaussian::SetPrecision | ( | Eigen::MatrixXd const & | newPrec | ) |
Set the precision matrix.
[in] | newprec | The new precision |
Definition at line 199 of file Gaussian.cpp.
References ComputeNormalization(), covPrec, muq::Modeling::GaussianBase::mean, mode, and sqrtCovPrec.
|
protected |
Definition at line 119 of file Gaussian.h.
Referenced by ApplyCovariance(), ApplyCovSqrt(), ApplyPrecision(), ApplyPrecSqrt(), ComputeNormalization(), Gaussian(), GetCovariance(), GetPrecision(), ResetHyperparameters(), SetCovariance(), and SetPrecision().
|
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().
|
protected |
The log determinant of the covariance matrix.
Definition at line 125 of file Gaussian.h.
Referenced by ComputeNormalization(), and LogDeterminant().
|
protected |
Have we specified the covariance or the precision.
Definition at line 113 of file Gaussian.h.
Referenced by ApplyCovariance(), ApplyCovSqrt(), ApplyPrecision(), ApplyPrecSqrt(), CheckInputTypes(), ComputeNormalization(), Gaussian(), GetCovariance(), GetMode(), GetPrecision(), SetCovariance(), and SetPrecision().
|
protected |
Definition at line 122 of file Gaussian.h.
Referenced by ApplyCovariance(), ApplyCovSqrt(), ApplyPrecision(), ApplyPrecSqrt(), ComputeNormalization(), Gaussian(), ResetHyperparameters(), SetCovariance(), and SetPrecision().