Public Member Functions | Protected Attributes | List of all members
stan::mcmc::base_hmc< Model, Hamiltonian, Integrator, BaseRNG > Class Template Referenceabstract

#include "stan/mcmc/hmc/base_hmc.hpp"

Inheritance diagram for stan::mcmc::base_hmc< Model, Hamiltonian, Integrator, BaseRNG >:
stan::mcmc::base_mcmc stan::mcmc::base_nuts< Model, Hamiltonian, Integrator, BaseRNG > stan::mcmc::base_nuts_classic< Model, Hamiltonian, Integrator, BaseRNG > stan::mcmc::base_static_hmc< Model, Hamiltonian, Integrator, BaseRNG > stan::mcmc::base_static_uniform< Model, Hamiltonian, Integrator, BaseRNG > stan::mcmc::base_xhmc< Model, Hamiltonian, Integrator, BaseRNG >

Public Member Functions

 base_hmc (const Model &model, BaseRNG &rng)
 
void write_sampler_stepsize (callbacks::writer &writer)
 
void write_sampler_metric (callbacks::writer &writer)
 
void write_sampler_state (callbacks::writer &writer)
 
void get_sampler_diagnostic_names (std::vector< std::string > &model_names, std::vector< std::string > &names)
 
void get_sampler_diagnostics (std::vector< double > &values)
 
void seed (const Eigen::VectorXd &q)
 
void init_hamiltonian (callbacks::logger &logger)
 
void init_stepsize (callbacks::logger &logger)
 
Hamiltonian< Model, BaseRNG >::PointType & z ()
 
virtual void set_nominal_stepsize (double e)
 
double get_nominal_stepsize ()
 
double get_current_stepsize ()
 
virtual void set_stepsize_jitter (double j)
 
double get_stepsize_jitter ()
 
void sample_stepsize ()
 
virtual sample transition (sample &init_sample, callbacks::logger &logger)=0
 
virtual void get_sampler_param_names (std::vector< std::string > &names)
 
virtual void get_sampler_params (std::vector< double > &values)
 

Protected Attributes

Hamiltonian< Model, BaseRNG >::PointType z_
 
Integrator< Hamiltonian< Model, BaseRNG > > integrator_
 
Hamiltonian< Model, BaseRNG > hamiltonian_
 
BaseRNG & rand_int_
 
boost::uniform_01< BaseRNG & > rand_uniform_
 
double nom_epsilon_
 
double epsilon_
 
double epsilon_jitter_
 

Detailed Description

template<class Model, template< class, class > class Hamiltonian, template< class > class Integrator, class BaseRNG>
class stan::mcmc::base_hmc< Model, Hamiltonian, Integrator, BaseRNG >

Definition at line 24 of file base_hmc.hpp.

Constructor & Destructor Documentation

template<class Model, template< class, class > class Hamiltonian, template< class > class Integrator, class BaseRNG>
stan::mcmc::base_hmc< Model, Hamiltonian, Integrator, BaseRNG >::base_hmc ( const Model model,
BaseRNG &  rng 
)
inline

Definition at line 26 of file base_hmc.hpp.

27  : base_mcmc(),
28  z_(model.num_params_r()),
29  integrator_(),
31  rand_int_(rng),
33  nom_epsilon_(0.1),
35  epsilon_jitter_(0.0) {}
Hamiltonian< Model, BaseRNG >::PointType z_
Definition: base_hmc.hpp:182
boost::uniform_01< BaseRNG & > rand_uniform_
Definition: base_hmc.hpp:189
Hamiltonian< Model, BaseRNG > hamiltonian_
Definition: base_hmc.hpp:184
const XML_Char XML_Content * model
Definition: expat.h:151
Integrator< Hamiltonian< Model, BaseRNG > > integrator_
Definition: base_hmc.hpp:183

Member Function Documentation

template<class Model, template< class, class > class Hamiltonian, template< class > class Integrator, class BaseRNG>
double stan::mcmc::base_hmc< Model, Hamiltonian, Integrator, BaseRNG >::get_current_stepsize ( )
inline

Definition at line 161 of file base_hmc.hpp.

