24 #include <queso/Defines.h>
26 #ifdef QUESO_HAVE_HDF5
33 #include <queso/Miscellaneous.h>
34 #include <queso/InfiniteDimensionalMeasureBase.h>
36 #include <queso/InfiniteDimensionalLikelihoodBase.h>
37 #include <queso/InfiniteDimensionalMCMCSampler.h>
38 #include <queso/InfiniteDimensionalMCMCSamplerOptions.h>
50 current_physical_state(prior.draw()),
51 proposed_physical_state(prior.draw()),
52 current_physical_mean(prior.draw()),
53 current_physical_var(prior.draw()),
61 #ifdef QUESO_MEMORY_DEBUGGING
62 std::cout <<
"Entering InfiniteDimensionalMCMCSampler class" << std::endl;
71 #ifdef QUESO_MEMORY_DEBUGGING
72 std::cout <<
"In InfiniteDimensionalMCMCSamplerOptions,"
73 <<
" finished scanning options" << std::endl;
79 if ((this->
m_env).subRank() == 0) {
81 (this->m_env).subIdString() +
82 "/test.txt").c_str());
83 queso_require_greater_equal_msg(irtrn, 0,
"unable to verify output path");
87 if ((this->
m_env).subRank() == 0) {
89 (this->m_env).subIdString() +
"/" +
91 H5F_ACC_TRUNC, H5P_DEFAULT, H5P_DEFAULT);
102 ((this->
m_env).fullComm()).Barrier();
129 std::cout.setf(std::ios::fixed, std::ios::floatfield);
130 std::cout.precision(15);
135 if ((this->
m_env).subRank() == 0) {
149 const double coeff = std::sqrt(1.0 - rwmh_step_sq);
169 double diff = this->
_llhd_val - proposed_llhd;
170 double alpha = std::min(1.0, std::exp(diff));
205 this->
_M2->add(1.0, *temp_ptr);
207 this->
_M2->add(-1.0, *temp_ptr);
234 if ((this->
m_env).subRank() == 0) {
238 hsize_t dims[ndims] = {0};
239 hsize_t maxdims[ndims] = {H5S_UNLIMITED};
241 hid_t file_space = H5Screate_simple(ndims, dims, maxdims);
245 hsize_t chunk_dims[ndims] = {1};
246 hid_t plist = H5Pcreate(H5P_DATASET_CREATE);
248 H5Pset_layout(plist, H5D_CHUNKED);
249 H5Pset_chunk(plist, ndims, chunk_dims);
252 hid_t dset = H5Dcreate(this->
_outfile, name.c_str(), H5T_NATIVE_DOUBLE,
253 file_space, H5P_DEFAULT, plist, H5P_DEFAULT);
259 H5Sclose(file_space);
271 if ((this->
m_env).subRank() == 0) {
275 hsize_t dims[ndims] = { 1 };
276 hid_t mem_space = H5Screate_simple(ndims, dims, NULL);
281 err = H5Dset_extent(dset, dims);
282 queso_require_greater_equal_msg(err, 0,
"H5DSet_extent(dset, dims) failed");
285 hid_t file_space = H5Dget_space(dset);
287 hsize_t count[1] = {1};
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");
296 H5Dwrite(dset, H5T_NATIVE_DOUBLE, mem_space, file_space, H5P_DEFAULT, &data);
299 H5Sclose(file_space);
317 std::stringstream basename;
319 basename << (this->
m_env).subIdString();
322 std::ostringstream curr_iter;
325 std::string sample_name(basename.str() +
"sample.e-s.");
326 sample_name += curr_iter.str();
329 std::string mean_name(basename.str() +
"mean.e-s.");
330 mean_name += curr_iter.str();
333 std::string var_name(basename.str() +
"var.e-s.");
334 var_name += curr_iter.str();
348 clone->_avg_acc_prob = 0.0;
375 #endif // QUESO_HAVE_HDF5
std::shared_ptr< T > Type
~InfiniteDimensionalMCMCSampler()
Destructor.
InfiniteDimensionalLikelihoodBase & llhd
This (virtual) class sets up the environment underlying the use of the QUESO library by an executable...
InfiniteDimensionalMeasureBase & prior
SharedPtr< FunctionBase >::Type current_physical_mean
std::string m_dataOutputFileName
Name of the HDF5 output file to store chain statistics.
InfiniteDimensionalMCMCSampler(const BaseEnvironment &env, InfiniteDimensionalMeasureBase &prior, InfiniteDimensionalLikelihoodBase &llhd, InfiniteDimensionalMCMCSamplerOptions *ov)
Construct an infinite dimensional MCMC chain with given prior, likelihood llhd, and options ov...
Abstract class representing the likelihood. Users must subclass this.
InfiniteDimensionalMCMCSamplerOptions * m_ov
double llhd_val() const
Get the current value of the llhd.
std::string m_dataOutputDirName
Name of the output dir to save infinite dimensional output files to.
double m_rwmh_step
The proposal step size.
std::string optionsInputFileName() const
Access to the attribute m_optionsInputFileName, which stores the name of the input file passed by the...
const RngBase * rngObject() const
Access to the RNG object.
SharedPtr< InfiniteDimensionalMCMCSampler >::Type clone_and_reset() const
Returns a pointer to new sampler, with all the moments reset.
void _metropolis_hastings()
hid_t _create_scalar_dataset(const std::string &name)
SharedPtr< FunctionBase >::Type _delta
unsigned int iteration() const
Returns the current iteration number.
int CheckFilePath(const char *path)
double acc_prob()
Returns current acceptance probability.
Abstract base class for infinite dimensional measures.
unsigned int m_save_freq
The frequency at which to save the state of the chain.
void step()
Do one iteration of the Markov chain.
void _append_scalar_dataset(hid_t dset, double data)
This class defines the options that specify the behaviour of the MCMC sampler.
SharedPtr< FunctionBase >::Type current_physical_state
double avg_acc_prob()
Returns current average acceptance probability.
SharedPtr< FunctionBase >::Type _M2
SharedPtr< FunctionBase >::Type current_physical_var
const BaseEnvironment & m_env
std::string m_prefix
The prefix to look for in the input file.
virtual SharedPtr< FunctionBase >::Type draw()=0
Draw from the measure, and then return a shared pointer to the draw.
virtual double uniformSample() const =0
Samples a value from a uniform distribution.
SharedPtr< FunctionBase >::Type proposed_physical_state
virtual double evaluate(FunctionBase &flow)=0
Evaluate the likelihood functional at flow. Subclasses must implement this method.
hid_t _L2_norm_samples_dset
std::ofstream * subDisplayFile() const
Access function for m_subDisplayFile (displays file on stream).