queso-0.56.1
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 actualValue (const V &domainVector, const V *domainDirection, V *gradVector, M *hessianMatrix, V *hessianEffect) const
 Actual value 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. 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::BaseGaussianLikelihood< V, M >
virtual void evaluateModel (const V &domainVector, const V *domainDirection, V &modelOutput, V *gradVector, M *hessianMatrix, V *hessianEffect) const =0
 Evaluates the user's model at the point domainVector. More...
 
 BaseGaussianLikelihood (const char *prefix, const VectorSet< V, M > &domainSet, const V &observations)
 Default constructor. More...
 
virtual ~BaseGaussianLikelihood ()
 Destructor. More...
 
- Public Member Functions inherited from QUESO::BaseScalarFunction< V, M >
 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...
 

Private Attributes

const GslBlockMatrixm_covariance
 

Additional Inherited Members

- Protected Attributes inherited from QUESO::BaseGaussianLikelihood< 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(), QUESO::GslMatrix::numRowsLocal(), and queso_error_msg.

36  : BaseGaussianLikelihood<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 }
GslMatrix & getBlock(unsigned int i) const
Return block i in the block diagonal matrix.
unsigned int numBlocks() const
Return the number of blocks in the block diagonal matrix.
unsigned int numRowsLocal() const
Returns the local row dimension of this matrix.
Definition: GslMatrix.C:275
#define queso_error_msg(msg)
Definition: asserts.h:47
const VectorSet< V, M > & domainSet() const
Access to the protected attribute m_domainSet: domain set of the scalar function. ...

Destructor.

Definition at line 51 of file GaussianLikelihoodBlockDiagonalCovarianceRandomCoefficients.C.

52 {
53 }

Member Function Documentation

template<class V , class M >
double QUESO::GaussianLikelihoodBlockDiagonalCovarianceRandomCoefficients< V, M >::actualValue ( const V &  domainVector,
const V *  domainDirection,
V *  gradVector,
M *  hessianMatrix,
V *  hessianEffect 
) const
virtual

Actual value of the scalar function.

Implements QUESO::BaseScalarFunction< V, M >.

Definition at line 57 of file GaussianLikelihoodBlockDiagonalCovarianceRandomCoefficients.C.

60 {
61  return std::exp(this->lnValue(domainVector, domainDirection, gradVector,
62  hessianMatrix, hessianEffect));
63 }
virtual double lnValue(const V &domainVector, const V *domainDirection, V *gradVector, M *hessianMatrix, V *hessianEffect) const
Logarithm of the value of the scalar function.
template<class V , class M >
double QUESO::GaussianLikelihoodBlockDiagonalCovarianceRandomCoefficients< V, M >::lnValue ( const V &  domainVector,
const V *  domainDirection,
V *  gradVector,
M *  hessianMatrix,
V *  hessianEffect 
) 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.

Implements QUESO::BaseScalarFunction< V, M >.

Definition at line 67 of file GaussianLikelihoodBlockDiagonalCovarianceRandomCoefficients.C.

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

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 Thu Dec 15 2016 13:23:14 for queso-0.56.1 by  doxygen 1.8.5