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

A class for handling Gaussian joint PDFs. More...

#include <GaussianJointPdf.h>

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

Public Member Functions

virtual void print (std::ostream &os) const
 Print method for informational and logging purposes. More...
 
Constructor/Destructor methods
 GaussianJointPdf (const char *prefix, const VectorSet< V, M > &domainSet, const V &lawExpVector, const V &lawVarVector)
 Constructor. More...
 
 GaussianJointPdf (const char *prefix, const VectorSet< V, M > &domainSet, const V &lawExpVector, const M &lawCovMatrix)
 Constructor. More...
 
 ~GaussianJointPdf ()
 Destructor. More...
 
Math methods
double actualValue (const V &domainVector, const V *domainDirection, V *gradVector, M *hessianMatrix, V *hessianEffect) const
 Actual value of the Gaussian PDF: More...
 
double lnValue (const V &domainVector, const V *domainDirection, V *gradVector, M *hessianMatrix, V *hessianEffect) const
 Logarithm of the value of the Gaussian PDF (scalar function). More...
 
virtual void distributionMean (V &meanVector) const
 Mean value of the underlying random variable. More...
 
virtual void distributionVariance (M &covMatrix) const
 Covariance matrix of the underlying random variable. More...
 
double computeLogOfNormalizationFactor (unsigned int numSamples, bool updateFactorInternally) const
 Computes the logarithm of the normalization factor. More...
 
void updateLawExpVector (const V &newLawExpVector)
 Updates the mean with the new value newLawExpVector. More...
 
void updateLawCovMatrix (const M &newLawCovMatrix)
 Updates the lower triangular matrix from Cholesky decomposition of the covariance matrix to the new value newLowerCholLawCovMatrix. More...
 
const M & lawCovMatrix () const
 Returns the covariance matrix; access to protected attribute m_lawCovMatrix. More...
 
const V & lawExpVector () const
 Access to the vector of mean values and private attribute: m_lawExpVector. More...
 
const V & lawVarVector () const
 Access to the vector of variance values and private attribute: m_lawVarVector. More...
 
- Public Member Functions inherited from QUESO::BaseJointPdf< V, M >
 BaseJointPdf (const char *prefix, const VectorSet< V, M > &domainSet)
 Default constructor. More...
 
virtual ~BaseJointPdf ()
 Destructor. More...
 
virtual void setNormalizationStyle (unsigned int value) const
 
void setLogOfNormalizationFactor (double value) const
 Sets a logarithmic value to be used in the normalization factor (stored in the protected attribute m_normalizationStyle.) 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
 Returns the logarithm of the function at domainVector. 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
 

Protected Attributes

V * m_lawExpVector
 
V * m_lawVarVector
 
bool m_diagonalCovMatrix
 
const M * m_lawCovMatrix
 
- Protected Attributes inherited from QUESO::BaseJointPdf< V, M >
unsigned int m_normalizationStyle
 
double m_logOfNormalizationFactor
 
- 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...
 

Additional Inherited Members

- Protected Member Functions inherited from QUESO::BaseJointPdf< V, M >
double commonComputeLogOfNormalizationFactor (unsigned int numSamples, bool updateFactorInternally) const
 Common method (to the derived classes) to compute the logarithm of the normalization factor. More...
 

Detailed Description

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

A class for handling Gaussian joint PDFs.

This class allows the mathematical definition of a Gaussian Joint PDF.

Definition at line 50 of file GaussianJointPdf.h.

Constructor & Destructor Documentation

template<class V , class M >
QUESO::GaussianJointPdf< V, M >::GaussianJointPdf ( const char *  prefix,
const VectorSet< V, M > &  domainSet,
const V &  lawExpVector,
const V &  lawVarVector 
)

Constructor.

Constructs a new object, given a prefix and the domain of the PDF, a vector of mean values, lawExpVector, and a vector of covariance values lawVarVector (an alternative representation for a diagonal covariance matrix).

Definition at line 33 of file GaussianJointPdf.C.

References QUESO::BaseEnvironment::displayVerbosity(), QUESO::GaussianJointPdf< V, M >::lawExpVector(), QUESO::GaussianJointPdf< V, M >::lawVarVector(), QUESO::BaseScalarFunction< V, M >::m_env, QUESO::BaseScalarFunction< V, M >::m_prefix, and QUESO::BaseEnvironment::subDisplayFile().

