Stan  2.10.0
probability, sampling & optimization
adapt_diag_e_static_hmc.hpp
Go to the documentation of this file.
1 #ifndef STAN_MCMC_HMC_STATIC_ADAPT_DIAG_E_STATIC_HMC_HPP
2 #define STAN_MCMC_HMC_STATIC_ADAPT_DIAG_E_STATIC_HMC_HPP
3 
7 
8 namespace stan {
9  namespace mcmc {
16  template <class Model, class BaseRNG>
17  class adapt_diag_e_static_hmc : public diag_e_static_hmc<Model, BaseRNG>,
18  public stepsize_var_adapter {
19  public:
20  adapt_diag_e_static_hmc(const Model& model, BaseRNG& rng)
21  : diag_e_static_hmc<Model, BaseRNG>(model, rng),
22  stepsize_var_adapter(model.num_params_r()) {}
23 
25 
26  sample
27  transition(sample& init_sample,
30  sample s
32  info_writer,
33  error_writer);
34 
35  if (this->adapt_flag_) {
37  s.accept_stat());
38  this->update_L_();
39 
40  bool update = this->var_adaptation_.learn_variance(this->z_.mInv,
41  this->z_.q);
42 
43  if (update) {
44  this->init_stepsize(info_writer, error_writer);
45  this->update_L_();
46 
47  this->stepsize_adaptation_.set_mu(log(10 * this->nom_epsilon_));
49  }
50  }
51  return s;
52  }
53 
57  }
58  };
59 
60  } // mcmc
61 } // stan
62 #endif
void complete_adaptation(double &epsilon)
double accept_stat() const
Definition: sample.hpp:41
Probability, optimization and sampling library.
void learn_stepsize(double &epsilon, double adapt_stat)
sample transition(sample &init_sample, interface_callbacks::writer::base_writer &info_writer, interface_callbacks::writer::base_writer &error_writer)
base_writer is an abstract base class defining the interface for Stan writer callbacks.
Definition: base_writer.hpp:20
virtual void disengage_adaptation()
Hamiltonian Monte Carlo implementation using the endpoint of trajectories with a static integration t...
adapt_diag_e_static_hmc(const Model &model, BaseRNG &rng)
sample transition(sample &init_sample, interface_callbacks::writer::base_writer &info_writer, interface_callbacks::writer::base_writer &error_writer)
Hamiltonian Monte Carlo implementation using the endpoint of trajectories with a static integration t...
bool learn_variance(Eigen::VectorXd &var, const Eigen::VectorXd &q)
void init_stepsize(interface_callbacks::writer::base_writer &info_writer, interface_callbacks::writer::base_writer &error_writer)
Definition: base_hmc.hpp:64

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