Class mcmc_stepper_rw (o2scl)

O2scl : Class List

template<class func_t, class data_t, class vec_t>
class mcmc_stepper_rw : public o2scl::mcmc_stepper_base<func_t, data_t, vec_t>

A simple random-walk stepper for MCMC.

This stepper performs a random walk. Given the parameter \( p_i \), the lower limit \( \ell_i \), the upper limit \( u_i \), a random number \( r_i \in [0,1) \), and the “step factor” \( s_i \), the new coordinate \( p_{\mathrm{new,i}} \) is

\[ p_{\mathrm{new,i}} = p_i + \frac{(2 r_i-1)}{s_i}(u_i - \ell_i) \]
The value of \( s_i \) is taken from
step_fac[k % step_fac.size()]
Thus larger values of step_fac result in smaller steps.

If the final point in parameter space, \( p_{\mathrm{new}} \), is out of bounds, then the value of func_ret is set to mcmc_stepper_base::mcmc_skip (which will lead to a rejection in mcmc_para_base or its children).

Public Functions

inline virtual const char *step_type()

Stepper type, “RW”.

inline mcmc_stepper_rw()
inline virtual ~mcmc_stepper_rw()
inline virtual void write_params(o2scl_hdf::hdf_file &hf)

Write stepper parameters to the HDF5 file.

inline virtual void step(size_t i_thread, size_t n_params, func_t &f, const vec_t &current, vec_t &next, double w_current, double &w_next, const vec_t &low, const vec_t &high, int &func_ret, bool &accept, data_t &dat, rng<> &r, int verbose)

Construct a step.

This function constructs next and w_next, the next point and log weight in parameter space. The objective function f is then evaluated at the new point, the return value is placed in func_ret, and the step acceptance or rejection is stored in accept.

Public Members

vec_t step_fac

The factor controlling the step size (default is a 1-element vector containing 2.0)