38  :
39  BaseJointPdf<V,M>(((std::string)(prefix)+"gau").c_str(),domainSet),
42  m_diagonalCovMatrix(true),
43  m_lawCovMatrix (m_domainSet.vectorSpace().newDiagMatrix(lawVarVector))
44 {
45 
46  if ((m_env.subDisplayFile()) && (m_env.displayVerbosity() >= 54)) {
47  *m_env.subDisplayFile() << "Entering GaussianJointPdf<V,M>::constructor() [1]"
48  << ": prefix = " << m_prefix
49  << std::endl;
50  }
51 
52  if ((m_env.subDisplayFile()) && (m_env.displayVerbosity() >= 55)) {
53  *m_env.subDisplayFile() << "In GaussianJointPdf<V,M>::constructor()"
54  //<< ", prefix = " << m_prefix
55  << ": meanVector = " << this->lawExpVector()
56  << ", Variances = " << this->lawVarVector()
57  << std::endl;
58  }
59 
60  if ((m_env.subDisplayFile()) && (m_env.displayVerbosity() >= 54)) {
61  *m_env.subDisplayFile() << "Leaving GaussianJointPdf<V,M>::constructor() [1]"
62  << ": prefix = " << m_prefix
63  << std::endl;
64  }
65 }
const V & lawExpVector() const
Access to the vector of mean values and private attribute: m_lawExpVector.
const VectorSet< V, M > & m_domainSet
Domain set of the scalar function.
const BaseEnvironment & m_env
const VectorSet< V, M > & domainSet() const
Access to the protected attribute m_domainSet: domain set of the scalar function. ...
const V & lawVarVector() const
Access to the vector of variance values and private attribute: m_lawVarVector.
std::ofstream * subDisplayFile() const
Access function for m_subDisplayFile (displays file on stream).
Definition: Environment.C:320
unsigned int displayVerbosity() const
Definition: Environment.C:450
template<class V , class M >
QUESO::GaussianJointPdf< V, M >::GaussianJointPdf ( const char *  prefix,
const VectorSet< V, M > &  domainSet,
const V &  lawExpVector,
const M &  lawCovMatrix 
)

Constructor.

Constructs a new object, given a prefix and the image set of the vector realizer, a vector of mean values, lawExpVector, and a covariance matrix, lawCovMatrix.

Definition at line 68 of file GaussianJointPdf.C.

References QUESO::BaseEnvironment::displayVerbosity(), QUESO::GaussianJointPdf< V, M >::lawExpVector(), QUESO::BaseScalarFunction< V, M >::m_env, QUESO::BaseScalarFunction< V, M >::m_prefix, and QUESO::BaseEnvironment::subDisplayFile().

73  :
74  BaseJointPdf<V,M>(((std::string)(prefix)+"gau").c_str(),domainSet),
76  m_lawVarVector (domainSet.vectorSpace().newVector(INFINITY)), // FIX ME
77  m_diagonalCovMatrix(false),
79 {
80  if ((m_env.subDisplayFile()) && (m_env.displayVerbosity() >= 54)) {
81  *m_env.subDisplayFile() << "Entering GaussianJointPdf<V,M>::constructor() [2]"
82  << ": prefix = " << m_prefix
83  << std::endl;
84  }
85 
86  if ((m_env.subDisplayFile()) && (m_env.displayVerbosity() >= 55)) {
87  *m_env.subDisplayFile() << "In GaussianJointPdf<V,M>::constructor()"
88  //<< ", prefix = " << m_prefix
89  << ": meanVector = " << this->lawExpVector()
90  << ", Covariance Matrix = " << lawCovMatrix
91  << std::endl;
92  }
93 
94  if ((m_env.subDisplayFile()) && (m_env.displayVerbosity() >= 54)) {
95  *m_env.subDisplayFile() << "Leaving GaussianJointPdf<V,M>::constructor() [2]"
96  << ": prefix = " << m_prefix
97  << std::endl;
98  }
99 }
const M & lawCovMatrix() const
Returns the covariance matrix; access to protected attribute m_lawCovMatrix.
const V & lawExpVector() const
Access to the vector of mean values and private attribute: m_lawExpVector.
const BaseEnvironment & m_env
const VectorSet< V, M > & domainSet() const
Access to the protected attribute m_domainSet: domain set of the scalar function. ...
std::ofstream * subDisplayFile() const
Access function for m_subDisplayFile (displays file on stream).
Definition: Environment.C:320
unsigned int displayVerbosity() const
Definition: Environment.C:450
template<class V , class M >
QUESO::GaussianJointPdf< V, M >::~GaussianJointPdf ( )