161  {
162  return this->epsilon_;
163  }
template<class Model, template< class, class > class Hamiltonian, template< class > class Integrator, class BaseRNG>
double stan::mcmc::base_hmc< Model, Hamiltonian, Integrator, BaseRNG >::get_nominal_stepsize ( )
inline

Definition at line 157 of file base_hmc.hpp.

Referenced by stan::mcmc::base_hmc< Model, diag_e_metric, expl_leapfrog, BaseRNG >::write_sampler_stepsize().

157  {
158  return this->nom_epsilon_;
159  }
template<class Model, template< class, class > class Hamiltonian, template< class > class Integrator, class BaseRNG>
void stan::mcmc::base_hmc< Model, Hamiltonian, Integrator, BaseRNG >::get_sampler_diagnostic_names ( std::vector< std::string > &  model_names,
std::vector< std::string > &  names 
)
inlinevirtual

Reimplemented from stan::mcmc::base_mcmc.

Definition at line 64 of file base_hmc.hpp.

65  {
66  z_.get_param_names(model_names, names);
67  }
Hamiltonian< Model, BaseRNG >::PointType z_
Definition: base_hmc.hpp:182
template<class Model, template< class, class > class Hamiltonian, template< class > class Integrator, class BaseRNG>
void stan::mcmc::base_hmc< Model, Hamiltonian, Integrator, BaseRNG >::get_sampler_diagnostics ( std::vector< double > &  values)
inlinevirtual

Reimplemented from stan::mcmc::base_mcmc.

Definition at line 69 of file base_hmc.hpp.

69  {
70  z_.get_params(values);
71  }
Hamiltonian< Model, BaseRNG >::PointType z_
Definition: base_hmc.hpp:182
virtual void stan::mcmc::base_mcmc::get_sampler_param_names ( std::vector< std::string > &  names)
inlinevirtualinherited

Reimplemented in stan::mcmc::base_xhmc< Model, Hamiltonian, Integrator, BaseRNG >, stan::mcmc::base_xhmc< Model, softabs_metric, impl_leapfrog, BaseRNG >, stan::mcmc::base_xhmc< mock_model, mock_hamiltonian, mock_integrator, rng_t >, stan::mcmc::base_xhmc< mock_model, divergent_hamiltonian, expl_leapfrog, rng_t >, stan::mcmc::base_xhmc< Model, unit_e_metric, expl_leapfrog, BaseRNG >, stan::mcmc::base_xhmc< Model, dense_e_metric, expl_leapfrog, BaseRNG >, stan::mcmc::base_xhmc< Model, diag_e_metric, expl_leapfrog, BaseRNG >, stan::mcmc::base_nuts< Model, Hamiltonian, Integrator, BaseRNG >, stan::mcmc::base_nuts< Model, softabs_metric, impl_leapfrog, BaseRNG >, stan::mcmc::base_nuts< stan_model, unit_e_metric, expl_leapfrog, boost::ecuyer1988 >, stan::mcmc::base_nuts< mock_model, mock_hamiltonian, mock_integrator, rng_t >, stan::mcmc::base_nuts< mock_model, divergent_hamiltonian, expl_leapfrog, rng_t >, stan::mcmc::base_nuts< Model, unit_e_metric, expl_leapfrog, BaseRNG >, stan::mcmc::base_nuts< Model, dense_e_metric, expl_leapfrog, BaseRNG >, stan::mcmc::base_nuts< Model, diag_e_metric, expl_leapfrog, BaseRNG >, stan::mcmc::base_nuts_classic< Model, Hamiltonian, Integrator, BaseRNG >, stan::mcmc::base_nuts_classic< mock_model, mock_hamiltonian, mock_integrator, rng_t >, stan::mcmc::base_nuts_classic< mock_model, divergent_hamiltonian, expl_leapfrog, rng_t >, stan::mcmc::base_nuts_classic< Model, unit_e_metric, expl_leapfrog, BaseRNG >, stan::mcmc::base_nuts_classic< Model, dense_e_metric, expl_leapfrog, BaseRNG >, stan::mcmc::base_nuts_classic< Model, diag_e_metric, expl_leapfrog, BaseRNG >, stan::mcmc::base_static_uniform< Model, Hamiltonian, Integrator, BaseRNG >, stan::mcmc::base_static_uniform< Model, softabs_metric, impl_leapfrog, BaseRNG >, stan::mcmc::base_static_uniform< mock_model, mock_hamiltonian, mock_integrator, rng_t >, stan::mcmc::base_static_uniform< Model, unit_e_metric, expl_leapfrog, BaseRNG >, stan::mcmc::base_static_uniform< Model, dense_e_metric, expl_leapfrog, BaseRNG >, stan::mcmc::base_static_uniform< Model, diag_e_metric, expl_leapfrog, BaseRNG >, mock_sampler, stan::mcmc::base_static_hmc< Model, Hamiltonian, Integrator, BaseRNG >, stan::mcmc::base_static_hmc< Model, softabs_metric, impl_leapfrog, BaseRNG >, stan::mcmc::base_static_hmc< mock_model, mock_hamiltonian, mock_integrator, rng_t >, stan::mcmc::base_static_hmc< Model, unit_e_metric, expl_leapfrog, BaseRNG >, stan::mcmc::base_static_hmc< Model, dense_e_metric, expl_leapfrog, BaseRNG >, stan::mcmc::base_static_hmc< Model, diag_e_metric, expl_leapfrog, BaseRNG >, mock_sampler, and stan::mcmc::mock_hmc.

