Stan  2.5.0
probability, sampling & optimization
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Pages
base_hmc.hpp
Go to the documentation of this file.
1 #ifndef STAN__MCMC__BASE__HMC__BETA
2 #define STAN__MCMC__BASE__HMC__BETA
3 
4 #include <math.h>
5 #include <stdexcept>
6 
7 #include <boost/math/special_functions/fpclassify.hpp>
8 #include <boost/random/variate_generator.hpp>
9 #include <boost/random/uniform_01.hpp>
10 
11 #include <stan/mcmc/base_mcmc.hpp>
13 
14 namespace stan {
15 
16  namespace mcmc {
17 
18  template <class M, class P, template<class, class> class H,
19  template<class, class> class I, class BaseRNG>
20  class base_hmc: public base_mcmc {
21 
22  public:
23 
24  base_hmc(M &m, BaseRNG& rng, std::ostream* o, std::ostream* e):
25  base_mcmc(o, e),
26  z_(m.num_params_r()),
27  integrator_(this->out_stream_),
28  hamiltonian_(m, this->err_stream_),
29  rand_int_(rng),
31  nom_epsilon_(0.1),
33  epsilon_jitter_(0.0)
34  {};
35 
36  void write_sampler_state(std::ostream* o) {
37  if(!o) return;
38  *o << "# Step size = " << get_nominal_stepsize() << std::endl;
39  z_.write_metric(o);
40  }
41 
42  void get_sampler_diagnostic_names(std::vector<std::string>& model_names,
43  std::vector<std::string>& names) {
44  z_.get_param_names(model_names, names);
45  };
46 
47  void get_sampler_diagnostics(std::vector<double>& values) {
48  z_.get_params(values);
49  };
50 
51  void seed(const Eigen::VectorXd& q) {
52  z_.q = q;
53  }
54 
55  void init_stepsize() {
56 
57  ps_point z_init(this->z_);
58 
59  this->hamiltonian_.sample_p(this->z_, this->rand_int_);
60  this->hamiltonian_.init(this->z_);
61 
62  double H0 = this->hamiltonian_.H(this->z_); // Guaranteed to be finite if randomly initialized
63 
64  this->integrator_.evolve(this->z_, this->hamiltonian_, this->nom_epsilon_);
65 
66  double h = this->hamiltonian_.H(this->z_);
67  if (boost::math::isnan(h)) h = std::numeric_limits<double>::infinity();
68 
69  double delta_H = H0 - h;
70 
71  int direction = delta_H > std::log(0.8) ? 1 : -1;
72 
73  while (1) {
74 
75  this->z_.ps_point::operator=(z_init);
76 
77  this->hamiltonian_.sample_p(this->z_, this->rand_int_);
78  this->hamiltonian_.init(this->z_);
79 
80  double H0 = this->hamiltonian_.H(this->z_);
81 
82  this->integrator_.evolve(this->z_, this->hamiltonian_, this->nom_epsilon_);
83 
84  double h = this->hamiltonian_.H(this->z_);
85  if (boost::math::isnan(h)) h = std::numeric_limits<double>::infinity();
86 
87  double delta_H = H0 - h;
88 
89  if ((direction == 1) && !(delta_H > std::log(0.8)))
90  break;
91  else if ((direction == -1) && !(delta_H < std::log(0.8)))
92  break;
93  else
94  this->nom_epsilon_ = ( (direction == 1)
95  ? 2.0 * this->nom_epsilon_
96  : 0.5 * this->nom_epsilon_ );
97 
98  if (this->nom_epsilon_ > 1e7)
99  throw std::runtime_error("Posterior is improper. Please check your model.");
100  if (this->nom_epsilon_ == 0)
101  throw std::runtime_error("No acceptably small step size could be found. Perhaps the posterior is not continuous?");
102 
103  }
104 
105  this->z_.ps_point::operator=(z_init);
106 
107  }
108 
109  P& z() { return z_; }
110 
111  virtual void set_nominal_stepsize(const double e) {
112  if(e > 0) nom_epsilon_ = e;
113  }
114 
115  double get_nominal_stepsize() { return this->nom_epsilon_; }
116 
117  double get_current_stepsize() { return this->epsilon_; }
118 
119  virtual void set_stepsize_jitter(const double j) {
120  if(j > 0 && j < 1) epsilon_jitter_ = j;
121  }
122 
123  double get_stepsize_jitter() { return this->epsilon_jitter_; }
124 
126  this->epsilon_ = this->nom_epsilon_;
127  if(this->epsilon_jitter_)
128  this->epsilon_ *= ( 1.0 + this->epsilon_jitter_ * (2.0 * this->rand_uniform_() - 1.0) );
129  }
130 
131  protected:
132 
133  P z_;
134  I<H<M, BaseRNG>, P> integrator_;
135  H<M, BaseRNG> hamiltonian_;
136 
137  BaseRNG& rand_int_;
138 
139  // Uniform(0, 1) RNG
140  boost::uniform_01<BaseRNG&> rand_uniform_;
141 
142  double nom_epsilon_;
143  double epsilon_;
145 
146  };
147 
148  } // mcmc
149 
150 } // stan
151 
152 #endif
bool isnan(const stan::agrad::var &v)
Checks if the given number is NaN.
void write_sampler_state(std::ostream *o)
Definition: base_hmc.hpp:36
H< M, BaseRNG > hamiltonian_
Definition: base_hmc.hpp:135
virtual void set_stepsize_jitter(const double j)
Definition: base_hmc.hpp:119
std::ostream * out_stream_
Definition: base_mcmc.hpp:44
virtual void set_nominal_stepsize(const double e)
Definition: base_hmc.hpp:111
void get_sampler_diagnostics(std::vector< double > &values)
Definition: base_hmc.hpp:47
std::ostream * err_stream_
Definition: base_mcmc.hpp:45
double get_stepsize_jitter()
Definition: base_hmc.hpp:123
double e()
Return the base of the natural logarithm.
Definition: constants.hpp:86
double get_nominal_stepsize()
Definition: base_hmc.hpp:115
void get_sampler_diagnostic_names(std::vector< std::string > &model_names, std::vector< std::string > &names)
Definition: base_hmc.hpp:42
base_hmc(M &m, BaseRNG &rng, std::ostream *o, std::ostream *e)
Definition: base_hmc.hpp:24
void seed(const Eigen::VectorXd &q)
Definition: base_hmc.hpp:51
boost::uniform_01< BaseRNG & > rand_uniform_
Definition: base_hmc.hpp:140
fvar< T > log(const fvar< T > &x)
Definition: log.hpp:15
double get_current_stepsize()
Definition: base_hmc.hpp:117
I< H< M, BaseRNG >, P > integrator_
Definition: base_hmc.hpp:134

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