Destructor.

Definition at line 102 of file GaussianJointPdf.C.

103 {
104  delete m_lawCovMatrix;
105  delete m_lawVarVector;
106  delete m_lawExpVector;
107 }

Member Function Documentation

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

Actual value of the Gaussian PDF:

This method calls lnValue() and applies the exponential to it.

Implements QUESO::BaseJointPdf< V, M >.

Definition at line 152 of file GaussianJointPdf.C.

References QUESO::queso_require_equal_to_msg.

158 {
159  if ((m_env.subDisplayFile()) && (m_env.displayVerbosity() >= 55)) {
160  *m_env.subDisplayFile() << "Entering GaussianJointPdf<V,M>::actualValue()"
161  << ", meanVector = " << *m_lawExpVector
162  << ", lawCovMatrix = " << *m_lawCovMatrix
163  << ": domainVector = " << domainVector
164  << std::endl;
165  }
166 
167  queso_require_equal_to_msg(domainVector.sizeLocal(), this->m_domainSet.vectorSpace().dimLocal(), "invalid input");
168 
169  queso_require_msg(!(hessianMatrix || hessianEffect), "incomplete code for gradVector, hessianMatrix and hessianEffect calculations");
170 
171  double returnValue = 0.;
172 
173  if (this->m_domainSet.contains(domainVector) == false) {
174  // What should the gradient be here?
175  returnValue = 0.;
176  }
177  else {
178  // Already normalised (so the gradient will have the normalisation constant
179  // in it)
180  returnValue = std::exp(this->lnValue(domainVector,domainDirection,gradVector,hessianMatrix,hessianEffect));
181 
182  if (gradVector) {
183  (*gradVector) *= returnValue;
184  }
185  }
186 
187  if ((m_env.subDisplayFile()) && (m_env.displayVerbosity() >= 55)) {
188  *m_env.subDisplayFile() << "Leaving GaussianJointPdf<V,M>::actualValue()"
189  << ", meanVector = " << *m_lawExpVector
190  << ", lawCovMatrix = " << *m_lawCovMatrix
191  << ": domainVector = " << domainVector
192  << ", returnValue = " << returnValue
193  << std::endl;
194  }
195 
196  return returnValue;
197 }
double lnValue(const V &domainVector, const V *domainDirection, V *gradVector, M *hessianMatrix, V *hessianEffect) const
Logarithm of the value of the Gaussian PDF (scalar function).
const VectorSet< V, M > & m_domainSet
Domain set of the scalar function.
const BaseEnvironment & m_env
MonteCarloSGOptions::MonteCarloSGOptions(const BaseEnvironment &env, const char *prefix, const McOptionsValues &alternativeOptionsValues queso_require_equal_to_msg)(m_env.optionsInputFileName(), std::string(""), std::string("this constructor is incompatible with the existence of an options input file"))
std::ofstream * subDisplayFile() const
Access function for m_subDisplayFile (displays file on stream).
Definition: Environment.C:320
unsigned int displayVerbosity() const
Definition: Environment.C:450
template<class V , class M >
double QUESO::GaussianJointPdf< V, M >::computeLogOfNormalizationFactor ( unsigned int  numSamples,
bool  updateFactorInternally 
) const
virtual

Computes the logarithm of the normalization factor.

This routine calls BaseJointPdf::commonComputeLogOfNormalizationFactor().

Implements QUESO::BaseJointPdf< V, M >.

Definition at line 329 of file GaussianJointPdf.C.

References QUESO::BaseJointPdf< V, M >::commonComputeLogOfNormalizationFactor().

