Stan  1.0
probability, sampling & optimization
Public Member Functions | List of all members
stan::mcmc::adaptive_hmc< BaseRNG > Class Template Reference

Adaptive Hamiltonian Monte Carlo (HMC) sampler. More...

#include <adaptive_hmc.hpp>

Inheritance diagram for stan::mcmc::adaptive_hmc< BaseRNG >:
stan::mcmc::hmc_base< boost::mt19937 > stan::mcmc::adaptive_sampler

Public Member Functions

 adaptive_hmc (stan::model::prob_grad &model, int L, double epsilon=-1, double epsilon_pm=0.0, bool epsilon_adapt=true, double delta=0.651, double gamma=0.05, BaseRNG rand_int=BaseRNG(std::time(0)), const std::vector< double > *params_r=0, const std::vector< int > *params_i=0)
 Construct an adaptive Hamiltonian Monte Carlo (HMC) sampler for the specified model, using the specified step size and number of leapfrog steps, with the specified random seed for randomization. More...
 
 ~adaptive_hmc ()
 Destructor. More...
 
virtual sample next_impl ()
 Returns the next sample. More...
 
virtual void write_sampler_param_names (std::ostream &o)
 Write out any sampler-specific parameter names for output. More...
 
virtual void write_sampler_params (std::ostream &o)
 Write out any sampler-specific parameters for output. More...
 
virtual void get_sampler_param_names (std::vector< std::string > &names)
 Get any sampler-specific parameter namess. More...
 
virtual void get_sampler_params (std::vector< double > &values)
 Get any sampler-specific parameters. More...
 
- Public Member Functions inherited from stan::mcmc::hmc_base< boost::mt19937 >
 hmc_base (stan::model::prob_grad &model, double epsilon=-1, double epsilon_pm=0.0, bool epsilon_adapt=true, double delta=0.651, double gamma=0.05, boost::mt19937 rand_int=boost::mt19937(std::time(0)), const std::vector< double > *params_r=0, const std::vector< int > *params_i=0)
 Construct a base HMC sampler. More...
 
virtual ~hmc_base ()
 
virtual void set_params (const std::vector< double > &x, const std::vector< int > &z)
 Sets the model real and integer parameters to the specified values. More...
 
void set_params_r (const std::vector< double > &x)
 Sets the model real parameters to the specified values and update gradients and log probability to match. More...
 
void set_params_i (const std::vector< int > &z)
 Sets the model integer parameters to the specified values and update gradients and log probability to match. More...
 
bool varying_epsilon ()
 
virtual void adapt_off ()
 Turn off parameter adaptation. More...
 
virtual void get_parameters (std::vector< double > &params)
 Write the step size into position 1 of the specified vector. More...
 
- Public Member Functions inherited from stan::mcmc::adaptive_sampler
 adaptive_sampler (bool adapt, std::ostream *error_msgs=0, std::ostream *output_msgs=0)
 Constructs an adaptive sampler with specified adaptation status. More...
 
virtual ~adaptive_sampler ()
 Destructor. More...
 
void set_error_stream (std::ostream &error_msgs)
 Set the stream into which errors will be written as the sampler runs. More...
 
void unset_error_stream ()
 Unset the stream into which errors are written to 0 so that error messages are ignored. More...
 
void set_output_stream (std::ostream &output_msgs)
 Set the stream into which output will be written as the sampler runs. More...
 
void unset_output_stream ()
 Unset the stream into which errors are written to 0 so that output messages are ignored. More...
 
sample next ()
 Returns the next sample from this sampler. More...
 
double mean_stat ()
 Returns the value of the statistic we are trying to coerce. More...
 
void set_mean_stat (double v)
 Sets the mean statistic to the specified value. More...
 
void update_mean_stat (double avg_eta, double adapt_stat)
 Updates the mean statistic given the specified adaptation statistic and weighting. More...
 
unsigned int nfevals ()
 Returns the number of times that the (possibly unnormalized) log probability function has been evaluated by this sampler. More...
 
void nfevals_plus_eq (int n)
 Add the specified number of evaluations to the number of function evaluations. More...
 
int n_steps ()
 Return the number of iterations for this sampler. More...
 
int n_adapt_steps ()
 Return how many iterations parameter adaptation has happened for. More...
 
virtual void adapt_on ()
 Turn on parameter adaptation. More...
 
bool adapting ()
 Return whether or not parameter adaptation is on. More...
 
virtual void write_adaptation_params (std::ostream &o)
 Use this method to write the adaptation parameters into the output. More...
 

Additional Inherited Members

- Protected Member Functions inherited from stan::mcmc::hmc_base< boost::mt19937 >
virtual void find_reasonable_parameters ()
 Search for a roughly reasonable (within a factor of 2) setting of the step size epsilon. More...
 
void adaptation_init (double epsilon_scale)
 
- Protected Attributes inherited from stan::mcmc::hmc_base< boost::mt19937 >
stan::model::prob_grad_model
 
double _epsilon
 
double _epsilon_pm
 
double _epsilon_last
 
bool _epsilon_adapt
 
const double _delta
 
const double _gamma
 
DualAverage _da
 
boost::mt19937 _rand_int
 
boost::variate_generator< boost::mt19937 &, boost::normal_distribution<> > _rand_unit_norm
 
boost::uniform_01< boost::mt19937 & > _rand_uniform_01
 
std::vector< double > _x
 
std::vector< int > _z
 
std::vector< double > _g
 
double _logp
 
