queso-0.56.1
GaussianLikelihoodBlockDiagonalCovarianceRandomCoefficients.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 <queso/GslVector.h>
26 #include <queso/GslMatrix.h>
27 #include <queso/VectorSet.h>
28 #include <queso/GaussianLikelihoodBlockDiagonalCovarianceRandomCoefficients.h>
29 
30 namespace QUESO {
31 
32 template<class V, class M>
34  const char * prefix, const VectorSet<V, M> & domainSet,
35  const V & observations, const GslBlockMatrix & covariance)
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 }
49 
50 template<class V, class M>
52 {
53 }
54 
55 template<class V, class M>
56 double
58  const V & domainVector, const V * domainDirection, V * gradVector,
59  M * hessianMatrix, V * hessianEffect) const
60 {
61  return std::exp(this->lnValue(domainVector, domainDirection, gradVector,
62  hessianMatrix, hessianEffect));
63 }
64 
65 template<class V, class M>
66 double
68  const V & domainVector, const V * domainDirection, V * gradVector,
69  M * hessianMatrix, V * hessianEffect) const
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 }
120 
121 } // End namespace QUESO
122 
Class for representing block matrices using GSL library.
GslMatrix & getBlock(unsigned int i) const
Return block i in the block diagonal matrix.
GaussianLikelihoodBlockDiagonalCovarianceRandomCoefficients(const char *prefix, const VectorSet< V, M > &domainSet, const V &observations, const GslBlockMatrix &covariance)
Default constructor.
virtual double actualValue(const V &domainVector, const V *domainDirection, V *gradVector, M *hessianMatrix, V *hessianEffect) const
Actual value of the scalar function.
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 double lnValue(const V &domainVector, const V *domainDirection, V *gradVector, M *hessianMatrix, V *hessianEffect) const
Logarithm of the value of the scalar function.
A class representing a Gaussian likelihood with block-diagonal covariance matrix. ...
A templated class for handling sets.
Definition: VectorSet.h:52
#define queso_error_msg(msg)
Definition: asserts.h:47
Base class for canned Gaussian likelihoods.

Generated on Thu Dec 15 2016 13:23:11 for queso-0.56.1 by  doxygen 1.8.5