queso-0.57.0
Public Member Functions | Private Attributes | List of all members
QUESO::GaussianLikelihoodBlockDiagonalCovarianceRandomCoefficients< V, M > Class Template Reference

A class representing a Gaussian likelihood with block-diagonal covariance matrix. More...

#include <GaussianLikelihoodBlockDiagonalCovarianceRandomCoefficients.h>

Inheritance diagram for QUESO::GaussianLikelihoodBlockDiagonalCovarianceRandomCoefficients< V, M >:
Inheritance graph
[legend]
Collaboration diagram for QUESO::GaussianLikelihoodBlockDiagonalCovarianceRandomCoefficients< V, M >:
Collaboration graph
[legend]

Public Member Functions

virtual double lnValue (const V &domainVector) const
 Logarithm of the value of the scalar function. More...
 
Constructor/Destructor methods.
 GaussianLikelihoodBlockDiagonalCovarianceRandomCoefficients (const char *prefix, const VectorSet< V, M > &domainSet, const V &observations, const GslBlockMatrix &covariance)
 Default constructor. More...
 
virtual ~GaussianLikelihoodBlockDiagonalCovarianceRandomCoefficients ()
 Destructor. More...
 
- Public Member Functions inherited from QUESO::LikelihoodBase< V, M >
virtual void evaluateModel (const V &domainVector, const V *domainDirection, V &modelOutput, V *gradVector, M *hessianMatrix, V *hessianEffect) const
 Deprecated. Evaluates the user's model at the point domainVector. More...
 
virtual void evaluateModel (const V &domainVector, V &modelOutput) const
 Evaluates the user's model at the point domainVector. More...
 
virtual double actualValue (const V &domainVector, const V *, V *, M *, V *) const
 Actual value of the scalar function. More...
 
 LikelihoodBase (const char *prefix, const VectorSet< V, M > &domainSet, const V &observations)
 Default constructor. More...
 
virtual ~LikelihoodBase ()=0
 Destructor, pure to make this class abstract. More...
 
- Public Member Functions inherited from QUESO::BaseScalarFunction< V, M >
void setFiniteDifferenceStepSize (double fdStepSize)
 Sets the step size for finite differencing gradients. More...
 
void setFiniteDifferenceStepSize (unsigned int i, double fdStepSize)
 
 BaseScalarFunction (const char *prefix, const VectorSet< V, M > &domainSet)
 Default constructor. More...
 
virtual ~BaseScalarFunction ()
 Destructor. More...
 
const VectorSet< V, M > & domainSet () const
 Access to the protected attribute m_domainSet: domain set of the scalar function. More...
 
virtual double lnValue (const V &domainVector, const V *domainDirection, V *gradVector, M *hessianMatrix, V *hessianEffect) const
 Logarithm of the value of the scalar function. Deprecated. More...
 
virtual double lnValue (const V &domainVector, V &gradVector) const
 Returns the logarithm of the function and its gradient at domainVector. More...
 
virtual double lnValue (const V &domainVector, V &gradVector, const V &domainDirection, V &hessianEffect) const
 

Private Attributes

const GslBlockMatrixm_covariance
 

Additional Inherited Members

- Protected Attributes inherited from QUESO::LikelihoodBase< V, M >
const V & m_observations
 
- Protected Attributes inherited from QUESO::BaseScalarFunction< V, M >
const BaseEnvironmentm_env
 
std::string m_prefix
 
const VectorSet< V, M > & m_domainSet
 Domain set of the scalar function. More...
 

Detailed Description

template<class V = GslVector, class M = GslMatrix>
class QUESO::GaussianLikelihoodBlockDiagonalCovarianceRandomCoefficients< V, M >

A class representing a Gaussian likelihood with block-diagonal covariance matrix.

Each block diagonal matrix has a multiplicative coefficient that is treated as a hyperparameter to be inferred during the sampling procedure.

Definition at line 47 of file GaussianLikelihoodBlockDiagonalCovarianceRandomCoefficients.h.

Constructor & Destructor Documentation

template<class V , class M >
QUESO::GaussianLikelihoodBlockDiagonalCovarianceRandomCoefficients< V, M >::GaussianLikelihoodBlockDiagonalCovarianceRandomCoefficients ( const char *  prefix,
const VectorSet< V, M > &  domainSet,
const V &  observations,
const GslBlockMatrix covariance 
)

Default constructor.

Instantiates a Gaussian likelihood function, given a prefix, its domain, a vector of observations and a block diagonal covariance matrix. The diagonal covariance matrix is of type GslBlockMatrix. Each block in the block diagonal matrix is an object of type GslMatrix.

Each block diagonal matrix has a multiplicative coefficient that is treated as a hyperparameter to be inferred during the sampling procedure.

Definition at line 33 of file GaussianLikelihoodBlockDiagonalCovarianceRandomCoefficients.C.

References QUESO::GslBlockMatrix::getBlock(), QUESO::GaussianLikelihoodBlockDiagonalCovarianceRandomCoefficients< V, M >::m_covariance, QUESO::GslBlockMatrix::numBlocks(), and QUESO::GslMatrix::numRowsLocal().

