queso-0.57.0
Public Member Functions | Private Member Functions | Private Attributes | List of all members
QUESO::InfiniteDimensionalMCMCSampler Class Reference

Class representing the infinite dimensional Markov chain Monte Carlo sampler. More...

#include <InfiniteDimensionalMCMCSampler.h>

Collaboration diagram for QUESO::InfiniteDimensionalMCMCSampler:
Collaboration graph
[legend]

Public Member Functions

 InfiniteDimensionalMCMCSampler (const BaseEnvironment &env, InfiniteDimensionalMeasureBase &prior, InfiniteDimensionalLikelihoodBase &llhd, InfiniteDimensionalMCMCSamplerOptions *ov)
 Construct an infinite dimensional MCMC chain with given prior, likelihood llhd, and options ov. More...
 
 ~InfiniteDimensionalMCMCSampler ()
 Destructor. More...
 
void step ()
 Do one iteration of the Markov chain. More...
 
double llhd_val () const
 Get the current value of the llhd. More...
 
double acc_prob ()
 Returns current acceptance probability. More...
 
double avg_acc_prob ()
 Returns current average acceptance probability. More...
 
unsigned int iteration () const
 Returns the current iteration number. More...
 
SharedPtr
< InfiniteDimensionalMCMCSampler >
::Type 
clone_and_reset () const
 Returns a pointer to new sampler, with all the moments reset. More...
 

Private Member Functions

void _propose ()
 
void _metropolis_hastings ()
 
void _update_moments ()
 
void _write_state ()
 
hid_t _create_scalar_dataset (const std::string &name)
 
void _append_scalar_dataset (hid_t dset, double data)
 

Private Attributes

unsigned int _iteration
 
double _llhd_val
 
double _acc_prob
 
double _avg_acc_prob
 
InfiniteDimensionalMeasureBaseprior
 
InfiniteDimensionalLikelihoodBasellhd
 
InfiniteDimensionalMCMCSamplerOptionsm_ov
 
const BaseEnvironmentm_env
 
SharedPtr< FunctionBase >::Type current_physical_state
 
SharedPtr< FunctionBase >::Type proposed_physical_state
 
SharedPtr< FunctionBase >::Type current_physical_mean
 
SharedPtr< FunctionBase >::Type current_physical_var
 
SharedPtr< FunctionBase >::Type _delta
 
SharedPtr< FunctionBase >::Type _M2
 
hid_t _outfile
 
hid_t _acc_dset
 
hid_t _avg_acc_dset
 
hid_t _neg_log_llhd_dset
 
hid_t _L2_norm_samples_dset
 
hid_t _L2_norm_mean_dset
 
hid_t _L2_norm_var_dset
 

Detailed Description

Class representing the infinite dimensional Markov chain Monte Carlo sampler.

Definition at line 51 of file InfiniteDimensionalMCMCSampler.h.

Constructor & Destructor Documentation

QUESO::InfiniteDimensionalMCMCSampler::InfiniteDimensionalMCMCSampler ( const BaseEnvironment env,
InfiniteDimensionalMeasureBase prior,
InfiniteDimensionalLikelihoodBase llhd,
InfiniteDimensionalMCMCSamplerOptions ov 
)

Construct an infinite dimensional MCMC chain with given prior, likelihood llhd, and options ov.

Definition at line 42 of file InfiniteDimensionalMCMCSampler.C.

References _acc_dset, _acc_prob, _avg_acc_dset, _avg_acc_prob, _create_scalar_dataset(), _delta, _iteration, _L2_norm_mean_dset, _L2_norm_samples_dset, _L2_norm_var_dset, _llhd_val, _M2, _neg_log_llhd_dset, _outfile, QUESO::CheckFilePath(), current_physical_mean, current_physical_state, current_physical_var, QUESO::InfiniteDimensionalLikelihoodBase::evaluate(), QUESO::InfiniteDimensionalMCMCSamplerOptions::m_dataOutputDirName, QUESO::InfiniteDimensionalMCMCSamplerOptions::m_dataOutputFileName, m_env, m_ov, QUESO::InfiniteDimensionalMCMCSamplerOptions::m_prefix, QUESO::BaseEnvironment::optionsInputFileName(), and QUESO::BaseEnvironment::subDisplayFile().

