Stan  2.5.0
probability, sampling & optimization
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Pages
base_static_hmc.hpp
Go to the documentation of this file.
1 #ifndef STAN__MCMC__BASE__STATIC__HMC__BETA
2 #define STAN__MCMC__BASE__STATIC__HMC__BETA
3 
4 #include <math.h>
5 #include <boost/math/special_functions/fpclassify.hpp>
8 
9 namespace stan {
10 
11  namespace mcmc {
12 
13  // Hamiltonian Monte Carlo
14  // with static integration time
15 
16  template <class M, class P, template<class, class> class H,
17  template<class, class> class I, class BaseRNG>
18  class base_static_hmc: public base_hmc<M, P, H, I, BaseRNG> {
19 
20  public:
21 
22  base_static_hmc(M &m, BaseRNG& rng, std::ostream* o, std::ostream* e):
23  base_hmc<M, P, H, I, BaseRNG>(m, rng, o, e), T_(1)
24  { update_L_(); }
25 
27 
28  sample transition(sample& init_sample) {
29 
30  this->sample_stepsize();
31 
32  this->seed(init_sample.cont_params());
33 
34  this->hamiltonian_.sample_p(this->z_, this->rand_int_);
35  this->hamiltonian_.init(this->z_);
36 
37  ps_point z_init(this->z_);
38 
39  double H0 = this->hamiltonian_.H(this->z_);
40 
41  for (int i = 0; i < L_; ++i)
42  this->integrator_.evolve(this->z_, this->hamiltonian_, this->epsilon_);
43 
44  double h = this->hamiltonian_.H(this->z_);
45  if (boost::math::isnan(h)) h = std::numeric_limits<double>::infinity();
46 
47  double acceptProb = std::exp(H0 - h);
48 
49  if (acceptProb < 1 && this->rand_uniform_() > acceptProb)
50  this->z_.ps_point::operator=(z_init);
51 
52  acceptProb = acceptProb > 1 ? 1 : acceptProb;
53 
54  return sample(this->z_.q, - this->hamiltonian_.V(this->z_), acceptProb);
55 
56  }
57 
58  void write_sampler_param_names(std::ostream& o) {
59  o << "stepsize__,int_time__,";
60  }
61 
62  void write_sampler_params(std::ostream& o) {
63  o << this->epsilon_ << "," << this->T_ << ",";
64  }
65 
66  void get_sampler_param_names(std::vector<std::string>& names) {
67  names.push_back("stepsize__");
68  names.push_back("int_time__");
69  }
70 
71  void get_sampler_params(std::vector<double>& values) {
72  values.push_back(this->epsilon_);
73  values.push_back(this->T_);
74  }
75 
76  void set_nominal_stepsize_and_T(const double e, const double t) {
77  if(e > 0 && t > 0) {
78  this->nom_epsilon_ = e;
79  T_ = t;
80  update_L_();
81  }
82  }
83 
84  void set_nominal_stepsize_and_L(const double e, const int l) {
85  if(e > 0 && l > 0) {
86  this->nom_epsilon_ = e;
87  L_ = l;
88  T_ = this->nom_epsilon_ * L_; }
89  }
90 
91  void set_T(const double t) {
92  if(t > 0) {
93  T_ = t;
94  update_L_();
95  }
96  }
97 
98  void set_nominal_stepsize(const double e) {
99  if(e > 0) {
100  this->nom_epsilon_ = e;
101  update_L_();
102  }
103  }
104 
105  double get_T() { return this->T_; }
106  int get_L() { return this->L_; }
107 
108  protected:
109 
110  double T_;
111  int L_;
112 
113  void update_L_() {
114  L_ = static_cast<int>(T_ / this->nom_epsilon_);
115  L_ = L_ < 1 ? 1 : L_;
116  }
117 
118  };
119 
120  } // mcmc
121 
122 } // stan
123 
124 
125 #endif
void write_sampler_param_names(std::ostream &o)
bool isnan(const stan::agrad::var &v)
Checks if the given number is NaN.
void write_sampler_params(std::ostream &o)
void set_nominal_stepsize_and_T(const double e, const double t)
H< M, BaseRNG > hamiltonian_
Definition: base_hmc.hpp:135
void set_nominal_stepsize_and_L(const double e, const int l)
void set_nominal_stepsize(const double e)
double cont_params(int k) const
Definition: sample.hpp:30
void set_T(const double t)
base_static_hmc(M &m, BaseRNG &rng, std::ostream *o, std::ostream *e)
double e()
Return the base of the natural logarithm.
Definition: constants.hpp:86
void seed(const Eigen::VectorXd &q)
Definition: base_hmc.hpp:51
boost::uniform_01< BaseRNG & > rand_uniform_
Definition: base_hmc.hpp:140
void get_sampler_param_names(std::vector< std::string > &names)
sample transition(sample &init_sample)
fvar< T > exp(const fvar< T > &x)
Definition: exp.hpp:16
void sample(stan::mcmc::base_mcmc *sampler, int num_warmup, int num_samples, int num_thin, int refresh, bool save, stan::io::mcmc_writer< Model, SampleRecorder, DiagnosticRecorder, MessageRecorder > &writer, stan::mcmc::sample &init_s, Model &model, RNG &base_rng, const std::string &prefix, const std::string &suffix, std::ostream &o, StartTransitionCallback &callback)
Definition: sample.hpp:13
void get_sampler_params(std::vector< double > &values)
I< H< M, BaseRNG >, P > integrator_
Definition: base_hmc.hpp:134

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