A Point/Value pair that is used to specify a value at a given point. More...
#include <StatisticalModel.h>
Public Types | |
typedef Representer< T > | RepresenterType |
typedef RepresenterType::DatasetPointerType | DatasetPointerType |
typedef RepresenterType::DatasetConstPointerType | DatasetConstPointerType |
typedef RepresenterType::ValueType | ValueType |
typedef RepresenterType::PointType | PointType |
typedef Domain< PointType > | DomainType |
typedef unsigned | PointIdType |
typedef std::pair< PointType, ValueType > | PointValuePairType |
typedef std::pair< unsigned, ValueType > | PointIdValuePairType |
typedef std::list < PointValuePairType > | PointValueListType |
typedef std::list < PointIdValuePairType > | PointIdValueListType |
Public Member Functions | |
virtual | ~StatisticalModel () |
General Info | |
unsigned int | GetNumberOfPrincipalComponents () const |
const ModelInfo & | GetModelInfo () const |
Sample from the model | |
| |
DatasetPointerType | DatasetToSample (DatasetConstPointerType dataset) const |
ValueType | EvaluateSampleAtPoint (DatasetConstPointerType sample, unsigned ptId) const |
ValueType | EvaluateSampleAtPoint (DatasetConstPointerType sample, const PointType &pt) const |
DatasetPointerType | DrawMean () const |
DatasetPointerType | DrawSample (const VectorType &coefficients, bool addNoise=false) const |
DatasetPointerType | DrawSample (bool addNoise=false) const |
DatasetPointerType | DrawPCABasisSample (unsigned componentNumber) const |
Point sampling and point information | |
ValueType | DrawMeanAtPoint (const PointType &point) const |
ValueType | DrawMeanAtPoint (unsigned pointId) const |
ValueType | DrawSampleAtPoint (const VectorType &coefficients, const PointType &point, bool addNoise=false) const |
ValueType | DrawSampleAtPoint (const VectorType &coefficients, unsigned pointId, bool addNoise=false) const |
MatrixType | GetJacobian (const PointType &pt) const |
MatrixType | GetCovarianceAtPoint (const PointType &pt1, const PointType &pt2) const |
MatrixType | GetCovarianceAtPoint (unsigned ptId1, unsigned ptId2) const |
Statistical Information from Dataset | |
MatrixType | GetCovarianceMatrix () const |
double | ComputeProbabilityOfDataset (DatasetConstPointerType dataset) const |
double | ComputeLogProbabilityOfDataset (DatasetConstPointerType dataset) const |
double | ComputeProbabilityOfCoefficients (const VectorType &coefficients) const |
double | ComputeLogProbabilityOfCoefficients (const VectorType &coefficients) const |
double | ComputeMahalanobisDistanceForDataset (DatasetConstPointerType dataset) const |
VectorType | ComputeCoefficientsForDataset (DatasetConstPointerType dataset) const |
VectorType | ComputeCoefficientsForSample (DatasetConstPointerType sample) const |
VectorType | ComputeCoefficientsForPointValues (const PointValueListType &pointValues, double pointValueNoiseVariance=0.0) const |
VectorType | ComputeCoefficientsForPointIDValues (const PointIdValueListType &pointValues, double pointValueNoiseVariance=0.0) const |
VectorType | RobustlyComputeCoefficientsForDataset (DatasetConstPointerType dataset, unsigned nIterations=100, unsigned nu=6, double sigma2=1) const |
Low level access | |
These methods provide a low level interface to the model content. They are of only limited use for an application. Prefer whenever possible the high level functions. | |
float | GetNoiseVariance () const |
const VectorType & | GetPCAVarianceVector () const |
const VectorType & | GetMeanVector () const |
const MatrixType & | GetPCABasisMatrix () const |
MatrixType | GetOrthonormalPCABasisMatrix () const |
VectorType | DrawSampleVector (const VectorType &coefficients, bool addNoise=false) const |
void | SetModelInfo (const ModelInfo &modelInfo) |
VectorType | ComputeCoefficientsForSampleVector (const VectorType &sample) const |
const RepresenterType * | GetRepresenter () const |
const DomainType & | GetDomain () const |
Loading and creating models | |
void | Delete () const |
void | Save (const std::string &filename) const |
void | Save (const H5::Group &modelRoot) const |
static StatisticalModel * | Create (const RepresenterType *representer, const VectorType &m, const MatrixType &orthonormalPCABasis, const VectorType &pcaVariance, double noiseVariance) |
static StatisticalModel * | Load (Representer< T > *representer, const std::string &filename, unsigned maxNumberOfPCAComponents=std::numeric_limits< unsigned >::max()) |
static StatisticalModel * | Load (Representer< T > *representer, const H5::Group &modelroot, unsigned maxNumberOfPCAComponents=std::numeric_limits< unsigned >::max()) |
A Point/Value pair that is used to specify a value at a given point.
Representation of a linear statistical model (PCA Model).
The statistical model class represents a statistical (PCA based) model. The implementation is based on the Probabilistic PCA, which includes the Standard PCA as a special case.
Mathematically, the statistical model is a generative model, where a sample is given as a linear combination of a dimensional latent variable
plus additive gaussian noise.
Here, is the mean,
is a linear mapping,
is a vector of latent variables (later also referred to as coefficients) and
is a noise component. The linear mapping
is referred to as the PCABasis. It is usually given as
where
where U is an orthonormal matrix and
is a scaling matrix, referred to as PCAvariance. Usually,
is the matrix of eigenvectors of the data covariance matrix and
the corresponding eigenvalues.
While all the matrices and vectors defined above could be obtained directly, the main goal of this class is to abstract from these technicalities, by providing a high level interface to shape models. In this high level view, the model represents a multivariate normal distribution over the types defined by the representer (which are typically either, surface meshes, point clouds, deformation fields or images). This class provides the method to sample from this probability distribution, and to compute the probability of given samples directly.
|
virtual |
Destructor
VectorType statismo::StatisticalModel< T >::ComputeCoefficientsForDataset | ( | DatasetConstPointerType | dataset | ) | const |
Converts the given dataset to a sample of the model and compute the latent variable coefficients of this sample under the model.
dataset | The dataset |
VectorType statismo::StatisticalModel< T >::ComputeCoefficientsForPointIDValues | ( | const PointIdValueListType & | pointValues, |
double | pointValueNoiseVariance = 0.0 |
||
) | const |
Same as ComputeCoefficientsForPointValues(const PointValueListType& pointValues), but used when the point ids, rather than the points are known.
pointValues | A list with (Point,Value) pairs, a list of (PointId, Value) is provided. |
pointValueNoiseVariance | The variance of estimated (gaussian) noise at the known points |
VectorType statismo::StatisticalModel< T >::ComputeCoefficientsForPointValues | ( | const PointValueListType & | pointValues, |
double | pointValueNoiseVariance = 0.0 |
||
) | const |
Returns the coefficients of the latent variables for the given values provided in the PointValueList. This is useful, when only a part of the dataset is given. The method is described in the paper
Probabilistic Modeling and Visualization of the Flexibility in Morphable Models, M. Luethi, T. Albrecht and T. Vetter, Mathematics of Surfaces, 2009
pointValues | A list with PointValuePairs . |
pointValueNoiseVariance | The variance of estimated (gaussian) noise at the known points |
VectorType statismo::StatisticalModel< T >::ComputeCoefficientsForSample | ( | DatasetConstPointerType | sample | ) | const |
Returns the coefficients of the latent variables for the given sample, i.e. the vectors of numbers , such that for the dataset
it holds that
VectorType statismo::StatisticalModel< T >::ComputeCoefficientsForSampleVector | ( | const VectorType & | sample | ) | const |
Computes the coefficients for the given sample vector. This is for library internal use only.
double statismo::StatisticalModel< T >::ComputeLogProbabilityOfCoefficients | ( | const VectorType & | coefficients | ) | const |
Returns the log probability of observing given coefficients.
dataset | The coefficients ![]() |
double statismo::StatisticalModel< T >::ComputeLogProbabilityOfDataset | ( | DatasetConstPointerType | dataset | ) | const |
Returns the log probability of observing a given dataset.
dataset | The dataset |
double statismo::StatisticalModel< T >::ComputeMahalanobisDistanceForDataset | ( | DatasetConstPointerType | dataset | ) | const |
Returns the mahalonoibs distance for the given dataset.
double statismo::StatisticalModel< T >::ComputeProbabilityOfCoefficients | ( | const VectorType & | coefficients | ) | const |
Returns the probability of observing the given coefficients under this model. If the coefficients define the dataset, the probability is
coefficients | The coefficients ![]() |
double statismo::StatisticalModel< T >::ComputeProbabilityOfDataset | ( | DatasetConstPointerType | dataset | ) | const |
Returns the probability of observing the given dataset under this model. If the coefficients define the dataset, the probability is
datatset | The dataset |
StatisticalModel< T >::DatasetPointerType statismo::StatisticalModel< T >::DatasetToSample | ( | DatasetConstPointerType | dataset | ) | const |
Returns a sample for the given Dataset by calling the representers internal functions. The resulting sample will have the same topology and alignment as all the other samples in the model.
The exact semantics of this call depends on the representer. It typically includes a rigid alignment of the dataset to the model but could even reparametrization or registration.
dataset |
|
inline |
Destroy the object. The same effect can be achieved by deleting the object in the usual way using the c++ delete keyword.
StatisticalModel< T >::DatasetPointerType statismo::StatisticalModel< T >::DrawMean | ( | ) | const |
StatisticalModel< T >::ValueType statismo::StatisticalModel< T >::DrawMeanAtPoint | ( | const PointType & | point | ) | const |
Returns the mean of the model, evaluated at the given point.
point | A point on the domain the model is defined |
StatisticalModel< T >::ValueType statismo::StatisticalModel< T >::DrawMeanAtPoint | ( | unsigned | pointId | ) | const |
Returns the mean of the model, evaluated at the given pointid.
pointid | The pointId of the point where it should be evaluated (as defined by the representer) |
StatisticalModel< T >::DatasetPointerType statismo::StatisticalModel< T >::DrawPCABasisSample | ( | unsigned | componentNumber | ) | const |
Draws the sample corresponding to the ith pca matrix
componentNumber | The number of the PCA Basis to be retrieved |
StatisticalModel< T >::DatasetPointerType statismo::StatisticalModel< T >::DrawSample | ( | const VectorType & | coefficients, |
bool | addNoise = false |
||
) | const |
Draws the sample with the given coefficients
coefficients | A coefficient vector. The size of the coefficient vector should be smaller than number of factors in the model. Otherwise an exception is thrown. |
addNoise | If true, the Gaussian noise assumed in the model is added to the sample |
StatisticalModel< T >::DatasetPointerType statismo::StatisticalModel< T >::DrawSample | ( | bool | addNoise = false | ) | const |
As StatisticalModel::DrawSample, but where the coefficients are chosen at random according to a standard normal distribution
addNoise | If true, the Gaussian noise assumed in the model is added to the sample |
StatisticalModel< T >::ValueType statismo::StatisticalModel< T >::DrawSampleAtPoint | ( | const VectorType & | coefficients, |
const PointType & | point, | ||
bool | addNoise = false |
||
) | const |
Returns the value of the sample defined by coefficients at the specified point. This method computes the value of the sample only for the given point, and is thus much more efficient that calling DrawSample, if only a few points are of interest.
coefficients | the coefficients of the sample |
the | point of the sample where it is evaluated |
addNoise | If true, the Gaussian noise assumed in the model is added to the sample |
StatisticalModel< T >::ValueType statismo::StatisticalModel< T >::DrawSampleAtPoint | ( | const VectorType & | coefficients, |
unsigned | pointId, | ||
bool | addNoise = false |
||
) | const |
Returns the value of the sample defined by coefficients at the specified pointID. This method computes the value of the sample only for the given point, and is thus much more efficient that calling DrawSample, if only a few points are of interest.
coefficients | the coefficients of the sample |
the | point of the sample where it is evaluated |
addNoise | If true, the Gaussian noise assumed in the model is added to the sample |
VectorType statismo::StatisticalModel< T >::DrawSampleVector | ( | const VectorType & | coefficients, |
bool | addNoise = false |
||
) | const |
Returns an instance for the given coefficients as a vector.
addNoise | If true, the Gaussian noise assumed in the model is added to the sample |
StatisticalModel< T >::ValueType statismo::StatisticalModel< T >::EvaluateSampleAtPoint | ( | DatasetConstPointerType | sample, |
unsigned | ptId | ||
) | const |
Returns the value of the given sample at the point specified with the ptId
sample | A sample |
ptId | the point id where to evaluate the sample |
StatisticalModel< T >::ValueType statismo::StatisticalModel< T >::EvaluateSampleAtPoint | ( | DatasetConstPointerType | sample, |
const PointType & | pt | ||
) | const |
Returns the value of the given sample corresponding to the given domain point
sample | A sample |
point | the (domain) point on which the sample should be evaluated. |
MatrixType statismo::StatisticalModel< T >::GetCovarianceAtPoint | ( | const PointType & | pt1, |
const PointType & | pt2 | ||
) | const |
Returns the variance in the model for point pt
pt | The point |
MatrixType statismo::StatisticalModel< T >::GetCovarianceAtPoint | ( | unsigned | ptId1, |
unsigned | ptId2 | ||
) | const |
Returns the variance in the model for point pt
pt | The point |
MatrixType statismo::StatisticalModel< T >::GetCovarianceMatrix | ( | ) | const |
Returns the covariance matrix for the model. If the model is defined on n points, in d dimensions, then this is a matrix of n
block matrices corresponding to the covariance at each point.
|
inline |
Return the domain of the statistical model
MatrixType statismo::StatisticalModel< T >::GetJacobian | ( | const PointType & | pt | ) | const |
Computes the jacobian of the Statistical model at a given point
pt | The point where the Jacobian is computed |
jacobian | Output parameter where the jacobian is stored. |
const VectorType & statismo::StatisticalModel< T >::GetMeanVector | ( | ) | const |
Returns a vector holding the mean. Assume the mean is defined on
points, the returned mean vector
has dimensionality
, i.e. the
components are stacked into the vector. The order of the components in the vector is undefined and depends on the representer.
const ModelInfo & statismo::StatisticalModel< T >::GetModelInfo | ( | ) | const |
float statismo::StatisticalModel< T >::GetNoiseVariance | ( | ) | const |
Returns the variance of the noise of the error term, that was set when the model was built.
unsigned int statismo::StatisticalModel< T >::GetNumberOfPrincipalComponents | ( | ) | const |
MatrixType statismo::StatisticalModel< T >::GetOrthonormalPCABasisMatrix | ( | ) | const |
Returns the PCA Matrix, but with its principal axis normalized to unit length.
const MatrixType & statismo::StatisticalModel< T >::GetPCABasisMatrix | ( | ) | const |
Returns a matrix with the PCA Basis as its columns. Assume the shapes are defined on
points, the returned matrix
has dimensionality
, i.e. the
components are stacked into the matrix. The order of the components in the matrix is undefined and depends on the representer.
const VectorType & statismo::StatisticalModel< T >::GetPCAVarianceVector | ( | ) | const |
Returns a vector where each element holds the variance of the corresponding principal component in data space
|
inline |
Return an instance of the representer
|
static |
Returns a new statistical model, which is loaded from the given HDF5 file
filename | The filename |
maxNumberOfPCAComponents | The maximal number of pca components that are loaded to create the model. |
|
static |
Returns a new statistical model, which is stored in the given HDF5 Group
modelroot | A h5 group where the model is saved |
maxNumberOfPCAComponents | The maximal number of pca components that are loaded to create the model. |
VectorType statismo::StatisticalModel< T >::RobustlyComputeCoefficientsForDataset | ( | DatasetConstPointerType | dataset, |
unsigned | nIterations = 100 , |
||
unsigned | nu = 6 , |
||
double | sigma2 = 1 |
||
) | const |
Computes the coefficients of the latent variables in a robust way. Instead of assuming Normally distributed noise on the data set points such as it is implicitely assumed in the ComputeCoefficientsForDataset, A student-t distribution is assumed for the noise. The solution is obtained using an EM algorithm.
dataset | The dataset |
nIterations | The number of iterations for the EM algorithm |
nu | The number of degrees of Freedom for the Student-t distribution defining the noise model |
sigma2 | The scale parameter of the t-distribution |
void statismo::StatisticalModel< T >::Save | ( | const std::string & | filename | ) | const |
Saves the statistical model to a HDF5 file
filename | The filename (preferred extension is .h5) |
void statismo::StatisticalModel< T >::Save | ( | const H5::Group & | modelRoot | ) | const |
Saves the statistical model to the given HDF5 group.
modelRoot | the group where to store the model |
void statismo::StatisticalModel< T >::SetModelInfo | ( | const ModelInfo & | modelInfo | ) |
Sets the model information. This is for library internal use only.