IMP logo
IMP Reference Guide  2.19.0
The Integrative Modeling Platform
IMP::isd::MultivariateFNormalSufficient Class Reference

MultivariateFNormalSufficient. More...

#include <IMP/isd/MultivariateFNormalSufficient.h>

+ Inheritance diagram for IMP::isd::MultivariateFNormalSufficient:

Detailed Description

MultivariateFNormalSufficient.

Probability density function and -log(p) of multivariate normal distribution of N M-variate observations.

\[ p(x_1,\cdots,x_N|\mu,F,\Sigma) = \left((2\pi\sigma^2)^M|\Sigma|\right)^{-N/2} J(F) \exp\left(-\frac{1}{2\sigma^2} \sum_{i=1}^N {}^t(F(\mu) - F(x_i))\Sigma^{-1}(F(\mu)-F(x_i)) \right) \]

which is implemented as

\[ p(x_1,\cdots,x_N|\mu,F,\Sigma) = ((2\pi\sigma^2)^M|\Sigma|)^{-N/2} J(F) \exp\left(-\frac{N}{2\sigma^2} {}^t\epsilon \Sigma^{-1} \epsilon\right) \exp\left(-\frac{1}{2\sigma^2} \text{tr}(W\Sigma^{-1})\right) \]

where

\[\epsilon = (F(\mu)- \overline{F(x)}) \quad \overline{F(x)} = \frac{1}{N} \sum_{i=1}^N F(x_i)\]

and

\[W=\sum_{i=1}^N(F(x_i)-\overline{F(x)}){}^t(F(x_i)-\overline{F(x)}) \]

\(\sigma\) is a multiplicative scaling factor that factors out of the \(\Sigma\) covariance matrix. It is set to 1 by default and its intent is to avoid inverting the \(\Sigma\) matrix unless necessary.

Set J(F) to 1 if you want the multivariate normal distribution. The distribution is normalized with respect to the matrix variable X. The Sufficient statistics are calculated at initialization.

Example: if F is the log function, the multivariate F-normal distribution is the multivariate lognormal distribution with mean \(\mu\) and standard deviation \(\Sigma\).

Note
This is an implementation of the matrix normal distribution for F(X), where rows of F(X) are independent and homoscedastic (they represent repetitions of the same experiment), but columns might be correlated, though the provided matrix.
For now, F must be monotonically increasing, so that J(F) > 0. The program will not check for that. Uses a Cholesky ( \(LDL^T\)) decomposition of \(\Sigma\), which is recomputed when needed.
All observations must be given, so if you have missing data you might want to do some imputation on it first.
References:
  • Multivariate Likelihood: Box and Tiao, "Bayesian Inference in Statistical Analysis", Addison-Wesley publishing company, 1973, pp 423.
  • Factorization in terms of sufficient statistics: Daniel Fink, "A Compendium of Conjugate Priors", online, May 1997, p.40-41.
  • Matrix calculations for derivatives: Petersen and Pedersen, "The Matrix Cookbook", 2008, matrixcookbook.com
  • Useful reading on missing data (for the case of varying Nobs): Little and Rubin, "Statistical Analysis with Missing Data", 2nd ed, Wiley, 2002, Chapters 6,7 and 11.

Definition at line 86 of file MultivariateFNormalSufficient.h.

Public Member Functions

 MultivariateFNormalSufficient (const Eigen::MatrixXd &FX, double JF, const Eigen::VectorXd &FM, const Eigen::MatrixXd &Sigma, double factor=1)
 Initialize with all observed data. More...
 
 MultivariateFNormalSufficient (const Eigen::VectorXd &Fbar, double JF, const Eigen::VectorXd &FM, int Nobs, const Eigen::MatrixXd &W, const Eigen::MatrixXd &Sigma, double factor=1)
 Initialize with sufficient statistics. More...
 
double density () const
 Probability density function. More...
 
double evaluate () const
 Energy (score) function, aka -log(p) More...
 