Referenced by clone_and_reset().

47  : prior(prior),
48  llhd(llhd),
49  m_env(env),
54  _delta(prior.draw()),
55  _M2(prior.draw())
56 {
57  if (ov != NULL) {
58  this->m_ov = ov;
59  }
60 
61 #ifdef QUESO_MEMORY_DEBUGGING
62  std::cout << "Entering InfiniteDimensionalMCMCSampler class" << std::endl;
63 #endif
64  if (m_env.subDisplayFile()) {
65  *m_env.subDisplayFile() << "Entering InfiniteDimensionalMCMCSampler::constructor()"
66  << ": prefix = " << this->m_ov->m_prefix
67  << ", m_env.optionsInputFileName() = " << this->m_env.optionsInputFileName()
68  << std::endl;
69  }
70 
71 #ifdef QUESO_MEMORY_DEBUGGING
72  std::cout << "In InfiniteDimensionalMCMCSamplerOptions,"
73  << " finished scanning options" << std::endl;
74 #endif
75 
76  // Verify parent directory exists (for cases when a user specifies a
77  // relative path for the desired output file).
78  // Only subprocess of subrank 0 creates the output file
79  if ((this->m_env).subRank() == 0) {
80  int irtrn = CheckFilePath((this->m_ov->m_dataOutputDirName +
81  (this->m_env).subIdString() +
82  "/test.txt").c_str());
83  queso_require_greater_equal_msg(irtrn, 0, "unable to verify output path");
84  }
85 
86  // Only subprocess of subrank 0 creates the output file
87  if ((this->m_env).subRank() == 0) {
88  this->_outfile = H5Fcreate((this->m_ov->m_dataOutputDirName +
89  (this->m_env).subIdString() + "/" +
90  this->m_ov->m_dataOutputFileName).c_str(),
91  H5F_ACC_TRUNC, H5P_DEFAULT, H5P_DEFAULT);
92  }
93 
94  this->_acc_dset = this->_create_scalar_dataset("acc");
95  this->_avg_acc_dset = this->_create_scalar_dataset("avg_acc");
96  this->_neg_log_llhd_dset = this->_create_scalar_dataset("neg_log_llhd");
97  this->_L2_norm_samples_dset = this->_create_scalar_dataset("L2_norm_samples");
98  this->_L2_norm_mean_dset = this->_create_scalar_dataset("L2_norm_mean");
99  this->_L2_norm_var_dset = this->_create_scalar_dataset("L2_norm_var");
100 
101  // Ensure that we created the path and the output files
102  ((this->m_env).fullComm()).Barrier();
103 
104  this->_iteration = 0;
105  this->_acc_prob = 0.0;
106  this->_avg_acc_prob = 0.0;
107 
108  // Zero out these guys. There's probably a better way of doing this than
109  // calling prior.draw() at the start, but creation of a Sampler object
110  // should only be done a O(1) times anyway.
111  this->_delta->zero();
112  this->_M2->zero();
113  this->current_physical_mean->zero();
114  this->current_physical_var->zero();
115 
116  // LibMeshFunction & p = libmesh_cast_ref<LibMeshFunction &>(*(this->current_physical_state));
117 
118  // Seed the sampler at the truth
119  // for (unsigned int ii = 0; ii < 513; ii++) {
120  // p.equation_systems->get_system<ExplicitSystem>("Function").solution->set(ii, std::cos(2.0 * boost::math::constants::pi<double>() * ii / 512.0));
121  // }
122  // p.equation_systems->get_system<ExplicitSystem>("Function").solution->close();
123  // std::string seedname("seedfn");
124  // p.save_function(seedname);
125 
126  // Initialise cost function
127  this->_llhd_val = this->llhd.evaluate(*(this->current_physical_state));
128 
129  std::cout.setf(std::ios::fixed, std::ios::floatfield);
130  std::cout.precision(15);
131 }
int CheckFilePath(const char *path)
SharedPtr< FunctionBase >::Type current_physical_state
std::string optionsInputFileName() const
Access to the attribute m_optionsInputFileName, which stores the name of the input file passed by the...
Definition: Environment.C:354
virtual double evaluate(FunctionBase &flow)=0
Evaluate the likelihood functional at flow. Subclasses must implement this method.
SharedPtr< FunctionBase >::Type current_physical_var
std::string m_dataOutputDirName
Name of the output dir to save infinite dimensional output files to.
virtual SharedPtr< FunctionBase >::Type draw()=0
Draw from the measure, and then return a shared pointer to the draw.
SharedPtr< FunctionBase >::Type proposed_physical_state
SharedPtr< FunctionBase >::Type current_physical_mean
InfiniteDimensionalMCMCSamplerOptions * m_ov
std::string m_prefix
The prefix to look for in the input file.
std::string m_dataOutputFileName
Name of the HDF5 output file to store chain statistics.
std::ofstream * subDisplayFile() const
Access function for m_subDisplayFile (displays file on stream).
Definition: Environment.C:320
hid_t _create_scalar_dataset(const std::string &name)
QUESO::InfiniteDimensionalMCMCSampler::~InfiniteDimensionalMCMCSampler ( )

