Defines an abstract Gaussian class.@seealso Gaussian. More...
#include <GaussianBase.h>
Defines an abstract Gaussian class.
@seealso Gaussian.
Definition at line 15 of file GaussianBase.h.
Public Member Functions | |
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::MatrixXd | ApplyCovariance (Eigen::Ref< const Eigen::MatrixXd > const &x) const =0 |
virtual Eigen::MatrixXd | ApplyPrecision (Eigen::Ref< const Eigen::MatrixXd > const &x) const =0 |
virtual Eigen::MatrixXd | ApplyCovSqrt (Eigen::Ref< const Eigen::MatrixXd > const &x) const =0 |
virtual Eigen::MatrixXd | ApplyPrecSqrt (Eigen::Ref< const Eigen::MatrixXd > const &x) const =0 |
virtual Eigen::VectorXd const & | GetMean () const |
virtual void | SetMean (Eigen::VectorXd const &newMu) |
virtual double | LogDeterminant () const |
virtual void | ResetHyperparameters (ref_vector< Eigen::VectorXd > const ¶ms) |
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 |
GaussianBase::GaussianBase | ( | unsigned int | dim | ) |
Basic constructor for Gaussians with not additional hyperparameters.
[in] | dim | The dimension of the Gaussian distribution. |
Definition at line 8 of file GaussianBase.cpp.
GaussianBase::GaussianBase | ( | unsigned int | dim, |
Eigen::VectorXi const & | hyperSizesIn | ||
) |
Basic constructor for Gaussians with hyperparameter inputs.
[in] | dim | The dimension of the Gaussian distribution |
[in] | hyperSizesIn | A vector of integers specify the size of any additional inputs (e.g., mean, standard deviation). |
Definition at line 12 of file GaussianBase.cpp.
GaussianBase::GaussianBase | ( | Eigen::VectorXd const & | meanIn | ) |
Construct a Gaussian with no hyperparameter inputs and a specified mean vector.
[in] | meanIn | A vector containing the mean of the Gaussian distribution. |
Definition at line 19 of file GaussianBase.cpp.
GaussianBase::GaussianBase | ( | Eigen::VectorXd const & | meanIn, |
Eigen::VectorXi const & | hyperSizesIn | ||
) |
Construct a Gaussian with hyperparameters and a specified mean vector.
[in] | meanIn | A vector containing the mean of the Gaussian distribution. |
[in] | hyperSizesIn | A vector of integers specify the size of any additional inputs (e.g., mean, standard deviation). |
Definition at line 24 of file GaussianBase.cpp.
|
virtualdefault |
|
pure virtual |
Applies the covariance matrix to a vector \(x\)
[in] | x | A reference to the input vector |
Implemented in muq::Modeling::Gaussian.
|
pure virtual |
Applies a matrix square root of the covariance to a vector x.
[in] | x | A reference to the input vector |
Implemented in muq::Modeling::Gaussian.
|
overridevirtual |
Reimplemented from muq::Modeling::Distribution.
Definition at line 69 of file GaussianBase.cpp.
References ApplyPrecision().
|
pure virtual |
Applies the precision matrix (inverse covariance) to a vector \(x\)
[in] | x | A reference to the input vector |
Implemented in muq::Modeling::Gaussian.
Referenced by ApplyLogDensityHessianImpl(), and GradLogDensityImpl().
|
pure virtual |
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 |
Implemented in muq::Modeling::Gaussian.
|
virtual |
Definition at line 48 of file GaussianBase.cpp.
References mean.
Referenced by muq::Modeling::Gaussian::Condition().
|
inlinevirtual |
Returns a vector containing the mean of the Gaussian.
Definition at line 91 of file GaussianBase.h.
References mean.
|
overridevirtual |
Compute the gradient of the log density with respect to either the distribution input or the hyperparameters. This should be overridden by child classes that have hyperparameters. If not overridden, this function will assert false if wrt>0.
[in] | wrt | Specifies the index of the variable we wish to take the gradient wrt. If wrt==0, then the gradient should be taken wrt the input variable. |
Reimplemented from muq::Modeling::Distribution.
Definition at line 58 of file GaussianBase.cpp.
References ApplyPrecision(), and mean.
|
overrideprotectedvirtual |
Compute the log density.
[inputs] | A vector of extra hyperparameter vectors. |
Reimplemented from muq::Modeling::Distribution.
Definition at line 30 of file GaussianBase.cpp.
References mean, and ResetHyperparameters().
|
inlinevirtual |
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 in muq::Modeling::Gaussian.
Definition at line 105 of file GaussianBase.h.
|
inlinevirtual |
Process the hyperparameter inputs. This should be overridden by children that have extra inputs.
[in] | params | A vector of hyperparameter vectors. |
Reimplemented in muq::Modeling::Gaussian.
Definition at line 111 of file GaussianBase.h.
Referenced by LogDensityImpl(), and SampleImpl().
|
overrideprotectedvirtual |
Sample the distribution.
Reimplemented from muq::Modeling::Distribution.
Definition at line 39 of file GaussianBase.cpp.
References mean, and ResetHyperparameters().
|
virtual |
Set the mean of this Gaussian. The new mean must be the same size as the old mean.
[in] | newMu | The new mean. |
Definition at line 52 of file GaussianBase.cpp.
References mean.
|
protected |
Definition at line 140 of file GaussianBase.h.
Referenced by muq::Modeling::Gaussian::Condition(), Dimension(), muq::Modeling::Gaussian::Gaussian(), muq::Modeling::Gaussian::GetCovariance(), GetMean(), muq::Modeling::Gaussian::GetPrecision(), GradLogDensityImpl(), LogDensityImpl(), muq::Modeling::Gaussian::ResetHyperparameters(), SampleImpl(), muq::Modeling::Gaussian::SetCovariance(), SetMean(), and muq::Modeling::Gaussian::SetPrecision().