Stan  1.0
probability, sampling & optimization
hmc.hpp
Go to the documentation of this file.
1 #ifndef __STAN__MCMC__HMC_HPP__
2 #define __STAN__MCMC__HMC_HPP__
3 
4 #include <ctime>
5 #include <cstddef>
6 #include <vector>
7 
8 #include <boost/random/mersenne_twister.hpp>
9 #include <boost/random/normal_distribution.hpp>
10 #include <boost/random/uniform_01.hpp>
11 #include <boost/random/variate_generator.hpp>
12 
13 #include <stan/model/prob_grad.hpp>
14 #include <stan/mcmc/sampler.hpp>
16 #include <stan/mcmc/util.hpp>
17 #include <stan/math/util.hpp>
18 
19 namespace stan {
20 
21  namespace mcmc {
22 
23 
34  class hmc : public adaptive_sampler {
35  private:
36  stan::model::prob_grad& _model;
37 
38  std::vector<double> _x;
39  std::vector<int> _z;
40  std::vector<double> _g;
41  double _logp;
42 
43  double _epsilon;
44  unsigned int _L;
45 
46  boost::mt19937 _rand_int;
47  boost::variate_generator<boost::mt19937&, boost::normal_distribution<> > _rand_unit_norm;
48  boost::uniform_01<boost::mt19937&> _rand_uniform_01;
49 
50  public:
51 
69  double epsilon,
70  int L,
71  unsigned int random_seed = static_cast<unsigned int>(std::time(0)))
72  : adaptive_sampler(false),
73  _model(model),
74  _x(model.num_params_r()),
75  _z(model.num_params_i()),
76  _g(model.num_params_r()),
77 
78  _epsilon(epsilon),
79  _L(L),
80 
81  _rand_int(random_seed),
82  _rand_unit_norm(_rand_int,
83  boost::normal_distribution<>()),
84  _rand_uniform_01(_rand_int) {
85 
86  model.init(_x,_z);
87  _logp = model.grad_log_prob(_x,_z,_g);
88  }
89 
95  virtual ~hmc() {
96  }
97 
110  virtual void set_params(const std::vector<double>& x,
111  const std::vector<int>& z) {
112  if (x.size() != _x.size() || z.size() != _z.size())
113  throw std::invalid_argument("x.size() or z.size() mismatch");
114  _x = x;
115  _z = z;
116  }
117 
129  virtual void set_params_r(const std::vector<double>& x) {
130  if (x.size() != _model.num_params_r())
131  throw std::invalid_argument ("x.size() must match the number of parameters of the model.");
132  _x = x;
133  _logp = _model.grad_log_prob(_x,_z,_g);
134  }
135 
147  virtual void set_params_i(const std::vector<int>& z) {
148  if (z.size() != _model.num_params_i())
149  throw std::invalid_argument ("z.size() must match the number of parameters of the model.");
150  _z = z;
151  _logp = _model.grad_log_prob(_x,_z,_g);
152  }
153 
159  virtual sample next_impl() {
160  // Gibbs for discrete
161  std::vector<double> probs;
162  for (size_t m = 0; m < _model.num_params_i(); ++m) {
163  probs.resize(0);
164  for (int k = _model.param_range_i_lower(m);
165  k < _model.param_range_i_upper(m);
166  ++k)
167  probs.push_back(_model.log_prob_star(m,k,_x,_z));
168  _z[m] = sample_unnorm_log(probs,_rand_uniform_01);
169  }
170 
171  // HMC for continuous
172  std::vector<double> m(_model.num_params_r());
173  for (size_t i = 0; i < m.size(); ++i)
174  m[i] = _rand_unit_norm();
175  double H = -(stan::math::dot_self(m) / 2.0) + _logp;
176 
177  std::vector<double> g_new(_g);
178  std::vector<double> x_new(_x);
179  //double epsilon_over_2 = _epsilon / 2.0;
180 
181  double logp_new = -1e100;
182  for (unsigned int l = 0; l < _L; ++l)
183  logp_new = leapfrog(_model, _z, x_new, m, g_new, _epsilon);
184  nfevals_plus_eq(_L);
185 
186  double H_new = -(stan::math::dot_self(m) / 2.0) + logp_new;
187  double dH = H_new - H;
188  if (_rand_uniform_01() < exp(dH)) {
189  _x = x_new;
190  _g = g_new;
191  _logp = logp_new;
192  }
193  mcmc::sample s(_x,_z,_logp);
194  return s;
195  }
196 
197  };
198 
199  }
200 
201 }
202 
203 #endif
An abstract base class for adaptive samplers.
void nfevals_plus_eq(int n)
Add the specified number of evaluations to the number of function evaluations.
Hamiltonian Monte Carlo (HMC) sampler.
Definition: hmc.hpp:34
virtual void set_params_r(const std::vector< double > &x)
Set the model real parameters to the specified values and update gradients and log probability to mat...
Definition: hmc.hpp:129
virtual void set_params(const std::vector< double > &x, const std::vector< int > &z)
Set the model real and integer parameters to the specified values.
Definition: hmc.hpp:110
virtual ~hmc()
Destructor.
Definition: hmc.hpp:95
virtual sample next_impl()
Return the next sample.
Definition: hmc.hpp:159
hmc(stan::model::prob_grad &model, double epsilon, int L, unsigned int random_seed=static_cast< unsigned int >(std::time(0)))
Constructs a Hamiltonian Monte Carlo (HMC) sampler for the specified model, using the specified step ...
Definition: hmc.hpp:68
virtual void set_params_i(const std::vector< int > &z)
Set the model real parameters to the specified values and update gradients and log probability to mat...
Definition: hmc.hpp:147
Representation of a MCMC sample.
Definition: sampler.hpp:16
The prob_grad class represents densities with fixed numbers of discrete and scalar parameters and the...
Definition: prob_grad.hpp:24
virtual void init(std::vector< double > &params_r, std::vector< int > &params_i)
Definition: prob_grad.hpp:80
virtual size_t num_params_r()
Definition: prob_grad.hpp:52
virtual size_t num_params_i()
Definition: prob_grad.hpp:56
virtual double log_prob_star(size_t idx, int val, std::vector< double > &params_r, std::vector< int > &params_i, std::ostream *output_stream=0)
Definition: prob_grad.hpp:144
virtual double grad_log_prob(std::vector< double > &params_r, std::vector< int > &params_i, std::vector< double > &gradient, std::ostream *output_stream=0)=0
int param_range_i_lower(size_t idx)
Definition: prob_grad.hpp:72
int param_range_i_upper(size_t idx)
Definition: prob_grad.hpp:76
Reimplementing boost functionality.
var exp(const var &a)
Return the exponentiation of the specified variable (cmath).
Definition: agrad.hpp:1716
double epsilon()
Return minimum positive number representable.
double dot_self(const Eigen::Matrix< double, R, C > &v)
Returns the dot product of the specified vector with itself.
Definition: matrix.hpp:846
double leapfrog(stan::model::prob_grad &model, std::vector< int > z, std::vector< double > &x, std::vector< double > &m, std::vector< double > &g, double epsilon, std::ostream *error_msgs=0, std::ostream *output_msgs=0)
Computes the log probability for a single leapfrog step in Hamiltonian Monte Carlo.
Definition: util.hpp:57
int sample_unnorm_log(std::vector< double > probs, boost::uniform_01< boost::mt19937 & > &rand_uniform_01)
Definition: util.hpp:103
Probability, optimization and sampling library.
Definition: agrad.cpp:6

     [ Stan Home Page ] © 2011–2012, Stan Development Team.