queso-0.53.0
Public Member Functions | Protected Attributes | List of all members
QUESO::BaseGaussianLikelihood< V, M > Class Template Referenceabstract

Base class for canned Gaussian likelihoods. More...

#include <GaussianLikelihood.h>

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

Public Member Functions

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...
 
Constructor/Destructor methods.
 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...
 
virtual double actualValue (const V &domainVector, const V *domainDirection, V *gradVector, M *hessianMatrix, V *hessianEffect) const =0
 Actual value of the scalar function. More...
 
virtual double lnValue (const V &domainVector, const V *domainDirection, V *gradVector, M *hessianMatrix, V *hessianEffect) const =0
 Logarithm of the value of the scalar function. More...
 

Protected Attributes

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::BaseGaussianLikelihood< V, M >

Base class for canned Gaussian likelihoods.

This class is an abstract base class for 'canned' Gaussian likelihoods. All this class does is add a pure virtual function called evaluateModel that the user will implement to interact with the forward code.

Definition at line 48 of file GaussianLikelihood.h.

Constructor & Destructor Documentation

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

Default constructor.

The vector of observations must be passed. This will be used when evaluating the likelihood functional

Definition at line 33 of file GaussianLikelihood.C.

36  : BaseScalarFunction<V, M>(prefix, domainSet),
37  m_observations(observations)
38 {
39 }
const VectorSet< V, M > & domainSet() const
Access to the protected attribute m_domainSet: domain set of the scalar function. ...
template<class V , class M >
QUESO::BaseGaussianLikelihood< V, M >::~BaseGaussianLikelihood ( )
virtual

Destructor.

Definition at line 42 of file GaussianLikelihood.C.

43 {
44 }

Member Function Documentation

template<class V = GslVector, class M = GslMatrix>
virtual void QUESO::BaseGaussianLikelihood< V, M >::evaluateModel ( const V &  domainVector,
const V *  domainDirection,
V &  modelOutput,
V *  gradVector,
M *  hessianMatrix,
V *  hessianEffect 
) const
pure virtual

Evaluates the user's model at the point domainVector.

This is pure virtual, so the user must implement this when subclassing a Gaussian likelihood class. Note that void is returned. The user will fill up the modelOutput vector with output from the model. This represents a vector of synthetic observations that will be to compare to actual observations when computing the likelihood functional.

The first n components of domainVector are the model parameters. The rest of domainVector contains the hyperparameters, if any. For example, in GaussianLikelihoodFullCovarainceRandomCoefficient, the last component of domainVector contains the multiplicative coefficient of the observational covariance matrix. In this case, the user need not concern themselves with this parameter as it is handled not in the model evaluation but by the likelihood evaluation.

Member Data Documentation

template<class V = GslVector, class M = GslMatrix>
const V& QUESO::BaseGaussianLikelihood< V, M >::m_observations
protected

Definition at line 86 of file GaussianLikelihood.h.


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

Generated on Thu Jun 11 2015 13:52:35 for queso-0.53.0 by  doxygen 1.8.5