Multivariate Analysis Package  0.3
A multivariate distributions analysis package
 All Classes Namespaces Functions Pages
Public Member Functions | Friends | List of all members
Multivariate::NormalDistribution Class Reference

Multivariate Normal Distribution More...

#include <NormalDist.h>

Inheritance diagram for Multivariate::NormalDistribution:
Multivariate::AbstarctDistribution

Public Member Functions

 NormalDistribution (unsigned int Dimension, const Eigen::VectorXd &mVect, const Eigen::MatrixXd &CovMatr)
 Construct a multivariate normal with the given parameters. More...
 
 NormalDistribution (unsigned int Dimension=1U)
 Construct a multivariate standard normal. More...
 
Eigen::MatrixXd ExtractSamples (unsigned int NumSamples) const
 Generates multiple simulations from the distribution. More...
 
Eigen::MatrixXd ExtractSamplesCDF (unsigned int NumSamples) const
 Extracts samples from the distribution and returns their marginal CDF. More...
 
Eigen::MatrixXd GetCorrelationMatrix () const
 Get the linear correlation matrix. More...
 
double GetCumulativeDesity (const Eigen::VectorXd &Coordinates) const
 Computes the cumulative density function of the distribution in correspondence of the supplied coordinates. More...
 
double GetDensity (const Eigen::VectorXd &Coordinates) const
 Computes the probability density function of the distribution in correspondence of the supplied coordinates. More...
 
const Eigen::VectorXd & GetMeanVector () const
 Get the mean vector of the distribution. More...
 
unsigned int GetNumSimul () const
 Get the maximum number of simulations that will be used to compute the CDF. More...
 
Eigen::VectorXd GetQuantile (double Prob) const
 Computes the inverse cumulative density function of the distribution in correspondence of the supplied probability. More...
 
bool GetUseGenz () const
 Get either the Genz algorithm is used or not. More...
 
const Eigen::MatrixXd & GetVarCovMatrix () const
 Get the Var-Cov matrix of the distribution. More...
 
bool SetDimension (unsigned int Dimension)
 Set the dimensionality of the distribution. More...
 
bool SetMeanVector (const Eigen::VectorXd &mVect)
 Set the mean vector. More...
 
bool SetMeanVector (const std::vector< double > &mVect)
 Set the expected values vector. More...
 
bool SetNumSimul (unsigned int a)
 Set the maximum number of simulations used in the CDF calculation. More...
 
void SetUseGenz (bool a=true)
 Set if Genz algorithm should be used. More...
 
bool SetVarCovMatrix (const Eigen::MatrixXd &CovMatr)
 Set the Var-Cov matrix of the distribution. More...
 
bool SetVarCovMatrix (const std::vector< double > &mVect, bool RowWise=true)
 Set the Var-Cov matrix of the distribution. More...
 
bool SetVarCovMatrix (const Eigen::MatrixXd &CorrelationMatrix, const Eigen::VectorXd &Variances)
 Set the Var-Cov matrix of the distribution. More...
 
Unsafe Methods

The methods in this group use unsafe memory access or return arrays allocated on the heap that must be manually deleted.

These functions are normally not compiled for safety reasons. To use them, the mvPackageUnsafeMethods symbol must be defined at compile time

virtual void SetMeanVector (double *mVect)
 Set the expected values vector. More...
 
virtual bool SetVarCovMatrix (double **mVect)
 Set the Var-Cov matrix of the distribution. More...
 
- Public Member Functions inherited from Multivariate::AbstarctDistribution
virtual Eigen::RowVectorXd ExtractSample () const
 Generates a single simulation from the distribution. More...
 
virtual Eigen::RowVectorXd ExtractSampleCDF () const
 Generates a single simulation from the distribution and returns its marginal CDF. More...
 
virtual std::vector< double > ExtractSampleCDFVect () const
 Generates a single simulation from the distribution and returns its marginal CDF. More...
 