36  : LikelihoodBase<V, M>(prefix, domainSet, observations),
37  m_covariance(covariance)
38 {
39  unsigned int totalDim = 0;
40 
41  for (unsigned int i = 0; i < this->m_covariance.numBlocks(); i++) {
42  totalDim += this->m_covariance.getBlock(i).numRowsLocal();
43  }
44 
45  if (totalDim != observations.sizeLocal()) {
46  queso_error_msg("Covariance matrix not same dimension as observation vector");
47  }
48 }
unsigned int numBlocks() const
Return the number of blocks in the block diagonal matrix.
const VectorSet< V, M > & domainSet() const
Access to the protected attribute m_domainSet: domain set of the scalar function. ...
GslMatrix & getBlock(unsigned int i) const
Return block i in the block diagonal matrix.
unsigned int numRowsLocal() const
Returns the local row dimension of this matrix.
Definition: GslMatrix.C:275

Destructor.

Definition at line 51 of file GaussianLikelihoodBlockDiagonalCovarianceRandomCoefficients.C.

52 {
53 }

Member Function Documentation

template<class V , class M >
double QUESO::GaussianLikelihoodBlockDiagonalCovarianceRandomCoefficients< V, M >::lnValue ( const V &  domainVector) const
virtual

Logarithm of the value of the scalar function.

The last n elements of domainVector, where n is the number of blocks in the block diagonal covariance matrix, are treated as hyperparameters and will be sample in a statistical inversion.

The user need not concern themselves with handling these in the model evaluation routine, since they are handled by the likelihood evaluation routine.

Reimplemented from QUESO::BaseScalarFunction< V, M >.

Definition at line 57 of file GaussianLikelihoodBlockDiagonalCovarianceRandomCoefficients.C.

58 {
59  V modelOutput(this->m_observations, 0, 0); // At least it's not a copy
60  V weightedMisfit(this->m_observations, 0, 0); // At least it's not a copy
61 
62  this->evaluateModel(domainVector, modelOutput);
63 
64  // Compute misfit G(x) - y
65  modelOutput -= this->m_observations;
66 
67  // Solve \Sigma u = G(x) - y for u
68  this->m_covariance.invertMultiply(modelOutput, weightedMisfit);
69 
70  // Deal with the multiplicative coefficients for each of the blocks
71  unsigned int numBlocks = this->m_covariance.numBlocks();
72  unsigned int offset = 0;
73 
74  // For each block...
75  double cov_norm_factor = 0.0;
76  for (unsigned int i = 0; i < this->m_covariance.numBlocks(); i++) {
77 
78  // ...find the right hyperparameter
79  unsigned int index = domainVector.sizeLocal() + (i - numBlocks);
80  double coefficient = domainVector[index];
81 
82  // ...divide the appropriate parts of the solution by the coefficient
83  unsigned int blockDim = this->m_covariance.getBlock(i).numRowsLocal();
84  for (unsigned int j = 0; j < blockDim; j++) {
85  // 'coefficient' is a variance, so we divide by it
86  modelOutput[offset+j] /= coefficient;
87  }
88 
89  // Keep track of the part of the covariance matrix that appears in the
90  // normalising constant because of the hyperparameter
91  double cov_determinant = this->m_covariance.getBlock(i).determinant();
92  cov_determinant = std::sqrt(cov_determinant);
93 
94  coefficient = std::sqrt(coefficient);
95  cov_norm_factor += std::log(std::pow(coefficient, blockDim) * cov_determinant);
96 
97  offset += blockDim;
98  }
99 
100  // Compute (G(x) - y)^T \Sigma^{-1} (G(x) - y)
101  modelOutput *= weightedMisfit;
102 
103  double norm2_squared = modelOutput.sumOfComponents(); // This is square of 2-norm
104 
105  return -0.5 * norm2_squared - cov_norm_factor;
106 }
void invertMultiply(const GslVector &b, GslVector &x) const
This function calculates the inverse of this matrix, multiplies it with vector b and stores the resul...
unsigned int numBlocks() const
Return the number of blocks in the block diagonal matrix.
double determinant() const
Calculates the determinant of this matrix.
Definition: GslMatrix.C:1038
GslMatrix & getBlock(unsigned int i) const
Return block i in the block diagonal matrix.
virtual void evaluateModel(const V &domainVector, const V *domainDirection, V &modelOutput, V *gradVector, M *hessianMatrix, V *hessianEffect) const
Deprecated. Evaluates the user&#39;s model at the point domainVector.
unsigned int numRowsLocal() const
Returns the local row dimension of this matrix.
Definition: GslMatrix.C:275

Member Data Documentation

template<class V = GslVector, class M = GslMatrix>
const GslBlockMatrix& QUESO::GaussianLikelihoodBlockDiagonalCovarianceRandomCoefficients< V, M >::m_covariance
private

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

Generated on Sat Apr 22 2017 14:04:40 for queso-0.57.0 by  doxygen 1.8.5