Definition at line 23 of file base_mcmc.hpp.

Referenced by stan::services::util::mcmc_writer::write_diagnostic_names(), and stan::services::util::mcmc_writer::write_sample_names().

23 {}
virtual void stan::mcmc::base_mcmc::get_sampler_params ( std::vector< double > &  values)
inlinevirtualinherited

Reimplemented in stan::mcmc::base_xhmc< Model, Hamiltonian, Integrator, BaseRNG >, stan::mcmc::base_xhmc< Model, softabs_metric, impl_leapfrog, BaseRNG >, stan::mcmc::base_xhmc< mock_model, mock_hamiltonian, mock_integrator, rng_t >, stan::mcmc::base_xhmc< mock_model, divergent_hamiltonian, expl_leapfrog, rng_t >, stan::mcmc::base_xhmc< Model, unit_e_metric, expl_leapfrog, BaseRNG >, stan::mcmc::base_xhmc< Model, dense_e_metric, expl_leapfrog, BaseRNG >, stan::mcmc::base_xhmc< Model, diag_e_metric, expl_leapfrog, BaseRNG >, stan::mcmc::base_nuts< Model, Hamiltonian, Integrator, BaseRNG >, stan::mcmc::base_nuts< Model, softabs_metric, impl_leapfrog, BaseRNG >, stan::mcmc::base_nuts< stan_model, unit_e_metric, expl_leapfrog, boost::ecuyer1988 >, stan::mcmc::base_nuts< mock_model, mock_hamiltonian, mock_integrator, rng_t >, stan::mcmc::base_nuts< mock_model, divergent_hamiltonian, expl_leapfrog, rng_t >, stan::mcmc::base_nuts< Model, unit_e_metric, expl_leapfrog, BaseRNG >, stan::mcmc::base_nuts< Model, dense_e_metric, expl_leapfrog, BaseRNG >, stan::mcmc::base_nuts< Model, diag_e_metric, expl_leapfrog, BaseRNG >, stan::mcmc::base_nuts_classic< Model, Hamiltonian, Integrator, BaseRNG >, stan::mcmc::base_nuts_classic< mock_model, mock_hamiltonian, mock_integrator, rng_t >, stan::mcmc::base_nuts_classic< mock_model, divergent_hamiltonian, expl_leapfrog, rng_t >, stan::mcmc::base_nuts_classic< Model, unit_e_metric, expl_leapfrog, BaseRNG >, stan::mcmc::base_nuts_classic< Model, dense_e_metric, expl_leapfrog, BaseRNG >, stan::mcmc::base_nuts_classic< Model, diag_e_metric, expl_leapfrog, BaseRNG >, stan::mcmc::base_static_uniform< Model, Hamiltonian, Integrator, BaseRNG >, stan::mcmc::base_static_uniform< Model, softabs_metric, impl_leapfrog, BaseRNG >, stan::mcmc::base_static_uniform< mock_model, mock_hamiltonian, mock_integrator, rng_t >, stan::mcmc::base_static_uniform< Model, unit_e_metric, expl_leapfrog, BaseRNG >, stan::mcmc::base_static_uniform< Model, dense_e_metric, expl_leapfrog, BaseRNG >, stan::mcmc::base_static_uniform< Model, diag_e_metric, expl_leapfrog, BaseRNG >, mock_sampler, stan::mcmc::base_static_hmc< Model, Hamiltonian, Integrator, BaseRNG >, stan::mcmc::base_static_hmc< Model, softabs_metric, impl_leapfrog, BaseRNG >, stan::mcmc::base_static_hmc< mock_model, mock_hamiltonian, mock_integrator, rng_t >, stan::mcmc::base_static_hmc< Model, unit_e_metric, expl_leapfrog, BaseRNG >, stan::mcmc::base_static_hmc< Model, dense_e_metric, expl_leapfrog, BaseRNG >, stan::mcmc::base_static_hmc< Model, diag_e_metric, expl_leapfrog, BaseRNG >, mock_sampler, and stan::mcmc::mock_hmc.

