queso-0.53.0
InfiniteDimensionalMCMCSampler.h
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/Defines.h>
26 
27 #ifdef QUESO_HAVE_HDF5
28 
29 #ifndef QUESO_INFINITE_SAMPLER_H
30 #define QUESO_INFINITE_SAMPLER_H
31 
32 // Boost includes
33 #include <boost/scoped_ptr.hpp>
34 #include <boost/shared_ptr.hpp>
35 
36 // Queso includes
37 #include <queso/Environment.h>
38 #include <queso/InfiniteDimensionalMeasureBase.h>
39 #include <queso/InfiniteDimensionalLikelihoodBase.h>
40 #include <queso/InfiniteDimensionalMCMCSamplerOptions.h>
41 
42 // GSL includes
43 #include <gsl/gsl_rng.h>
44 
45 // HDF5 includes
46 #include <hdf5.h>
47 
48 namespace QUESO {
49 
57 class InfiniteDimensionalMCMCSampler
58 {
59 public:
61  InfiniteDimensionalMCMCSampler(
62  const BaseEnvironment& env,
63  InfiniteDimensionalMeasureBase & prior,
64  InfiniteDimensionalLikelihoodBase & llhd,
65  InfiniteDimensionalMCMCSamplerOptions * ov);
66 
68  ~InfiniteDimensionalMCMCSampler();
69 
71  void step();
72 
74  double llhd_val() const;
75 
77  double acc_prob();
78 
80  double avg_acc_prob();
81 
83  unsigned int iteration() const;
84 
86  boost::shared_ptr<InfiniteDimensionalMCMCSampler> clone_and_reset() const;
87 
88 private:
89  // Current iteration
90  unsigned int _iteration;
91 
92  // The current value of the negative log-likelihood functional
93  double _llhd_val;
94 
95  // The current acceptance probability
96  double _acc_prob;
97 
98  // The current running mean of the acceptance probability
99  double _avg_acc_prob;
100 
101  // The prior measure from which to draw
102  InfiniteDimensionalMeasureBase & prior;
103 
104  // The negative log-likelihood functional. Operates on functions.
105  // Should be const?
106  InfiniteDimensionalLikelihoodBase & llhd;
107 
108  // Aggregate options object
109  InfiniteDimensionalMCMCSamplerOptions * m_ov;
110 
111  // The QUESO environment
112  const BaseEnvironment& m_env;
113 
114  // Pointer to the current physical state
115  boost::shared_ptr<FunctionBase> current_physical_state;
116 
117  // Pointer to the current proposed state
118  boost::shared_ptr<FunctionBase> proposed_physical_state;
119 
120  // Pointer to the current physical mean
121  boost::shared_ptr<FunctionBase> current_physical_mean;
122 
123  // Pointer to the current physical variance
124  boost::shared_ptr<FunctionBase> current_physical_var;
125 
126  // Stores the differences from the mean
127  boost::shared_ptr<FunctionBase> _delta;
128 
129  // Stores a running sum-of-squares (kinda)
130  boost::shared_ptr<FunctionBase> _M2;
131 
132  // A pointer to the random number generator to use.
133  // Should probably use the one in the queso environment.
134  gsl_rng *r;
135 
136  // HDF5 file identifier
137  hid_t _outfile;
138 
139  // HDF5 datasets
140  hid_t _acc_dset;
141  hid_t _avg_acc_dset;
142  hid_t _neg_log_llhd_dset;
143  hid_t _L2_norm_samples_dset;
144  hid_t _L2_norm_mean_dset;
145  hid_t _L2_norm_var_dset;
146 
150  void _propose();
151 
152  // Do a metropolis random walk step
153  void _metropolis_hastings();
154 
155  // Increments the current iteration and updates the running mean and variance.
156  void _update_moments();
157 
158  // Write the current state of the chain to disk
159  void _write_state();
160 
161  // Creates a scalar dataset called \c name in \c _outfile file and returns it
162  hid_t _create_scalar_dataset(const std::string & name);
163 
164  // Appends to a scalar dataset in the hdf5 file
165  void _append_scalar_dataset(hid_t dset, double data);
166 };
167 
168 } // End namespace QUESO
169 
170 #endif // QUESO_INFINITE_SAMPLER_H
171 
172 #endif // QUESO_HAVE_HDF5

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