Destructor.

Definition at line 133 of file InfiniteDimensionalMCMCSampler.C.

References _acc_dset, _avg_acc_dset, _L2_norm_mean_dset, _L2_norm_samples_dset, _L2_norm_var_dset, _neg_log_llhd_dset, _outfile, and m_env.

Member Function Documentation

void QUESO::InfiniteDimensionalMCMCSampler::_append_scalar_dataset ( hid_t  dset,
double  data 
)
private

Definition at line 268 of file InfiniteDimensionalMCMCSampler.C.

References _iteration, m_env, m_ov, and QUESO::InfiniteDimensionalMCMCSamplerOptions::m_save_freq.

Referenced by _write_state().

269 {
270  // Only subprocess with rank 0 manipulates the output file
271  if ((this->m_env).subRank() == 0) {
272  int err;
273  // Create a memory dataspace for data to append
274  const int ndims = 1;
275  hsize_t dims[ndims] = { 1 }; // Only writing one double
276  hid_t mem_space = H5Screate_simple(ndims, dims, NULL);
277 
278  // Extend the dataset
279  // Set dims to be the *new* dimension of the extended dataset
280  dims[0] = this->_iteration / this->m_ov->m_save_freq;
281  err = H5Dset_extent(dset, dims);
282  queso_require_greater_equal_msg(err, 0, "H5DSet_extent(dset, dims) failed");
283 
284  // Select hyperslab on file dataset
285  hid_t file_space = H5Dget_space(dset);
286  hsize_t start[1] = {(this->_iteration / this->m_ov->m_save_freq) - 1};
287  hsize_t count[1] = {1};
288 
289  err = H5Sselect_hyperslab(file_space, H5S_SELECT_SET, start, NULL, count, NULL);
290  queso_require_greater_equal_msg(err, 0, "H5Sselect_hyperslab(file_space, H5S_SELECT_SET, start, NULL, count, NULL) failed");
291 
292  // hsize_t size[1];
293  // size[0] = this->_iteration / this->m_ov->m_save_freq;
294 
295  // Write the data
296  H5Dwrite(dset, H5T_NATIVE_DOUBLE, mem_space, file_space, H5P_DEFAULT, &data);
297 
298  // Close a bunch of stuff
299  H5Sclose(file_space);
300  H5Sclose(mem_space);
301  }
302 }
unsigned int m_save_freq
The frequency at which to save the state of the chain.
InfiniteDimensionalMCMCSamplerOptions * m_ov
hid_t QUESO::InfiniteDimensionalMCMCSampler::_create_scalar_dataset ( const std::string &  name)
private