Definition at line 25 of file base_mcmc.hpp.

Referenced by stan::services::util::mcmc_writer::write_diagnostic_params(), and stan::services::util::mcmc_writer::write_sample_params().

25 {}
template<class Model, template< class, class > class Hamiltonian, template< class > class Integrator, class BaseRNG>
double stan::mcmc::base_hmc< Model, Hamiltonian, Integrator, BaseRNG >::get_stepsize_jitter ( )
inline

Definition at line 170 of file base_hmc.hpp.

170  {
171  return this->epsilon_jitter_;
172  }
template<class Model, template< class, class > class Hamiltonian, template< class > class Integrator, class BaseRNG>
void stan::mcmc::base_hmc< Model, Hamiltonian, Integrator, BaseRNG >::init_hamiltonian ( callbacks::logger logger)
inline

Definition at line 78 of file base_hmc.hpp.

78  {
79  this->hamiltonian_.init(this->z_, logger);
80  }
Hamiltonian< Model, BaseRNG >::PointType z_
Definition: base_hmc.hpp:182
Hamiltonian< Model, BaseRNG > hamiltonian_
Definition: base_hmc.hpp:184
template<class Model, template< class, class > class Hamiltonian, template< class > class Integrator, class BaseRNG>
void stan::mcmc::base_hmc< Model, Hamiltonian, Integrator, BaseRNG >::init_stepsize ( callbacks::logger logger)
inline

Definition at line 83 of file base_hmc.hpp.