double evaluate_derivative_factor () const
 derivative wrt scalar factor More...
 
Eigen::VectorXd evaluate_derivative_FM () const
 gradient of the energy wrt the mean F(M) More...
 
Eigen::MatrixXd evaluate_derivative_Sigma () const
 gradient of the energy wrt the variance-covariance matrix Sigma More...
 
Eigen::MatrixXd evaluate_second_derivative_FM_FM () const
 second derivative wrt FM and FM More...
 
Eigen::MatrixXd evaluate_second_derivative_FM_Sigma (unsigned l) const
 second derivative wrt FM(l) and Sigma More...
 
Eigen::MatrixXd evaluate_second_derivative_Sigma_Sigma (unsigned k, unsigned l) const
 second derivative wrt Sigma and Sigma(k,l) More...
 
Eigen::VectorXd get_epsilon () const
 
double get_factor () const
 
Eigen::VectorXd get_Fbar () const
 
Eigen::VectorXd get_FM () const
 
Eigen::MatrixXd get_FX () const
 
double get_jacobian () const
 
double get_log_generalized_variance () const
 return log |^2 | More...
 
double get_mean_square_residuals () const
 return transpose(epsilon)*P*epsilon More...
 
double get_minus_exponent () const
 return minus exponent More...
 
double get_minus_log_jacobian () const
 
double get_minus_log_normalization () const
 return minus log normalization More...
 
Eigen::MatrixXd get_Sigma () const
 
double get_Sigma_condition_number () const
 return Sigma's condition number More...
 
Eigen::VectorXd get_Sigma_eigenvalues () const
 return Sigma's eigenvalues from smallest to biggest More...
 
virtual std::string get_type_name () const override
 
virtual ::IMP::VersionInfo get_version_info () const override
 Get information about the module and version of the object. More...
 
Eigen::MatrixXd get_W () const
 
void reset_flags ()
 if you want to force a recomputation of all stored variables More...
 
void set_factor (double f)
 
void set_Fbar (const Eigen::VectorXd &f)
 
void set_FM (const Eigen::VectorXd &f)
 
void set_FX (const Eigen::MatrixXd &f)
 
void set_jacobian (double f)
 
void set_minus_log_jacobian (double f)
 
void set_Sigma (const Eigen::MatrixXd &f)
 
void set_use_cg (bool use, double tol)
 use conjugate gradients (default false) More...
 
void set_W (const Eigen::MatrixXd &f)
 
Eigen::MatrixXd solve (Eigen::MatrixXd B) const
 Solve for Sigma*X = B, yielding X. More...
 
- Public Member Functions inherited from IMP::Object
virtual void clear_caches ()
 
CheckLevel get_check_level () const
 
LogLevel get_log_level () const
 
void set_check_level (CheckLevel l)
 
void set_log_level (LogLevel l)
 Set the logging level used in this object. More...
 
void set_was_used (bool tf) const
 
void show (std::ostream &out=std::cout) const
 
const std::string & get_name () const
 
void set_name (std::string name)
 

Additional Inherited Members

- Protected Member Functions inherited from IMP::Object
 Object (std::string name)
 Construct an object with the given name. More...
 
virtual void do_destroy ()
 

Constructor & Destructor Documentation

IMP::isd::MultivariateFNormalSufficient::MultivariateFNormalSufficient ( const Eigen::MatrixXd &  FX,
double  JF,
const Eigen::VectorXd &  FM,
const Eigen::MatrixXd &  Sigma,
double  factor = 1 
)

Initialize with all observed data.

Parameters
[in]FXF(X) matrix of observations with M columns and N rows.
[in]JFJ(F) determinant of Jacobian of F with respect to observation matrix X.
[in]FMF(M) mean vector \(F(\mu)\) of size M.
[in]Sigma: MxM variance-covariance matrix \(\Sigma\).
[in]factor: multiplicative factor (default 1)
IMP::isd::MultivariateFNormalSufficient::MultivariateFNormalSufficient ( const Eigen::VectorXd &  Fbar,
double  JF,
const Eigen::VectorXd &  FM,
int  Nobs,
const Eigen::MatrixXd &  W,
const Eigen::MatrixXd &  Sigma,
double  factor = 1 
)