Definition at line 231 of file InfiniteDimensionalMCMCSampler.C.

References _outfile, and m_env.

Referenced by InfiniteDimensionalMCMCSampler().

232 {
233  // Only subprocess with rank 0 manipulates the output file
234  if ((this->m_env).subRank() == 0) {
235  // Create a 1D dataspace. Unlimited size. Initially set to 0. We will
236  // extend it later
237  const int ndims = 1;
238  hsize_t dims[ndims] = {0}; // dataset dimensions at creation
239  hsize_t maxdims[ndims] = {H5S_UNLIMITED};
240 
241  hid_t file_space = H5Screate_simple(ndims, dims, maxdims);
242 
243  // Create dataset creation property list. Unlimited datasets must be
244  // chunked. Choosing the chunk size is an issue, here we set it to 1
245  hsize_t chunk_dims[ndims] = {1};
246  hid_t plist = H5Pcreate(H5P_DATASET_CREATE);
247 
248  H5Pset_layout(plist, H5D_CHUNKED);
249  H5Pset_chunk(plist, ndims, chunk_dims);
250 
251  // Create the dataset
252  hid_t dset = H5Dcreate(this->_outfile, name.c_str(), H5T_NATIVE_DOUBLE,
253  file_space, H5P_DEFAULT, plist, H5P_DEFAULT);
254 
255  // We don't need the property list anymore. We also don't need the file
256  // dataspace anymore because we'll extend it later, making this one
257  // invalild anyway.
258  H5Pclose(plist);
259  H5Sclose(file_space);
260 
261  return dset;
262  }
263 
264  hid_t dummy = -1;
265  return dummy;
266 }
void QUESO::InfiniteDimensionalMCMCSampler::_metropolis_hastings ( )
private

Definition at line 158 of file InfiniteDimensionalMCMCSampler.C.

References _acc_prob, _llhd_val, current_physical_state, QUESO::InfiniteDimensionalLikelihoodBase::evaluate(), llhd, m_env, proposed_physical_state, QUESO::BaseEnvironment::rngObject(), and QUESO::RngBase::uniformSample().

Referenced by step().

159 {
160  // Downcast since we know we're dealing with a libmesh function
161  // LibMeshFunction & p = static_cast<LibMeshFunction &>(*(this->proposed_physical_state));
162  // LibMeshFunction & q = static_cast<LibMeshFunction &>(*(this->current_physical_state));
163  // Evaluate the likelihood at the proposed state
164  double proposed_llhd = this->llhd.evaluate(*(this->proposed_physical_state));
165  // double current_llhd = this->llhd.evaluate(q);
166  // std::cout << "llhd of current is: " << current_llhd << std::endl;
167  // std::cout << "llhd of proposal is: " << proposed_llhd << std::endl;
168 
169  double diff = this->_llhd_val - proposed_llhd;
170  double alpha = std::min(1.0, std::exp(diff));
171  double rand = this->m_env.rngObject()->uniformSample();
172  if (rand < alpha) {
173  // Accept
174  // std::cout << "accepted" << std::endl;
176  this->proposed_physical_state->zero();
177  this->_acc_prob = alpha;
178  this->_llhd_val = proposed_llhd;
179  }
180  else {
181  // std::cout << "rejected" << std::endl;
182  this->_acc_prob = 0.0;
183  }
184 }
SharedPtr< FunctionBase >::Type current_physical_state
virtual double evaluate(FunctionBase &flow)=0
Evaluate the likelihood functional at flow. Subclasses must implement this method.
const RngBase * rngObject() const
Access to the RNG object.
Definition: Environment.C:471
SharedPtr< FunctionBase >::Type proposed_physical_state
virtual double uniformSample() const =0
Samples a value from a uniform distribution.
void QUESO::InfiniteDimensionalMCMCSampler::_propose ( )
private

Make a proposal from the prior using a standard random walk

Definition at line 146 of file InfiniteDimensionalMCMCSampler.C.