- Protected Attributes inherited from stan::mcmc::adaptive_sampler
bool _adapt
 
unsigned int _n_steps
 
int _n_adapt_steps
 
unsigned int _nfevals
 
double _mean_stat
 
std::ostream * _error_msgs
 
std::ostream * _output_msgs
 

Detailed Description

template<class BaseRNG = boost::mt19937>
class stan::mcmc::adaptive_hmc< BaseRNG >

Adaptive Hamiltonian Monte Carlo (HMC) sampler.

adaptive_hmc automatically adapts the step size, epsilon, to try to coerce the average acceptance probability to some value, delta.

The HMC sampler requires a probability model with the ability to compute gradients. This is provided through an instance of stan::model::prob_grad.

Samples from the sampler are returned through the superclass stan::mcmc::sample and adaptation is handled through the superclass stan::mcmc::adaptive_sampler.

Definition at line 41 of file adaptive_hmc.hpp.

Constructor & Destructor Documentation

◆ adaptive_hmc()

template<class BaseRNG = boost::mt19937>
stan::mcmc::adaptive_hmc< BaseRNG >::adaptive_hmc ( stan::model::prob_grad model,
int  L,
double  epsilon = -1,
double  epsilon_pm = 0.0,
bool  epsilon_adapt = true,
double  delta = 0.651,
double  gamma = 0.05,
BaseRNG  rand_int = BaseRNG(std::time(0)),
const std::vector< double > *  params_r = 0,
const std::vector< int > *  params_i = 0 
)
inline

Construct an adaptive Hamiltonian Monte Carlo (HMC) sampler for the specified model, using the specified step size and number of leapfrog steps, with the specified random seed for randomization.

If the same seed is used twice, the series of samples should be the same. This property is most helpful for testing. If no random seed is specified, the std::time(0) function is called from the ctime library.

Parameters
modelProbability model with gradients.
LNumber of leapfrog steps per simulation.
deltaTarget value of E[acceptance probability]. Optional; defaults to the value of 0.651, which has some theoretical justification.
epsilonHamiltonian dynamics simulation step size. Optional; if not specified or set < 0, find_reasonable_parameters() will be called to initialize epsilon.
epsilon_pm
epsilon_adapt
deltaTarget value of E[acceptance probability]. Optional; defaults to the value of 0.651, which has some theoretical justification.
gammaRegularization parameter. See stan::mcmc::DualAverage.
rand_intSeed for random number generator; optional, if not specified, generate new seen based on system time.

Definition at line 77 of file adaptive_hmc.hpp.

◆ ~adaptive_hmc()

template<class BaseRNG = boost::mt19937>
stan::mcmc::adaptive_hmc< BaseRNG >::~adaptive_hmc ( )
inline

Destructor.

The implementation for this class is a no-op.

Definition at line 104 of file adaptive_hmc.hpp.

Member Function Documentation

◆ get_sampler_param_names()

template<class BaseRNG = boost::mt19937>
virtual void stan::mcmc::adaptive_hmc< BaseRNG >::get_sampler_param_names ( std::vector< std::string > &  names)
inlinevirtual

Get any sampler-specific parameter namess.

Parameters
[out]namesOutput vector to which param names are written.

Reimplemented from stan::mcmc::adaptive_sampler.

Definition at line 185 of file adaptive_hmc.hpp.

◆ get_sampler_params()

template<class BaseRNG = boost::mt19937>
virtual void stan::mcmc::adaptive_hmc< BaseRNG >::get_sampler_params ( std::vector< double > &  values)
inlinevirtual

Get any sampler-specific parameters.

Parameters
[out]valuesOutput vector to which params are written. All values are casted to type double. This function should match get_sampler_param_names.

Reimplemented from stan::mcmc::adaptive_sampler.

Definition at line 191 of file adaptive_hmc.hpp.

◆ next_impl()

template<class BaseRNG = boost::mt19937>
virtual sample stan::mcmc::adaptive_hmc< BaseRNG >::next_impl ( )
inlinevirtual

Returns the next sample.

Returns
The next sample.

Implements stan::mcmc::adaptive_sampler.

Definition at line 113 of file adaptive_hmc.hpp.

◆ write_sampler_param_names()

template<class BaseRNG = boost::mt19937>
virtual void stan::mcmc::adaptive_hmc< BaseRNG >::write_sampler_param_names ( std::ostream &  o)
inlinevirtual

Write out any sampler-specific parameter names for output.

This method must match write_sampler_params() in terms of number of parameters written.

Params should be writte starting with a comma, then the first parameter, then a comma, then the second parameter, ending on the final parameter.

The base class implementation is a no-op.

Parameters
[out]oOutput stream to which param names are written.

Reimplemented from stan::mcmc::adaptive_sampler.

Definition at line 175 of file adaptive_hmc.hpp.

◆ write_sampler_params()

template<class BaseRNG = boost::mt19937>
virtual void stan::mcmc::adaptive_hmc< BaseRNG >::write_sampler_params ( std::ostream &  o)
inlinevirtual

Write out any sampler-specific parameters for output.

This method must match write_sampler_param_names() in terms of number of parameters written.

Params should be writte starting with a comma, then the first parameter, then a comma, then the second parameter, ending on the final parameter.

The base class implementation is a no-op.

Parameters
[out]oOutput stream to which params are written.

Reimplemented from stan::mcmc::adaptive_sampler.

Definition at line 180 of file adaptive_hmc.hpp.


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

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