83  {
84  ps_point z_init(this->z_);
85 
86  // Skip initialization for extreme step sizes
87  if (this->nom_epsilon_ == 0 || this->nom_epsilon_ > 1e7)
88  return;
89 
90  this->hamiltonian_.sample_p(this->z_, this->rand_int_);
91  this->hamiltonian_.init(this->z_, logger);
92 
93  // Guaranteed to be finite if randomly initialized
94  double H0 = this->hamiltonian_.H(this->z_);
95 
96  this->integrator_.evolve(this->z_, this->hamiltonian_,
97  this->nom_epsilon_,
98  logger);
99 
100  double h = this->hamiltonian_.H(this->z_);
101  if (boost::math::isnan(h))
102  h = std::numeric_limits<double>::infinity();
103 
104  double delta_H = H0 - h;
105 
106  int direction = delta_H > std::log(0.8) ? 1 : -1;
107 
108  while (1) {
109  this->z_.ps_point::operator=(z_init);
110 
111  this->hamiltonian_.sample_p(this->z_, this->rand_int_);
112  this->hamiltonian_.init(this->z_, logger);
113 
114  double H0 = this->hamiltonian_.H(this->z_);
115 
116  this->integrator_.evolve(this->z_, this->hamiltonian_,
117  this->nom_epsilon_,
118  logger);
119 
120  double h = this->hamiltonian_.H(this->z_);
121  if (boost::math::isnan(h))
122  h = std::numeric_limits<double>::infinity();
123 
124  double delta_H = H0 - h;
125 
126  if ((direction == 1) && !(delta_H > std::log(0.8)))
127  break;
128  else if ((direction == -1) && !(delta_H < std::log(0.8)))
129  break;
130  else
131  this->nom_epsilon_
132  = direction == 1
133  ? 2.0 * this->nom_epsilon_
134  : 0.5 * this->nom_epsilon_;
135 
136  if (this->nom_epsilon_ > 1e7)
137  throw std::runtime_error("Posterior is improper. "
138  "Please check your model.");
139  if (this->nom_epsilon_ == 0)
140  throw std::runtime_error("No acceptably small step size could "
141  "be found. Perhaps the posterior is "
142  "not continuous?");
143  }
144 
145  this->z_.ps_point::operator=(z_init);
146  }
Hamiltonian< Model, BaseRNG >::PointType z_
Definition: base_hmc.hpp:182
bool isnan(const stan::math::var &v)
Definition: boost_isnan.hpp:20
Hamiltonian< Model, BaseRNG > hamiltonian_
Definition: base_hmc.hpp:184
Integrator< Hamiltonian< Model, BaseRNG > > integrator_
Definition: base_hmc.hpp:183
template<class Model, template< class, class > class Hamiltonian, template< class > class Integrator, class BaseRNG>
void stan::mcmc::base_hmc< Model, Hamiltonian, Integrator, BaseRNG >::sample_stepsize ( )
inline
template<class Model, template< class, class > class Hamiltonian, template< class > class Integrator, class BaseRNG>
void stan::mcmc::base_hmc< Model, Hamiltonian, Integrator, BaseRNG >::seed ( const Eigen::VectorXd &  q)
inline
template<class Model, template< class, class > class Hamiltonian, template< class > class Integrator, class BaseRNG>
virtual void stan::mcmc::base_hmc< Model, Hamiltonian, Integrator, BaseRNG >::set_nominal_stepsize ( double  e)
inlinevirtual
template<class Model, template< class, class > class Hamiltonian, template< class > class Integrator, class BaseRNG>
virtual void stan::mcmc::base_hmc< Model, Hamiltonian, Integrator, BaseRNG >::set_stepsize_jitter ( double  j)
inlinevirtual

Definition at line 165 of file base_hmc.hpp.

Referenced by TEST().

165  {
166  if (j > 0 && j < 1)
167  epsilon_jitter_ = j;
168  }
const double j
Definition: BetheBloch.cxx:29
virtual sample stan::mcmc::base_mcmc::transition ( sample init_sample,
callbacks::logger logger 
)
pure virtualinherited