Initialize with sufficient statistics.

Parameters
[in]FbarM-dimensional vector of mean observations.
[in]JFJ(F) determinant of Jacobian of F with respect to observation matrix X.
[in]FMF(M) : M-dimensional true mean vector \(\mu\).
[in]Nobsnumber of observations for each variable.
[in]WMxM matrix of sample variance-covariances.
[in]SigmaMxM variance-covariance matrix Sigma.
[in]factormultiplicative factor (default 1)

Member Function Documentation

double IMP::isd::MultivariateFNormalSufficient::density ( ) const

Probability density function.

double IMP::isd::MultivariateFNormalSufficient::evaluate ( ) const

Energy (score) function, aka -log(p)

double IMP::isd::MultivariateFNormalSufficient::evaluate_derivative_factor ( ) const

derivative wrt scalar factor

Eigen::VectorXd IMP::isd::MultivariateFNormalSufficient::evaluate_derivative_FM ( ) const

gradient of the energy wrt the mean F(M)

Eigen::MatrixXd IMP::isd::MultivariateFNormalSufficient::evaluate_derivative_Sigma ( ) const

gradient of the energy wrt the variance-covariance matrix Sigma

Eigen::MatrixXd IMP::isd::MultivariateFNormalSufficient::evaluate_second_derivative_FM_FM ( ) const

second derivative wrt FM and FM

Eigen::MatrixXd IMP::isd::MultivariateFNormalSufficient::evaluate_second_derivative_FM_Sigma ( unsigned  l) const

second derivative wrt FM(l) and Sigma

row and column indices in the matrix returned are for Sigma

Eigen::MatrixXd IMP::isd::MultivariateFNormalSufficient::evaluate_second_derivative_Sigma_Sigma ( unsigned  k,
unsigned  l 
) const

second derivative wrt Sigma and Sigma(k,l)

double IMP::isd::MultivariateFNormalSufficient::get_log_generalized_variance ( ) const

return log |^2 |

double IMP::isd::MultivariateFNormalSufficient::get_mean_square_residuals ( ) const

return transpose(epsilon)*P*epsilon

double IMP::isd::MultivariateFNormalSufficient::get_minus_exponent ( ) const

return minus exponent

\[-\frac{1}{2\sigma^2} \sum_{i=1}^N {}^t(F(\mu) - F(x_i))\Sigma^{-1}(F(\mu)-F(x_i)) \]

double IMP::isd::MultivariateFNormalSufficient::get_minus_log_normalization ( ) const

return minus log normalization

\[\frac{N}{2}\left(\log(2\pi\sigma^2) + \log |\Sigma|\right) -\log J(F) \]

double IMP::isd::MultivariateFNormalSufficient::get_Sigma_condition_number ( ) const

return Sigma's condition number

Eigen::VectorXd IMP::isd::MultivariateFNormalSufficient::get_Sigma_eigenvalues ( ) const

return Sigma's eigenvalues from smallest to biggest

virtual ::IMP::VersionInfo IMP::isd::MultivariateFNormalSufficient::get_version_info ( ) const
overridevirtual

Get information about the module and version of the object.

Reimplemented from IMP::Object.

Definition at line 225 of file MultivariateFNormalSufficient.h.

void IMP::isd::MultivariateFNormalSufficient::reset_flags ( )

if you want to force a recomputation of all stored variables

void IMP::isd::MultivariateFNormalSufficient::set_use_cg ( bool  use,
double  tol 
)

use conjugate gradients (default false)

Eigen::MatrixXd IMP::isd::MultivariateFNormalSufficient::solve ( Eigen::MatrixXd  B) const

Solve for Sigma*X = B, yielding X.


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