Stan  1.0
probability, sampling & optimization
adaptive_hmc.hpp
Go to the documentation of this file.
1 #ifndef __STAN__MCMC__ADAPTIVE_HMC_H__
2 #define __STAN__MCMC__ADAPTIVE_HMC_H__
3 
4 #include <ctime>
5 #include <cstddef>
6 #include <iostream>
7 #include <vector>
8 
9 #include <boost/random/mersenne_twister.hpp>
10 #include <boost/random/normal_distribution.hpp>
11 #include <boost/random/uniform_01.hpp>
12 #include <boost/random/variate_generator.hpp>
13 
14 #include <stan/math/util.hpp>
17 #include <stan/mcmc/hmc_base.hpp>
18 #include <stan/mcmc/util.hpp>
19 #include <stan/model/prob_grad.hpp>
20 
21 namespace stan {
22 
23  namespace mcmc {
24 
40  template <class BaseRNG = boost::mt19937>
41  class adaptive_hmc : public hmc_base<BaseRNG> {
42  private:
43 
44  unsigned int _L; // fixed number of Hamiltonian simulation steps
45 
46  public:
47 
78  int L,
79  double epsilon=-1,
80  double epsilon_pm = 0.0,
81  bool epsilon_adapt = true,
82  double delta = 0.651,
83  double gamma = 0.05,
84  BaseRNG rand_int = BaseRNG(std::time(0)),
85  const std::vector<double>* params_r = 0,
86  const std::vector<int>* params_i = 0)
87  : hmc_base<BaseRNG>(model,
88  epsilon,
89  epsilon_pm,
90  epsilon_adapt,
91  delta,
92  gamma,
93  rand_int,
94  params_r,
95  params_i),
96  _L(L) {
97  this->adaptation_init(1.0); // target is just epsilon
98  }
99 
100 
105  }
106 
107 
113  virtual sample next_impl() {
114  // Gibbs for discrete
115  // std::vector<double> probs;
116  // for (size_t m = 0; m < this->_model.num_params_i(); ++m) {
117  // probs.resize(0);
118  // for (int k = this->_model.param_range_i_lower(m);
119  // k < this->_model.param_range_i_upper(m);
120  // ++k)
121  // probs.push_back(this->_model.log_prob_star(m,k,this->_x,this->_z));
122  // this->_z[m] = sample_unnorm_log(probs,this->_rand_uniform_01);
123  // }
124  // HMC for continuous
125  std::vector<double> m(this->_model.num_params_r());
126  for (size_t i = 0; i < m.size(); ++i)
127  m[i] = this->_rand_unit_norm();
128  double H = -(stan::math::dot_self(m) / 2.0) + this->_logp;
129 
130  std::vector<double> g_new(this->_g);
131  std::vector<double> x_new(this->_x);
132  double logp_new = -1e100;
133  double epsilon = this->_epsilon;
134  // only vary epsilon after done adapting
135  if (!this->adapting() && this->varying_epsilon()) {
136  double low = epsilon * (1.0 - this->_epsilon_pm);
137  double high = epsilon * (1.0 + this->_epsilon_pm);
138  double range = high - low;
139  epsilon = low + (range * this->_rand_uniform_01());
140  }
141  this->_epsilon_last = epsilon;
142  for (unsigned int l = 0; l < _L; ++l)
143  logp_new = leapfrog(this->_model, this->_z, x_new, m, g_new, epsilon,
144  this->_error_msgs, this->_output_msgs);
145  this->nfevals_plus_eq(_L);
146 
147  double H_new = -(stan::math::dot_self(m) / 2.0) + logp_new;
148  double dH = H_new - H;
149  if (this->_rand_uniform_01() < exp(dH)) {
150  this->_x = x_new;
151  this->_g = g_new;
152  this->_logp = logp_new;
153  }
154 
155  // Now we just have to update epsilon, if adaptation is on.
156  double adapt_stat = stan::math::min(1, exp(dH));
157  if (adapt_stat != adapt_stat)
158  adapt_stat = 0;
159  if (this->adapting()) {
160  double adapt_g = adapt_stat - this->_delta;
161  std::vector<double> gvec(1, -adapt_g);
162  std::vector<double> result; // FIXME: update directly to _epsilon?
163  this->_da.update(gvec, result);
164  this->_epsilon = exp(result[0]);
165  }
166  std::vector<double> result;
167  this->_da.xbar(result);
168  // fprintf(stderr, "xbar = %f\n", exp(result[0]));
169  double avg_eta = 1.0 / this->n_steps();
170  this->update_mean_stat(avg_eta,adapt_stat);
171 
172  return mcmc::sample(this->_x, this->_z, this->_logp);
173  }
174 
175  virtual void write_sampler_param_names(std::ostream& o) {
176  if (this->_epsilon_adapt || this->varying_epsilon())
177  o << "stepsize__,";
178  }
179 
180  virtual void write_sampler_params(std::ostream& o) {
181  if (this->_epsilon_adapt || this->varying_epsilon())
182  o << this->_epsilon_last << ',';
183  }
184 
185  virtual void get_sampler_param_names(std::vector<std::string>& names) {
186  names.clear();
187  if (this->_epsilon_adapt || this->varying_epsilon())
188  names.push_back("stepsize__");
189  }
190 
191  virtual void get_sampler_params(std::vector<double>& values) {
192  values.clear();
193  if (this->_epsilon_adapt || this->varying_epsilon())
194  values.push_back(this->_epsilon_last);
195  }
196  };
197 
198  }
199 
200 }
201 
202 #endif
void update(const std::vector< double > &g, std::vector< double > &xk)
Produces the next iterate xk given the current gradient g.
Definition: dualaverage.hpp:53
void xbar(std::vector< double > &xbar)
Get the exponentially weighted moving average of all previous iterates.
Definition: dualaverage.hpp:86
Adaptive Hamiltonian Monte Carlo (HMC) sampler.
adaptive_hmc(stan::model::prob_grad &model, int L, double epsilon=-1, double epsilon_pm=0.0, bool epsilon_adapt=true, double delta=0.651, double gamma=0.05, BaseRNG rand_int=BaseRNG(std::time(0)), const std::vector< double > *params_r=0, const std::vector< int > *params_i=0)
Construct an adaptive Hamiltonian Monte Carlo (HMC) sampler for the specified model,...
virtual void get_sampler_param_names(std::vector< std::string > &names)
Get any sampler-specific parameter namess.
virtual void write_sampler_param_names(std::ostream &o)
Write out any sampler-specific parameter names for output.
virtual void write_sampler_params(std::ostream &o)
Write out any sampler-specific parameters for output.
virtual sample next_impl()
Returns the next sample.
virtual void get_sampler_params(std::vector< double > &values)
Get any sampler-specific parameters.
bool adapting()
Return whether or not parameter adaptation is on.
int n_steps()
Return the number of iterations for this sampler.
void update_mean_stat(double avg_eta, double adapt_stat)
Updates the mean statistic given the specified adaptation statistic and weighting.
void nfevals_plus_eq(int n)
Add the specified number of evaluations to the number of function evaluations.
boost::uniform_01< boost::mt19937 & > _rand_uniform_01
Definition: hmc_base.hpp:44
stan::model::prob_grad & _model
Definition: hmc_base.hpp:27
void adaptation_init(double epsilon_scale)
Definition: hmc_base.hpp:90
boost::variate_generator< boost::mt19937 &, boost::normal_distribution<> > _rand_unit_norm
Definition: hmc_base.hpp:41
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 size_t num_params_r()
Definition: prob_grad.hpp:52
var exp(const var &a)
Return the exponentiation of the specified variable (cmath).
Definition: agrad.hpp:1716
double epsilon()
Return minimum positive number representable.
int min(const std::vector< int > &x)
Returns the minimum coefficient in the specified column vector.
Definition: matrix.hpp:917
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
Probability, optimization and sampling library.
Definition: agrad.cpp:6

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