queso-0.53.0
GaussianLikelihoodBlockDiagonalCovariance.C
Go to the documentation of this file.
1 //-----------------------------------------------------------------------bl-
2 //--------------------------------------------------------------------------
3 //
4 // QUESO - a library to support the Quantification of Uncertainty
5 // for Estimation, Simulation and Optimization
6 //
7 // Copyright (C) 2008-2015 The PECOS Development Team
8 //
9 // This library is free software; you can redistribute it and/or
10 // modify it under the terms of the Version 2.1 GNU Lesser General
11 // Public License as published by the Free Software Foundation.
12 //
13 // This library is distributed in the hope that it will be useful,
14 // but WITHOUT ANY WARRANTY; without even the implied warranty of
15 // MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the GNU
16 // Lesser General Public License for more details.
17 //
18 // You should have received a copy of the GNU Lesser General Public
19 // License along with this library; if not, write to the Free Software
20 // Foundation, Inc. 51 Franklin Street, Fifth Floor,
21 // Boston, MA 02110-1301 USA
22 //
23 //-----------------------------------------------------------------------el-
24 
25 #include <cmath>
26 
27 #include <queso/GslVector.h>
28 #include <queso/GslMatrix.h>
29 #include <queso/VectorSet.h>
30 #include <queso/GaussianLikelihoodBlockDiagonalCovariance.h>
31 
32 namespace QUESO {
33 
34 template<class V, class M>
36  const char * prefix, const VectorSet<V, M> & domainSet,
37  const V & observations, const GslBlockMatrix & covariance)
38  : BaseGaussianLikelihood<V, M>(prefix, domainSet, observations),
39  m_covarianceCoefficients(covariance.numBlocks(), 1.0),
40  m_covariance(covariance)
41 {
42  unsigned int totalDim = 0;
43 
44  for (unsigned int i = 0; i < this->m_covariance.numBlocks(); i++) {
45  totalDim += this->m_covariance.getBlock(i).numRowsLocal();
46  }
47 
48  if (totalDim != observations.sizeLocal()) {
49  queso_error_msg("Covariance matrix not same dimension as observation vector");
50  }
51 }
52 
53 template<class V, class M>
55 {
56 }
57 
58 template<class V, class M>
59 double &
61  unsigned int i)
62 {
63  return this->m_covarianceCoefficients[i];
64 }
65 
66 template<class V, class M>
67 const double &
69  unsigned int i) const
70 {
71  return this->m_covarianceCoefficients[i];
72 }
73 
74 template<class V, class M>
75 double
77  const V & domainVector, const V * domainDirection, V * gradVector,
78  M * hessianMatrix, V * hessianEffect) const
79 {
80  return std::exp(this->lnValue(domainVector, domainDirection, gradVector,
81  hessianMatrix, hessianEffect));
82 }
83 
84 template<class V, class M>
85 double
87  const V & domainVector, const V * domainDirection, V * gradVector,
88  M * hessianMatrix, V * hessianEffect) const
89 {
90  V modelOutput(this->m_observations, 0, 0); // At least it's not a copy
91  V weightedMisfit(this->m_observations, 0, 0); // At least it's not a copy
92 
93  this->evaluateModel(domainVector, domainDirection, modelOutput, gradVector,
94  hessianMatrix, hessianEffect);
95 
96  // Compute misfit G(x) - y
97  modelOutput -= this->m_observations;
98 
99  // Solve \Sigma u = G(x) - y for u
100  this->m_covariance.invertMultiply(modelOutput, weightedMisfit);
101 
102  // Deal with the multiplicative coefficients for each of the blocks
103  unsigned int offset = 0;
104 
105  // For each block...
106  for (unsigned int i = 0; i < this->m_covariance.numBlocks(); i++) {
107  // ...divide the appropriate parts of the solution by the coefficient
108  unsigned int blockDim = this->m_covariance.getBlock(i).numRowsLocal();
109  for (unsigned int j = 0; j < blockDim; j++) {
110  // coefficient is a variance, so we divide by it
111  modelOutput[offset+j] /= this->m_covarianceCoefficients[i];
112  }
113  offset += blockDim;
114  }
115 
116  // Compute (G(x) - y)^T \Sigma^{-1} (G(x) - y)
117  modelOutput *= weightedMisfit;
118 
119  double norm2_squared = modelOutput.sumOfComponents(); // This is square of 2-norm
120 
121  return -0.5 * norm2_squared;
122 }
123 
124 } // End namespace QUESO
125 
A class representing a Gaussian likelihood with block-diagonal covariance matrix. ...
A templated class for handling sets.
Definition: VectorSet.h:52
virtual double lnValue(const V &domainVector, const V *domainDirection, V *gradVector, M *hessianMatrix, V *hessianEffect) const
Logarithm of the value of the scalar function.
#define queso_error_msg(msg)
Definition: asserts.h:47
GslMatrix & getBlock(unsigned int i) const
Return block i in the block diagonal matrix.
Class for representing block matrices using GSL library.
unsigned int numRowsLocal() const
Returns the local row dimension of this matrix.
Definition: GslMatrix.C:300
virtual double actualValue(const V &domainVector, const V *domainDirection, V *gradVector, M *hessianMatrix, V *hessianEffect) const
Actual value of the scalar function.
GaussianLikelihoodBlockDiagonalCovariance(const char *prefix, const VectorSet< V, M > &domainSet, const V &observations, const GslBlockMatrix &covariance)
Default constructor.
const double & getBlockCoefficient(unsigned int i) const
Get (const) multiplicative coefficient for block i.
Base class for canned Gaussian likelihoods.
unsigned int numBlocks() const
Return the number of blocks in the block diagonal matrix.
double & blockCoefficient(unsigned int i)
Get (non-const) multiplicative coefficient for block i.

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