Stan  2.5.0
probability, sampling & optimization
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Pages
adapt_unit_e_nuts.hpp
Go to the documentation of this file.
1 #ifndef STAN__MCMC__ADAPT__UNIT__E__NUTS__BETA
2 #define STAN__MCMC__ADAPT__UNIT__E__NUTS__BETA
3 
6 
7 namespace stan {
8 
9  namespace mcmc {
10 
11  // The No-U-Turn Sampler (NUTS) on a
12  // Euclidean manifold with unit metric
13  // and adaptive stepsize
14 
15  template <typename M, class BaseRNG>
16  class adapt_unit_e_nuts: public unit_e_nuts<M, BaseRNG>,
17  public stepsize_adapter {
18 
19  public:
20 
21  adapt_unit_e_nuts(M &m, BaseRNG& rng,
22  std::ostream* o = &std::cout, std::ostream* e = 0):
23  unit_e_nuts<M, BaseRNG>(m, rng, o, e) {};
24 
26 
27  sample transition(sample& init_sample) {
28 
30 
31  if (this->adapt_flag_)
33 
34  return s;
35 
36  }
37 
41  }
42 
43  };
44 
45  } // mcmc
46 
47 } // stan
48 
49 
50 #endif
void complete_adaptation(double &epsilon)
double accept_stat() const
Definition: sample.hpp:46
adapt_unit_e_nuts(M &m, BaseRNG &rng, std::ostream *o=&std::cout, std::ostream *e=0)
void learn_stepsize(double &epsilon, double adapt_stat)
virtual void disengage_adaptation()
sample transition(sample &init_sample)
double e()
Return the base of the natural logarithm.
Definition: constants.hpp:86
stepsize_adaptation stepsize_adaptation_

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