References current_physical_state, QUESO::InfiniteDimensionalMeasureBase::draw(), m_ov, QUESO::InfiniteDimensionalMCMCSamplerOptions::m_rwmh_step, prior, and proposed_physical_state.

Referenced by step().

147 {
148  const double rwmh_step_sq = (this->m_ov->m_rwmh_step * this->m_ov->m_rwmh_step);
149  const double coeff = std::sqrt(1.0 - rwmh_step_sq);
150 
152 
153  this->proposed_physical_state->zero();
154  this->proposed_physical_state->add(coeff, *(this->current_physical_state));
155  this->proposed_physical_state->add(this->m_ov->m_rwmh_step, *p);
156 }
SharedPtr< FunctionBase >::Type current_physical_state
virtual SharedPtr< FunctionBase >::Type draw()=0
Draw from the measure, and then return a shared pointer to the draw.
SharedPtr< FunctionBase >::Type proposed_physical_state
InfiniteDimensionalMCMCSamplerOptions * m_ov
Definition of a shared pointer.
void QUESO::InfiniteDimensionalMCMCSampler::_update_moments ( )
private

Definition at line 186 of file InfiniteDimensionalMCMCSampler.C.

References _avg_acc_prob, _delta, _iteration, _M2, acc_prob(), avg_acc_prob(), current_physical_mean, current_physical_state, current_physical_var, and iteration().

Referenced by step().

187 {
188  // Increment the current iteration number and update the running mean and
189  // variance.
190  this->_iteration += 1;
191 
192  // Make _delta contain the difference from the mean
193  this->_delta->zero();
194  this->_delta->add(1.0, *(this->current_physical_state));
195  this->_delta->add(-1.0, *(this->current_physical_mean));
196 
197  // Scale _delta to update the mean field
198  this->current_physical_mean->add(1.0 / this->iteration(), *(this->_delta));
199 
200  // Update running sum-of-squares
201  SharedPtr<FunctionBase>::Type temp_ptr(this->_delta->zero_clone());
202  // LibMeshFunction & temp = static_cast<LibMeshFunction &>(*temp_ptr);
203 
204  temp_ptr->pointwise_mult(*(this->_delta), *(this->current_physical_state));
205  this->_M2->add(1.0, *temp_ptr);
206  temp_ptr->pointwise_mult(*(this->_delta), *(this->current_physical_mean));
207  this->_M2->add(-1.0, *temp_ptr);
208 
209  if (this->iteration() > 1) {
210  this->current_physical_var->zero();
211  this->current_physical_var->add(1.0 / (this->iteration() - 1), *(this->_M2));
212  }
213 
214  // Update acceptance rate
215  double delta_acc = this->acc_prob() - this->avg_acc_prob();
216  this->_avg_acc_prob += delta_acc / this->iteration();
217 }
SharedPtr< FunctionBase >::Type current_physical_state
double avg_acc_prob()
Returns current average acceptance probability.
SharedPtr< FunctionBase >::Type current_physical_var
SharedPtr< FunctionBase >::Type current_physical_mean
Definition of a shared pointer.
unsigned int iteration() const
Returns the current iteration number.
double acc_prob()
Returns current acceptance probability.
void QUESO::InfiniteDimensionalMCMCSampler::_write_state ( )
private

Definition at line 304 of file InfiniteDimensionalMCMCSampler.C.

References _acc_dset, _append_scalar_dataset(), _avg_acc_dset, _L2_norm_mean_dset, _L2_norm_samples_dset, _L2_norm_var_dset, _llhd_val, _neg_log_llhd_dset, acc_prob(), avg_acc_prob(), current_physical_mean, current_physical_state, current_physical_var, iteration(), QUESO::InfiniteDimensionalMCMCSamplerOptions::m_dataOutputDirName, m_env, and m_ov.

Referenced by step().