330 {
331  double value = 0.;
332 
333  if ((m_env.subDisplayFile()) && (m_env.displayVerbosity() >= 2)) {
334  *m_env.subDisplayFile() << "Entering GaussianJointPdf<V,M>::computeLogOfNormalizationFactor()"
335  << std::endl;
336  }
337  value = BaseJointPdf<V,M>::commonComputeLogOfNormalizationFactor(numSamples, updateFactorInternally);
338  if ((m_env.subDisplayFile()) && (m_env.displayVerbosity() >= 2)) {
339  *m_env.subDisplayFile() << "Leaving GaussianJointPdf<V,M>::computeLogOfNormalizationFactor()"
340  << ", m_logOfNormalizationFactor = " << m_logOfNormalizationFactor
341  << std::endl;
342  }
343 
344  return value;
345 }
const BaseEnvironment & m_env
double commonComputeLogOfNormalizationFactor(unsigned int numSamples, bool updateFactorInternally) const
Common method (to the derived classes) to compute the logarithm of the normalization factor...
Definition: JointPdf.C:115
double m_logOfNormalizationFactor
Definition: JointPdf.h:122
std::ofstream * subDisplayFile() const
Access function for m_subDisplayFile (displays file on stream).
Definition: Environment.C:320
unsigned int displayVerbosity() const
Definition: Environment.C:450
template<class V , class M >
void QUESO::GaussianJointPdf< V, M >::distributionMean ( V &  meanVector) const
virtual

Mean value of the underlying random variable.

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

Definition at line 301 of file GaussianJointPdf.C.

302 {
303  meanVector = this->lawExpVector();
304 }
const V & lawExpVector() const
Access to the vector of mean values and private attribute: m_lawExpVector.
template<class V , class M >
void QUESO::GaussianJointPdf< V, M >::distributionVariance ( M &  covMatrix) const
virtual

Covariance matrix of the underlying random variable.

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

Definition at line 308 of file GaussianJointPdf.C.

309 {
310  queso_assert_equal_to (covMatrix.numCols(), covMatrix.numRowsGlobal());
311 
312  if (m_diagonalCovMatrix) {
313  covMatrix.zeroLower();
314  covMatrix.zeroUpper();
315 
316  unsigned int n_comp = this->lawVarVector().sizeLocal();
317  queso_assert_equal_to (n_comp, covMatrix.numCols());
318 
319  for (unsigned int i = 0; i < n_comp; ++i) {
320  covMatrix(i,i) = this->lawVarVector()[i];
321  }
322  } else {
323  covMatrix = *this->m_lawCovMatrix;
324  }
325 }
const V & lawVarVector() const
Access to the vector of variance values and private attribute: m_lawVarVector.
template<class V , class M >
const M & QUESO::GaussianJointPdf< V, M >::lawCovMatrix ( ) const

Returns the covariance matrix; access to protected attribute m_lawCovMatrix.

Definition at line 369 of file GaussianJointPdf.C.

Referenced by QUESO::ScaledCovMatrixTKGroup< V, M >::setPreComputingPosition(), and QUESO::MetropolisAdjustedLangevinTK< V, M >::setPreComputingPosition().

370 {
371  return *m_lawCovMatrix;
372 }
template<class V , class M >
const V & QUESO::GaussianJointPdf< V, M >::lawExpVector ( ) const

Access to the vector of mean values and private attribute: m_lawExpVector.

Definition at line 111 of file GaussianJointPdf.C.

Referenced by QUESO::GaussianJointPdf< V, M >::GaussianJointPdf().

112 {
113  return *m_lawExpVector;
114 }
template<class V , class M >
const V & QUESO::GaussianJointPdf< V, M >::lawVarVector ( ) const

Access to the vector of variance values and private attribute: m_lawVarVector.

Definition at line 118 of file GaussianJointPdf.C.

Referenced by QUESO::GaussianJointPdf< V, M >::GaussianJointPdf().

119 {
120  return *m_lawVarVector;
121 }
template<class V , class M >
double QUESO::GaussianJointPdf< V, M >::lnValue ( const V &  domainVector,
const V *  domainDirection,
V *  gradVector,
M *  hessianMatrix,
V *  hessianEffect 
) const
virtual

Logarithm of the value of the Gaussian PDF (scalar function).

The ln(value) comes from a summation of the Gaussian density:

\[ lnValue =- \sum_i \frac{1}{\sqrt{|covMatrix|} \sqrt{2 \pi}} exp(-\frac{(domainVector_i - lawExpVector_i)* covMatrix^{-1}* (domainVector_i - lawExpVector_i) }{2}, \]

where the $ covMatrix $ may recovered via this->lawVarVector(), in case of diagonal matrices or via this->m_lawCovMatrix, otherwise.

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

Definition at line 201 of file GaussianJointPdf.C.

207 {
208  if ((m_env.subDisplayFile()) && (m_env.displayVerbosity() >= 55)) {
209  *m_env.subDisplayFile() << "Entering GaussianJointPdf<V,M>::lnValue()"
210  << ", meanVector = " << *m_lawExpVector
211  << ", lawCovMatrix = " << *m_lawCovMatrix
212  << ": domainVector = " << domainVector
213  << std::endl;
214  }
215 
216  queso_require_msg(!(domainDirection || hessianMatrix || hessianEffect), "incomplete code for gradVector, hessianMatrix and hessianEffect calculations");
217 
218  double returnValue = 0.;
219 
220  double lnDeterminant = 0.;
221  if (this->m_domainSet.contains(domainVector) == false) {
222  // What should the gradient be here?
223  returnValue = -INFINITY;
224  }
225  else {
226  V diffVec(domainVector - this->lawExpVector());
227  if (m_diagonalCovMatrix) {
228  returnValue = ((diffVec*diffVec)/this->lawVarVector()).sumOfComponents();
229 
230  // Compute the gradient of log of the pdf.
231  // The log of a Gaussian pdf is:
232  // f(x) = - 1/2 (x - \mu)^T \Sigma^{-1} (x - \mu)
233  // Therefore
234  // \frac{df}{dx}(x) = - (x - \mu)^T \Sigma^{-1}
235  // = - (\Sigma^{-1}^T (x - \mu))^T
236  // = - (\Sigma^{-1} (x - \mu))^T
237  // = - \Sigma^{-1} (x - \mu) (row/column vector doesn't matter)
238  //
239  // So if \Sigma is diagonal we have a component-wise product of two
240  // vectors (x - \mu) and the diagonal elements of \Sigma^{-1}
241  if (gradVector) {
242  (*gradVector) = diffVec; // Copy
243  (*gradVector) /= this->lawVarVector();
244  (*gradVector) *= -1.0;
245  }
246 
247  if (m_normalizationStyle == 0) {
248  unsigned int iMax = this->lawVarVector().sizeLocal();
249  for (unsigned int i = 0; i < iMax; ++i) {
250  lnDeterminant += std::log(this->lawVarVector()[i]);
251  }
252  }
253  }
254  else {
255  V tmpVec = this->m_lawCovMatrix->invertMultiply(diffVec);
256  returnValue = (diffVec*tmpVec).sumOfComponents();
257 
258  // Compute the gradient of log of the pdf.
259  // The log of a Gaussian pdf is:
260  // f(x) = - 1/2 (x - \mu)^T \Sigma^{-1} (x - \mu)
261  // Therefore
262  // \frac{df}{dx}(x) = - (x - \mu)^T \Sigma^{-1}
263  // = - (\Sigma^{-1}^T (x - \mu))^T
264  // = - (\Sigma^{-1} (x - \mu))^T
265  // = - \Sigma^{-1} (x - \mu) (row/column vector doesn't matter)
266  if (gradVector) {
267  (*gradVector) = tmpVec;
268  (*gradVector) *= -1.0;
269  }
270 
271  if (m_normalizationStyle == 0) {
272  lnDeterminant = this->m_lawCovMatrix->lnDeterminant();
273  }
274  }
275  if (m_normalizationStyle == 0) {
276  returnValue += ((double) this->lawVarVector().sizeLocal()) * std::log(2*M_PI); // normalization of pdf
277  returnValue += lnDeterminant; // normalization of pdf
278  }
279  returnValue *= -0.5;
280  }
281  returnValue += m_logOfNormalizationFactor;
282 
283  if ((m_env.subDisplayFile()) && (m_env.displayVerbosity() >= 99)) {
284  *m_env.subDisplayFile() << "Leaving GaussianJointPdf<V,M>::lnValue()"
285  << ", m_normalizationStyle = " << m_normalizationStyle
286  << ", m_diagonalCovMatrix = " << m_diagonalCovMatrix
287  << ", m_logOfNormalizationFactor = " << m_logOfNormalizationFactor
288  << ", lnDeterminant = " << lnDeterminant
289  << ", meanVector = " << *m_lawExpVector
290  << ", lawCovMatrix = " << *m_lawCovMatrix
291  << ": domainVector = " << domainVector
292  << ", returnValue = " << returnValue
293  << std::endl;
294  }
295 
296  return returnValue;
297 }
const V & lawExpVector() const
Access to the vector of mean values and private attribute: m_lawExpVector.
const VectorSet< V, M > & m_domainSet
Domain set of the scalar function.
const BaseEnvironment & m_env
unsigned int m_normalizationStyle
Definition: JointPdf.h:121
const V & lawVarVector() const
Access to the vector of variance values and private attribute: m_lawVarVector.
double m_logOfNormalizationFactor
Definition: JointPdf.h:122
std::ofstream * subDisplayFile() const
Access function for m_subDisplayFile (displays file on stream).
Definition: Environment.C:320
unsigned int displayVerbosity() const
Definition: Environment.C:450
template<class V , class M >
void QUESO::GaussianJointPdf< V, M >::print ( std::ostream &  os) const
virtual

Print method for informational and logging purposes.

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

Definition at line 125 of file GaussianJointPdf.C.

126 {
127  // Print m_env?
128 
129  os << "Start printing GaussianJointPdf<V, M>" << std::endl;
130  os << "m_prefix:" << std::endl;
131  os << this->m_prefix << std::endl;
132  os << "m_domainSet:" << std::endl;
133  os << this->m_domainSet << std::endl;
134  os << "m_normalizationStyle:" << std::endl;
135  os << this->m_normalizationStyle << std::endl;
136  os << "m_logOfNormalizationFactor:" << std::endl;
137  os << this->m_logOfNormalizationFactor << std::endl;
138  os << "Mean:" << std::endl;
139  os << this->lawExpVector() << std::endl;
140  os << "Variance vector:" << std::endl;
141  os << this->lawVarVector() << std::endl;
142  os << "Covariance matrix:" << std::endl;
143  os << this->lawCovMatrix() << std::endl;
144  os << "Diagonal covariance?" << std::endl;
145  os << this->m_diagonalCovMatrix << std::endl;
146  os << "End printing GaussianJointPdf<V, M>" << std::endl;
147 }
const M & lawCovMatrix() const
Returns the covariance matrix; access to protected attribute m_lawCovMatrix.
const V & lawExpVector() const
Access to the vector of mean values and private attribute: m_lawExpVector.
const VectorSet< V, M > & m_domainSet
Domain set of the scalar function.
unsigned int m_normalizationStyle
Definition: JointPdf.h:121
const V & lawVarVector() const
Access to the vector of variance values and private attribute: m_lawVarVector.
double m_logOfNormalizationFactor
Definition: JointPdf.h:122
template<class V , class M >
void QUESO::GaussianJointPdf< V, M >::updateLawCovMatrix ( const M &  newLawCovMatrix)

Updates the lower triangular matrix from Cholesky decomposition of the covariance matrix to the new value newLowerCholLawCovMatrix.

This method deletes old expected values (allocated at construction or last call to this method).

Definition at line 359 of file GaussianJointPdf.C.

360 {
361  // delete old expected values (allocated at construction or last call to this function)
362  delete m_lawCovMatrix;
363  m_lawCovMatrix = new M(newLawCovMatrix);
364  return;
365 }
template<class V , class M >
void QUESO::GaussianJointPdf< V, M >::updateLawExpVector ( const V &  newLawExpVector)

Updates the mean with the new value newLawExpVector.

This method deletes old expected values (allocated at construction or last call to this method).

Definition at line 349 of file GaussianJointPdf.C.

350 {
351  // delete old expected values (allocated at construction or last call to this function)
352  delete m_lawExpVector;
353  m_lawExpVector = new V(newLawExpVector);
354  return;
355 }

Member Data Documentation

template<class V = GslVector, class M = GslMatrix>
bool QUESO::GaussianJointPdf< V, M >::m_diagonalCovMatrix
protected

Definition at line 126 of file GaussianJointPdf.h.

template<class V = GslVector, class M = GslMatrix>
const M* QUESO::GaussianJointPdf< V, M >::m_lawCovMatrix
protected

Definition at line 127 of file GaussianJointPdf.h.

template<class V = GslVector, class M = GslMatrix>
V* QUESO::GaussianJointPdf< V, M >::m_lawExpVector
protected

Definition at line 124 of file GaussianJointPdf.h.

template<class V = GslVector, class M = GslMatrix>
V* QUESO::GaussianJointPdf< V, M >::m_lawVarVector
protected

Definition at line 125 of file GaussianJointPdf.h.


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