virtual std::map< unsigned int,
std::vector< double > > 
ExtractSamplesCDFMap (unsigned int NumSamples) const
 Extracts samples from the distribution and returns their marginal CDF. More...
 
virtual std::map< unsigned int,
std::vector< double > > 
ExtractSamplesMap (unsigned int NumSamples) const
 Generates multiple simulations from the distribution. More...
 
virtual std::vector< double > ExtractSampleVector () const
 Generates a single simulation from the distribution. More...
 
virtual double GetCumulativeDesity (const std::vector< double > &Coordinates) const
 Computes the cumulative density function of the distribution in correspondence of the supplied coordinates. More...
 
unsigned int GetCurrentSeed () const
 Get the random number generator seed. More...
 
virtual double GetDensity (const std::vector< double > &Coordinates) const
 Computes the probability density function of the distribution in correspondence of the supplied coordinates. More...
 
unsigned int GetDimension () const
 Get the dimensionality of the distribution. More...
 
virtual std::vector< double > GetQuantileVector (double Prob) const
 Computes the inverse cumulative density function of the distribution in correspondence of the supplied probability. More...
 
bool IsValid () const
 Check if the distribution is valid. More...
 
void SetRandomSeed (unsigned int NewSeed)
 Set the random number generator seed. More...
 
virtual double * GetQuantileArray (double Prob)
 Computes the inverse cumulative density function of the distribution in correspondence of the supplied probability. More...
 
virtual double GetCumulativeDesity (double *Coordinates) const
 Computes the cumulative density function of the distribution in correspondence of the supplied coordinates. More...
 
virtual double GetDensity (double *Coordinates) const
 Computes the probability density function of the distribution in correspondence of the supplied coordinates. More...
 
virtual double * ExtractSampleArray () const
 Generates a single simulation from the distribution. More...
 
virtual double ** ExtractSamplesMatix (unsigned int NumSamples) const
 Generates a single simulation from the distribution. More...
 
virtual double * ExtractSampleCDFArray () const
 Generates a single simulation from the distribution and returns its marginal CDF. More...
 
virtual double ** ExtractSamplesCDFMatix (unsigned int NumSamples) const
 Generates a single simulation from the distribution and returns its marginal CDF. More...
 

Friends

template<class F , class T >
void boost::math::tools::detail::handle_zero_derivative (F f, T &last_f0, const T &f0, T &delta, T &result, T &guess, const T &min, const T &max)
 
template<class F , class T >
boost::math::tools::newton_raphson_iterate (F f, T guess, T min, T max, int digits)
 
template<class F , class T >
boost::math::tools::newton_raphson_iterate (F f, T guess, T min, T max, int digits, boost::uintmax_t &max_iter)
 

Detailed Description

Multivariate Normal Distribution

This class provides the functionality of calculating the probability density value, cumulative probability density value, inverse cumulative probability density and generate random samples from a multivariate normal.

Defining:

