Skip to content

Commit a52f42d

Browse files
committed
Likelihood evaluation can now do marginalization
If we detect the user has supplied a marinal parameter pdf, then we instead will use the user-supplied integration rule to integrate over the parameter space, calling the user-overrideen evaluateModel() method that accepts marginal parameters as well as inference parameters.
1 parent 0d077b4 commit a52f42d

File tree

1 file changed

+46
-3
lines changed

1 file changed

+46
-3
lines changed

src/stats/src/LikelihoodBase.C

Lines changed: 46 additions & 3 deletions
Original file line numberDiff line numberDiff line change
@@ -98,11 +98,54 @@ template<class V, class M>
9898
double
9999
LikelihoodBase<V, M>::lnValue(const V & domainVector) const
100100
{
101-
V modelOutput(this->m_observations, 0, 0); // At least it's not a copy
101+
V modelOutput(this->m_observations, 0, 0);
102102

103-
this->evaluateModel(domainVector, modelOutput);
103+
double lnLikelihood_value = 0.0;
104104

105-
return this->lnLikelihood(domainVector, modelOutput);
105+
// If we're not marginalizing, then we just need f(m)
106+
if( !m_marg_param_pdf )
107+
{
108+
this->evaluateModel(domainVector, modelOutput);
109+
110+
// Note modelOutput made be modified in lnLikelihood()
111+
lnLikelihood_value = this->lnLikelihood(domainVector, modelOutput);
112+
}
113+
// Otherwise we're integrating over marginal parameter space
114+
else
115+
{
116+
queso_assert(m_marg_integration);
117+
118+
const std::vector< typename QUESO::SharedPtr< V >::Type > & x =
119+
this->m_marg_integration->positions();
120+
121+
const std::vector< double > & w = this->m_marg_integration->weights();
122+
123+
unsigned int n_qpoints = x.size();
124+
125+
for( unsigned int q = 0; q < n_qpoints; q++ )
126+
{
127+
this->evaluateModel(domainVector, *(x[q]), modelOutput);
128+
129+
// Note modelOutput made be modified in lnLikelihood()
130+
double ln_pi_m_q = this->lnLikelihood(domainVector, modelOutput);
131+
132+
double ln_pi_q = 0.0;
133+
134+
// If we're not treating the marginal parameters pdf as
135+
// the weighting function in the quadrature evaluation,
136+
// then we need to also evaluate it.
137+
if(!m_marg_pdf_is_weight_func)
138+
this->m_marg_param_pdf->pdf().lnValue( *(x[q]) );
139+
140+
/*! \todo [PB]: We might want to play games with a log(\sum) identity
141+
if precision starts to become an issue. */
142+
lnLikelihood_value += std::exp(ln_pi_m_q + ln_pi_q)*w[q];
143+
}
144+
145+
lnLikelihood_value = std::log(lnLikelihood_value);
146+
}
147+
148+
return lnLikelihood_value;
106149
}
107150

108151
} // End namespace QUESO

0 commit comments

Comments
 (0)