305 {
306  this->_append_scalar_dataset(this->_acc_dset, this->acc_prob());
307  this->_append_scalar_dataset(this->_avg_acc_dset, this->avg_acc_prob());
310  this->_append_scalar_dataset(this->_L2_norm_mean_dset, this->current_physical_mean->L2_norm());
311  this->_append_scalar_dataset(this->_L2_norm_var_dset, this->current_physical_var->L2_norm());
312 
313  // Now to write the fields. It appears to be a pain in the arse to write a
314  // method to spit this out to HDF5 format. Also, this won't scale to
315  // non-uniform finite element meshes. Therefore, I'm going to spit out an
316  // ExodusII file for each sample, average and variance. Got disk space?
317  std::stringstream basename;
318  basename << this->m_ov->m_dataOutputDirName;
319  basename << (this->m_env).subIdString();
320  basename << "/"; // Sigh
321 
322  std::ostringstream curr_iter;
323  curr_iter << this->iteration();
324 
325  std::string sample_name(basename.str() + "sample.e-s.");
326  sample_name += curr_iter.str();
327  this->current_physical_state->save_function(sample_name, this->iteration());
328 
329  std::string mean_name(basename.str() + "mean.e-s.");
330  mean_name += curr_iter.str();
331  this->current_physical_mean->save_function(mean_name, this->iteration());
332 
333  std::string var_name(basename.str() + "var.e-s.");
334  var_name += curr_iter.str();
335  this->current_physical_var->save_function(var_name, this->iteration());
336 }
SharedPtr< FunctionBase >::Type current_physical_state
double avg_acc_prob()
Returns current average acceptance probability.
SharedPtr< FunctionBase >::Type current_physical_var
std::string m_dataOutputDirName
Name of the output dir to save infinite dimensional output files to.
SharedPtr< FunctionBase >::Type current_physical_mean
InfiniteDimensionalMCMCSamplerOptions * m_ov
unsigned int iteration() const
Returns the current iteration number.
double acc_prob()
Returns current acceptance probability.
double QUESO::InfiniteDimensionalMCMCSampler::acc_prob ( )

Returns current acceptance probability.

Definition at line 353 of file InfiniteDimensionalMCMCSampler.C.

References _acc_prob.

Referenced by _update_moments(), and _write_state().

354 {
355  return this->_acc_prob;
356 }
double QUESO::InfiniteDimensionalMCMCSampler::avg_acc_prob ( )

Returns current average acceptance probability.

Definition at line 358 of file InfiniteDimensionalMCMCSampler.C.

References _avg_acc_prob.

Referenced by _update_moments(), and _write_state().

359 {
360  return this->_avg_acc_prob;
361 }
SharedPtr< InfiniteDimensionalMCMCSampler >::Type QUESO::InfiniteDimensionalMCMCSampler::clone_and_reset ( ) const

Returns a pointer to new sampler, with all the moments reset.

Definition at line 338 of file InfiniteDimensionalMCMCSampler.C.

References _acc_prob, _llhd_val, current_physical_state, InfiniteDimensionalMCMCSampler(), llhd, m_env, m_ov, prior, and proposed_physical_state.

339 {
340  // Set up a clone
342 
343  // Copy the state.
344  clone->current_physical_state = this->current_physical_state;
345  clone->proposed_physical_state = this->proposed_physical_state;
346  clone->_llhd_val = this->_llhd_val;
347  clone->_acc_prob = this->_acc_prob;
348  clone->_avg_acc_prob = 0.0;
349 
350  return clone;
351 }
SharedPtr< FunctionBase >::Type current_physical_state
SharedPtr< FunctionBase >::Type proposed_physical_state
InfiniteDimensionalMCMCSamplerOptions * m_ov
Definition of a shared pointer.
InfiniteDimensionalMCMCSampler(const BaseEnvironment &env, InfiniteDimensionalMeasureBase &prior, InfiniteDimensionalLikelihoodBase &llhd, InfiniteDimensionalMCMCSamplerOptions *ov)
Construct an infinite dimensional MCMC chain with given prior, likelihood llhd, and options ov...
unsigned int QUESO::InfiniteDimensionalMCMCSampler::iteration ( ) const