The multivariate normal distribution funtion is defined as: \( f(\textbf{x})=((2\pi)^{-\frac{k}{2}} |\boldsymbol{\Sigma}|^{-\frac{1}{2}} e^{(-\frac{1}{2} (\textbf{x}-\boldsymbol{\mu})' \boldsymbol{\Sigma}^{-1} (\textbf{x}-\boldsymbol{\mu}))} \)

The analytical process for computing pdf and simulate from the distribution is based upon the mvtnorm package for R by Alan Genz, Frank Bretz, Tetsuhisa Miwa, Xuefei Mi, Friedrich Leisch, Fabian Scheipl, Bjoern Bornkamp, Torsten Hothorn The algorithm for cdf calculation is based on A. Genz (1992)

If you construct multiple instances of this class, to avoid the generated samples to be the same, you should supply a different seed. To do so, for example, you can call MyDistribution.SetRandomSeed(MyDistribution.GetCurrentSeed()+1U);

Please refer to the Examples page for usage examples.

Remarks
This class is re-entrant
Date
November 2013
License
This program is free software: you can redistribute it and/or modify it under the terms of the GNU Lesser General Public License as published by the Free Software Foundation, either version 3 of the License, or (at your option) any later version.

This program is distributed in the hope that it will be useful, but WITHOUT ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the GNU Lesser General Public License for more details.

Here, you can find a copy of the GNU Lesser General Public License. Alternatively, see gnu.org.

Constructor & Destructor Documentation

NormalDistribution::NormalDistribution ( unsigned int  Dimension,
const Eigen::VectorXd &  mVect,
const Eigen::MatrixXd &  CovMatr 
)

Construct a multivariate normal with the given parameters.

Parameters
DimensionThe dimensionality of the multivariate normal (supports also univariate gaussian distributions in case this is 1)
mVectThe column vector of expected values
CovMatrThe variance-covariance

Construct a multivariate normal distribution.

In case:

  • The Dimension is 0
  • The mean vector has a number of elements different from the Dimension
  • The variance-covariance is not square
  • The variance-covariance is not symmetric
  • The variance-covariance is not semi-positive definite
  • The variance-covariance has a number of rows different from the Dimension

The class will be considered invalid (it can be checked using IsValid()) and won't produce any result until the problem is fixed.

NormalDistribution::NormalDistribution ( unsigned int  Dimension = 1U)

Construct a multivariate standard normal.

Parameters
DimensionThe dimensionality of the multivariate normal (supports also univariate Gaussian distributions in case this is 1)

Construct a multivariate normal distribution with all mean values set to 0, unitary variances and null covariates.

In case The Dimension is 0 the class will be considered invalid (it can be checked using IsValid()) and won't produce any result until the problem is fixed.

If dimension is unspecified, a univariate standard normal is constructed

Member Function Documentation

Eigen::MatrixXd NormalDistribution::ExtractSamples ( unsigned int  NumSamples) const
virtual

Generates multiple simulations from the distribution.

Parameters
NumSamplesThe number of simulation to run
Returns
A matrix with columns equal to the dimensionality of the distribution and rows equal to the number of simulations

This function generates NumSamples simulation from the current distribution and returns them in matrix form.

If NumSamples is 0 or the distribution is invalid, a null matrix is returned

Implements Multivariate::AbstarctDistribution.

Eigen::MatrixXd NormalDistribution::ExtractSamplesCDF ( unsigned int  NumSamples) const
virtual

Extracts samples from the distribution and returns their marginal CDF.

Parameters
NumSamplesThe number of simulation to run
Returns
A matrix with columns equal to the dimensionality of the distribution and rows equal to the number of simulations

This function generates NumSamples simulation from the current distribution, computes the marginal cumulative density function for each of them and returns them in matrix form.

If NumSamples is 0 or the distribution is invalid, a null matrix is returned

Implements Multivariate::AbstarctDistribution.

Eigen::MatrixXd NormalDistribution::GetCorrelationMatrix ( ) const

Get the linear correlation matrix.

Returns
The linear correlation matrix associated with the current variance-covariance matrix of the distribution
See Also
GetVarCovMatrix()
double NormalDistribution::GetCumulativeDesity ( const Eigen::VectorXd &  Coordinates) const
virtual

Computes the cumulative density function of the distribution in correspondence of the supplied coordinates.

Parameters
CoordinatesA vector containing the coordinates of the point for which the cdf should be computed
Returns
The value of the cumulative density function

This function computes the cumulative density function of the current distribution associated with the given coordinates.

If the number of elements in Coordinates is different from the dimensionality of the distribution or the distribution is invalid, -1 is returned.

Implements Multivariate::AbstarctDistribution.

double NormalDistribution::GetDensity ( const Eigen::VectorXd &  Coordinates) const
virtual

Computes the probability density function of the distribution in correspondence of the supplied coordinates.

Parameters
CoordinatesA vector containing the coordinates of the point for which the pdf should be computed
Returns
The value of the probability density function

This function computes the probability density function of the current distribution associated with the given coordinates.

If the number of elements in Coordinates is different from the dimensionality of the distribution or the distribution is invalid, -1 is returned

Implements Multivariate::AbstarctDistribution.

const Eigen::VectorXd& Multivariate::NormalDistribution::GetMeanVector ( ) const
inline

Get the mean vector of the distribution.

Returns
The current mean vector of the distribution
See Also
SetMeanVector(const Eigen::VectorXd&)
unsigned int Multivariate::NormalDistribution::GetNumSimul ( ) const
inline

Get the maximum number of simulations that will be used to compute the CDF.

Returns
The maximum number of simulations
See Also
SetNumSimul();
Eigen::VectorXd NormalDistribution::GetQuantile ( double  Prob) const
virtual

Computes the inverse cumulative density function of the distribution in correspondence of the supplied probability.

Parameters
ProbThe probability for which the corresponding quantile must be found
Returns
A vector containing the coordinates of the quantile

This function computes the inverse cumulative density function of the current distribution associated with the given probability.

The solution is unique only in the univariate case.
Generally the system of equations \( F^{-1}(Coordinates_1 \cdots Coordinates_k)=Prob \) has k-1 degrees of freedom, where k is the dimensionality of the distribution.
The additional restriction imposed to get to an unique solution is that each coordinate has equal distance from it's mean.

If the probability supplied is greater than 1, less than 0 or the distribution is invalid, an empty vector is returned.

Implements Multivariate::AbstarctDistribution.

bool Multivariate::NormalDistribution::GetUseGenz ( ) const
inline

Get either the Genz algorithm is used or not.

Returns
A bool determining if Genz algorithm will be used for calculating the CDF
See Also
SetUseGenz();
const Eigen::MatrixXd& Multivariate::NormalDistribution::GetVarCovMatrix ( ) const
inline

Get the Var-Cov matrix of the distribution.

Returns
The current variance-covariance matrix of the distribution
See Also
SetVarCovMatrix(const Eigen::MatrixXd&)
SetVarCovMatrix(const std::vector<double>&,bool)
SetVarCovMatrix(const Eigen::MatrixXd&, const Eigen::VectorXd&)
bool NormalDistribution::SetDimension ( unsigned int  Dimension)
virtual

Set the dimensionality of the distribution.

Parameters
Dimensionthe new dimensionality of the distribution
Returns
A boolean determining if the dimensionality was changed successfully

This function will try to change the dimensionality of the distribution (e.g. 2 for bivariate, 3 for trivariate, etc.)

All the components of the mean vector will be reset to 0 and the variance covariance matrix will default to an identity matrix

If the argument passed is 0 the dimensionality will not be changed and the function will return false

See Also
GetDimension()

Implements Multivariate::AbstarctDistribution.

bool NormalDistribution::SetMeanVector ( const Eigen::VectorXd &  mVect)

Set the mean vector.

Parameters
mVectthe vector of new values for the mean vector
Returns
A boolean determining if the mean vector was changed successfully

This function attempts to set the mean vector to the new values supplied.

If the dimension of the vector is different from the dimension of the distribution the mean vector is not changed and false is returned

See Also
GetMeanVector()
bool NormalDistribution::SetMeanVector ( const std::vector< double > &  mVect)

Set the expected values vector.

Parameters
mVectthe vector of new values for the mean vector

This is an overloaded version of SetMeanVector(const Eigen::VectorXd&)

void NormalDistribution::SetMeanVector ( double *  mVect)
virtual

Set the expected values vector.

Parameters
mVectan array containing the new values for the expected values vector

This is an overloaded version of SetMeanVector(const Eigen::VectorXd&)

Warning
This function will search for a number of elements equal to the dimensionality of the distribution in the array. This may mean accessing unallocated memory blocks if the supplied array is not big enough
bool NormalDistribution::SetNumSimul ( unsigned int  a)

Set the maximum number of simulations used in the CDF calculation.

void Multivariate::NormalDistribution::SetUseGenz ( bool  a = true)
inline

Set if Genz algorithm should be used.

Parameters
aDeterminses if Genz algorithm will be used

If UseGenz is true, the algorithm used for computing the cumulative density function and the quantiles is the one described in Alan Genz - "Numerical Computation of Multivariate Normal Probabilities", Journal of Computational and Graphical Statistics, 1(1992), pp. 141-149

If set to false, full monte-carlo simulation will be used (resulting in much slower performance)

Remarks
By default Genz algorithm is used in CDF and quantiles calculations.
Note
The GetQuantile() function needs to calculate the CDF several times, using full monte-carlo may lead to very long execution times of that function
See Also
GetUseGenz()
SetNumSimul()
bool NormalDistribution::SetVarCovMatrix ( const Eigen::MatrixXd &  CovMatr)

Set the Var-Cov matrix of the distribution.

Parameters
CovMatrthe new variance covariance matrix of the distribution
Returns
A boolean determining if the variance covariance matrix of the distribution was changed successfully

This function tries to set the variance covariance matrix of the distribution to the new one.

In case:

  • The variance-covariance is not square
  • The variance-covariance is not symmetric
  • The variance-covariance is not semi-positive definite
  • The variance-covariance has a number of rows different from the Dimension

The variance covariance matrix of the distribution will not be changed and this function will return false

See Also
GetVarCovMatrix()
bool NormalDistribution::SetVarCovMatrix ( const std::vector< double > &  mVect,
bool  RowWise = true 
)

Set the Var-Cov matrix of the distribution.

Parameters
mVecta vector containing the elements of the new variance covariance matrix of the distribution
RowWiseif it's set to true (the default) the matrix will be filled by row. If it's false it will be filled by columns
Returns
A boolean determining if the variance covariance matrix of the distribution was changed successfully

This function tries to set the variance covariance matrix of the distribution to the new one.

Constructs square a matrix with number of rows equal to the dimensionality of the distribution, it is then filled with the values supplied in order according to the RowWise parameter

In case:

  • The vector size is different from the square of the distribution dimensionality
  • The variance-covariance is not symmetric
  • The variance-covariance is not semi-positive definite
  • The variance-covariance has a number of rows different from the Dimension

The variance covariance matrix of the distribution will not be changed and this function will return false

See Also
GetVarCovMatrix()
bool NormalDistribution::SetVarCovMatrix ( const Eigen::MatrixXd &  CorrelationMatrix,
const Eigen::VectorXd &  Variances 
)

Set the Var-Cov matrix of the distribution.

Parameters
CorrelationMatrixThe correlation coefficients matrix
Variancesthe vector of variances
Returns
A boolean determining if the variance covariance matrix of the distribution was changed successfully

This function tries to set the variance covariance matrix of the distribution according to the correlation matrix supplied.

Constructs the variance covariance matrix of the distribution using Variances as diagonal elements and covariates corresponding to the linear correlation coefficient supplied

In case:

  • The size of the vector of variances is different from the distribution dimensionality
  • The correlation matrix is not symmetric
  • The correlation matrix has elements different from 1 on the diagonal
  • The correlation matrix has a number of rows different from the Dimension
  • The correlation matrix is not squared
  • The correlation matrix has off-diagonal elements greater than 1 in absolute value

The variance covariance matrix of the distribution will not be changed and this function will return false

See Also
GetVarCovMatrix()
bool NormalDistribution::SetVarCovMatrix ( double **  mVect)
virtual

Set the Var-Cov matrix of the distribution.

Parameters
mVecta matrix containing the new values for the variance covariance matrix of the distribution

This is an overloaded version of SetVarCovMatrix(const Eigen::MatrixXd&)

Warning
This function will search for a number of elements, both in the row and the column dimension, equal to the dimensionality of the distribution in the matrix. This may mean accessing unallocated memory blocks if the supplied matrix is not big enough