Implemented in mock_sampler, stan::mcmc::base_xhmc< Model, Hamiltonian, Integrator, BaseRNG >, stan::mcmc::base_xhmc< Model, softabs_metric, impl_leapfrog, BaseRNG >, stan::mcmc::base_xhmc< mock_model, mock_hamiltonian, mock_integrator, rng_t >, stan::mcmc::base_xhmc< mock_model, divergent_hamiltonian, expl_leapfrog, rng_t >, stan::mcmc::base_xhmc< Model, unit_e_metric, expl_leapfrog, BaseRNG >, stan::mcmc::base_xhmc< Model, dense_e_metric, expl_leapfrog, BaseRNG >, stan::mcmc::base_xhmc< Model, diag_e_metric, expl_leapfrog, BaseRNG >, stan::mcmc::base_nuts< Model, Hamiltonian, Integrator, BaseRNG >, stan::mcmc::base_nuts< Model, softabs_metric, impl_leapfrog, BaseRNG >, stan::mcmc::base_nuts< stan_model, unit_e_metric, expl_leapfrog, boost::ecuyer1988 >, stan::mcmc::base_nuts< mock_model, mock_hamiltonian, mock_integrator, rng_t >, stan::mcmc::base_nuts< mock_model, divergent_hamiltonian, expl_leapfrog, rng_t >, stan::mcmc::base_nuts< Model, unit_e_metric, expl_leapfrog, BaseRNG >, stan::mcmc::base_nuts< Model, dense_e_metric, expl_leapfrog, BaseRNG >, stan::mcmc::base_nuts< Model, diag_e_metric, expl_leapfrog, BaseRNG >, stan::mcmc::base_nuts_classic< Model, Hamiltonian, Integrator, BaseRNG >, stan::mcmc::base_nuts_classic< mock_model, mock_hamiltonian, mock_integrator, rng_t >, stan::mcmc::base_nuts_classic< mock_model, divergent_hamiltonian, expl_leapfrog, rng_t >, stan::mcmc::base_nuts_classic< Model, unit_e_metric, expl_leapfrog, BaseRNG >, stan::mcmc::base_nuts_classic< Model, dense_e_metric, expl_leapfrog, BaseRNG >, stan::mcmc::base_nuts_classic< Model, diag_e_metric, expl_leapfrog, BaseRNG >, stan::mcmc::base_static_hmc< Model, Hamiltonian, Integrator, BaseRNG >, stan::mcmc::base_static_hmc< Model, softabs_metric, impl_leapfrog, BaseRNG >, stan::mcmc::base_static_hmc< mock_model, mock_hamiltonian, mock_integrator, rng_t >, stan::mcmc::base_static_hmc< Model, unit_e_metric, expl_leapfrog, BaseRNG >, stan::mcmc::base_static_hmc< Model, dense_e_metric, expl_leapfrog, BaseRNG >, stan::mcmc::base_static_hmc< Model, diag_e_metric, expl_leapfrog, BaseRNG >, stan::mcmc::base_static_uniform< Model, Hamiltonian, Integrator, BaseRNG >, stan::mcmc::base_static_uniform< Model, softabs_metric, impl_leapfrog, BaseRNG >, stan::mcmc::base_static_uniform< mock_model, mock_hamiltonian, mock_integrator, rng_t >, stan::mcmc::base_static_uniform< Model, unit_e_metric, expl_leapfrog, BaseRNG >, stan::mcmc::base_static_uniform< Model, dense_e_metric, expl_leapfrog, BaseRNG >, stan::mcmc::base_static_uniform< Model, diag_e_metric, expl_leapfrog, BaseRNG >, mock_sampler, stan::mcmc::adapt_diag_e_nuts_classic< Model, BaseRNG >, stan::mcmc::adapt_dense_e_static_hmc< Model, BaseRNG >, stan::mcmc::adapt_diag_e_static_hmc< Model, BaseRNG >, stan::mcmc::adapt_dense_e_static_uniform< Model, BaseRNG >, stan::mcmc::adapt_diag_e_static_uniform< Model, BaseRNG >, stan::mcmc::adapt_unit_e_static_uniform< Model, BaseRNG >, stan::mcmc::adapt_dense_e_nuts< Model, BaseRNG >, stan::mcmc::adapt_diag_e_nuts< Model, BaseRNG >, stan::mcmc::adapt_dense_e_nuts_classic< Model, BaseRNG >, stan::mcmc::adapt_unit_e_nuts_classic< Model, BaseRNG >, stan::mcmc::adapt_unit_e_static_hmc< Model, BaseRNG >, stan::mcmc::adapt_softabs_static_uniform< Model, BaseRNG >, stan::mcmc::adapt_dense_e_xhmc< Model, BaseRNG >, stan::mcmc::adapt_diag_e_xhmc< Model, BaseRNG >, stan::mcmc::adapt_unit_e_nuts< Model, BaseRNG >, stan::mcmc::adapt_softabs_static_hmc< Model, BaseRNG >, stan::mcmc::adapt_unit_e_xhmc< Model, BaseRNG >, stan::mcmc::mock_hmc, stan::mcmc::adapt_unit_e_nuts< stan_model, boost::ecuyer1988 >, stan::mcmc::adapt_softabs_nuts< Model, BaseRNG >, stan::mcmc::adapt_softabs_xhmc< Model, BaseRNG >, and stan::mcmc::fixed_param_sampler.

Referenced by stan::services::util::generate_transitions(), and stan::mcmc::base_mcmc::~base_mcmc().

template<class Model, template< class, class > class Hamiltonian, template< class > class Integrator, class BaseRNG>
void stan::mcmc::base_hmc< Model, Hamiltonian, Integrator, BaseRNG >::write_sampler_metric ( callbacks::writer writer)
inline

write elements of mass matrix

Definition at line 51 of file base_hmc.hpp.

Referenced by stan::mcmc::base_hmc< Model, diag_e_metric, expl_leapfrog, BaseRNG >::write_sampler_state().