Returns the current iteration number.

Definition at line 368 of file InfiniteDimensionalMCMCSampler.C.

References _iteration.

Referenced by _update_moments(), and _write_state().

369 {
370  return this->_iteration;
371 }
double QUESO::InfiniteDimensionalMCMCSampler::llhd_val ( ) const

Get the current value of the llhd.

Definition at line 363 of file InfiniteDimensionalMCMCSampler.C.

References _llhd_val.

364 {
365  return this->_llhd_val;
366 }
void QUESO::InfiniteDimensionalMCMCSampler::step ( )

Do one iteration of the Markov chain.

Definition at line 219 of file InfiniteDimensionalMCMCSampler.C.

References _iteration, _metropolis_hastings(), _propose(), _update_moments(), _write_state(), m_ov, and QUESO::InfiniteDimensionalMCMCSamplerOptions::m_save_freq.

220 {
221  this->_propose();
222  this->_metropolis_hastings();
223  this->_update_moments();
224 
225  // We never save the 0th iteration
226  if (this->_iteration % this->m_ov->m_save_freq == 0) {
227  this->_write_state();
228  }
229 }
unsigned int m_save_freq
The frequency at which to save the state of the chain.
InfiniteDimensionalMCMCSamplerOptions * m_ov

Member Data Documentation

hid_t QUESO::InfiniteDimensionalMCMCSampler::_acc_dset
private
double QUESO::InfiniteDimensionalMCMCSampler::_acc_prob
private
hid_t QUESO::InfiniteDimensionalMCMCSampler::_avg_acc_dset
private
double QUESO::InfiniteDimensionalMCMCSampler::_avg_acc_prob
private
SharedPtr<FunctionBase>::Type QUESO::InfiniteDimensionalMCMCSampler::_delta
private
unsigned int QUESO::InfiniteDimensionalMCMCSampler::_iteration
private
hid_t QUESO::InfiniteDimensionalMCMCSampler::_L2_norm_mean_dset
private
hid_t QUESO::InfiniteDimensionalMCMCSampler::_L2_norm_samples_dset
private
hid_t QUESO::InfiniteDimensionalMCMCSampler::_L2_norm_var_dset
private
double QUESO::InfiniteDimensionalMCMCSampler::_llhd_val
private
SharedPtr<FunctionBase>::Type QUESO::InfiniteDimensionalMCMCSampler::_M2
private
hid_t QUESO::InfiniteDimensionalMCMCSampler::_neg_log_llhd_dset
private
hid_t QUESO::InfiniteDimensionalMCMCSampler::_outfile
private
SharedPtr<FunctionBase>::Type QUESO::InfiniteDimensionalMCMCSampler::current_physical_mean
private
SharedPtr<FunctionBase>::Type QUESO::InfiniteDimensionalMCMCSampler::current_physical_state
private
SharedPtr<FunctionBase>::Type QUESO::InfiniteDimensionalMCMCSampler::current_physical_var
private
InfiniteDimensionalLikelihoodBase& QUESO::InfiniteDimensionalMCMCSampler::llhd
private

Definition at line 100 of file InfiniteDimensionalMCMCSampler.h.

Referenced by _metropolis_hastings(), and clone_and_reset().

const BaseEnvironment& QUESO::InfiniteDimensionalMCMCSampler::m_env
private
InfiniteDimensionalMCMCSamplerOptions* QUESO::InfiniteDimensionalMCMCSampler::m_ov
private
InfiniteDimensionalMeasureBase& QUESO::InfiniteDimensionalMCMCSampler::prior
private

Definition at line 96 of file InfiniteDimensionalMCMCSampler.h.

Referenced by _propose(), and clone_and_reset().

SharedPtr<FunctionBase>::Type QUESO::InfiniteDimensionalMCMCSampler::proposed_physical_state
private

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

Generated on Sat Apr 22 2017 14:04:39 for queso-0.57.0 by  doxygen 1.8.5