51  {
52  z_.write_metric(writer);
53  }
Hamiltonian< Model, BaseRNG >::PointType z_
Definition: base_hmc.hpp:182
template<class Model, template< class, class > class Hamiltonian, template< class > class Integrator, class BaseRNG>
void stan::mcmc::base_hmc< Model, Hamiltonian, Integrator, BaseRNG >::write_sampler_state ( callbacks::writer writer)
inlinevirtual

write stepsize and elements of mass matrix

Reimplemented from stan::mcmc::base_mcmc.

Definition at line 59 of file base_hmc.hpp.

Referenced by TEST().

59  {
60  write_sampler_stepsize(writer);
61  write_sampler_metric(writer);
62  }
void write_sampler_metric(callbacks::writer &writer)
Definition: base_hmc.hpp:51
void write_sampler_stepsize(callbacks::writer &writer)
Definition: base_hmc.hpp:41
template<class Model, template< class, class > class Hamiltonian, template< class > class Integrator, class BaseRNG>
void stan::mcmc::base_hmc< Model, Hamiltonian, Integrator, BaseRNG >::write_sampler_stepsize ( callbacks::writer writer)
inline

format and write stepsize

Definition at line 41 of file base_hmc.hpp.

Referenced by stan::mcmc::base_hmc< Model, diag_e_metric, expl_leapfrog, BaseRNG >::write_sampler_state().

41  {
42  std::stringstream nominal_stepsize;
43  nominal_stepsize << "Step size = " << get_nominal_stepsize();
44  writer(nominal_stepsize.str());
45  }
double get_nominal_stepsize()
Definition: base_hmc.hpp:157
template<class Model, template< class, class > class Hamiltonian, template< class > class Integrator, class BaseRNG>
Hamiltonian<Model, BaseRNG>::PointType& stan::mcmc::base_hmc< Model, Hamiltonian, Integrator, BaseRNG >::z ( )
inline

Definition at line 148 of file base_hmc.hpp.

Referenced by TEST().

148  {
149  return z_;
150  }
Hamiltonian< Model, BaseRNG >::PointType z_
Definition: base_hmc.hpp:182

Member Data Documentation

template<class Model, template< class, class > class Hamiltonian, template< class > class Integrator, class BaseRNG>
double stan::mcmc::base_hmc< Model, Hamiltonian, Integrator, BaseRNG >::epsilon_
protected
template<class Model, template< class, class > class Hamiltonian, template< class > class Integrator, class BaseRNG>
double stan::mcmc::base_hmc< Model, Hamiltonian, Integrator, BaseRNG >::epsilon_jitter_
protected
template<class Model, template< class, class > class Hamiltonian, template< class > class Integrator, class BaseRNG>
Hamiltonian<Model, BaseRNG> stan::mcmc::base_hmc< Model, Hamiltonian, Integrator, BaseRNG >::hamiltonian_
protected
template<class Model, template< class, class > class Hamiltonian, template< class > class Integrator, class BaseRNG>
Integrator<Hamiltonian<Model, BaseRNG> > stan::mcmc::base_hmc< Model, Hamiltonian, Integrator, BaseRNG >::integrator_
protected
template<class Model, template< class, class > class Hamiltonian, template< class > class Integrator, class BaseRNG>
double stan::mcmc::base_hmc< Model, Hamiltonian, Integrator, BaseRNG >::nom_epsilon_
protected
template<class Model, template< class, class > class Hamiltonian, template< class > class Integrator, class BaseRNG>
BaseRNG& stan::mcmc::base_hmc< Model, Hamiltonian, Integrator, BaseRNG >::rand_int_
protected
template<class Model, template< class, class > class Hamiltonian, template< class > class Integrator, class BaseRNG>
boost::uniform_01<BaseRNG&> stan::mcmc::base_hmc< Model, Hamiltonian, Integrator, BaseRNG >::rand_uniform_
protected
template<class Model, template< class, class > class Hamiltonian, template< class > class Integrator, class BaseRNG>
Hamiltonian<Model, BaseRNG>::PointType stan::mcmc::base_hmc< Model, Hamiltonian, Integrator, BaseRNG >::z_
protected

The documentation for this class was generated from the following file: