|
template<class ExprType , class ConfigType , class Sampler > |
MCMCResult | base_mcmc (const ExprType &expr, const ConfigType &config, Sampler f) |
|
double | hamiltonian (double potential, double kinetic) |
|
template<class ADExprType , class MatType > |
double | reset_autodiff (ADExprType &ad_expr, Eigen::MatrixBase< MatType > &adjoints, Eigen::MatrixBase< MatType > &tp_adjoints) |
|
template<class ADExprType , class MatType , class MomentumHandlerType > |
double | leapfrog (ADExprType &ad_expr, Eigen::MatrixBase< MatType > &theta, Eigen::MatrixBase< MatType > &theta_adj, Eigen::MatrixBase< MatType > &tp_adj, Eigen::MatrixBase< MatType > &r, const MomentumHandlerType &m_handler, double epsilon, bool reuse_adj) |
|
template<class MatType1 , class MatType2 , class MatType3 > |
bool | check_entropy (const MatType1 &rho, const MatType2 &p_beg_scaled, const MatType3 &p_end_scaled) |
|
template<class InputType , class UniformDistType , class GenType , class MomentumHandlerType > |
TreeOutput | build_tree (size_t n_params, InputType &input, uint8_t depth, UniformDistType &unif_sampler, GenType &gen, MomentumHandlerType &momentum_handler, double *tree_cache) |
|
template<class ADExprType , class MatType , class GenType , class MomentumHandlerType > |
double | find_reasonable_epsilon (double eps, ADExprType &ad_expr, MatType &theta, MatType &theta_adj, MatType &tp_adj, GenType &gen, MomentumHandlerType &momentum_handler) |
|
template<class ProgramType , class OffsetPackType , class MCMCResultType , class NUTSConfigType = NUTSConfig<>> |
void | nuts_ (ProgramType &program, const NUTSConfigType &config, const OffsetPackType &pack, MCMCResultType &res) |
|
template<class ProgramType , class OffsetPackType , class MCMCResultType > |
void | mh_ (const ProgramType &program, const MHConfig &config, const OffsetPackType &pack, MCMCResultType &res) |
|
size_t | random_seed () |
|
template<class UniformDistType , class GenType > |
bool | accept_or_reject (double p, UniformDistType &&unif_sampler, GenType &&gen) |
|
template<class InputType , class UniformDistType , class GenType , class MomentumHandlerType >
TreeOutput ppl::mcmc::build_tree |
( |
size_t |
n_params, |
|
|
InputType & |
input, |
|
|
uint8_t |
depth, |
|
|
UniformDistType & |
unif_sampler, |
|
|
GenType & |
gen, |
|
|
MomentumHandlerType & |
momentum_handler, |
|
|
double * |
tree_cache |
|
) |
| |
Building binary tree for sampling candidates. Helper function to obtain the forward/backward-most position and momentum. Accept/reject policy is based on UniformDistType parameter and GenType
Note that the caller, i.e. nuts(), MUST have theta_adj already pre-computed (theta_adj is a member of input and input will be an instance of TreeInput).
- Parameters
-
n_params | number of (continuous) parameters |
input | TreeInput-like input object |
depth | current depth of building tree |
unif_sampler | an object like std::uniform_distribution(0,1) used for metropolis acceptance |
gen | rng device |
momentum_handler | MomentumHandler-like object to compute kinetic energy and momentum |
tree_cache | pointer to cache memory that will be used by build_tree. The array of doubles must be of size n_params * 7 * max_depth. |
template<class ADExprType , class MatType , class MomentumHandlerType >
double ppl::mcmc::leapfrog |
( |
ADExprType & |
ad_expr, |
|
|
Eigen::MatrixBase< MatType > & |
theta, |
|
|
Eigen::MatrixBase< MatType > & |
theta_adj, |
|
|
Eigen::MatrixBase< MatType > & |
tp_adj, |
|
|
Eigen::MatrixBase< MatType > & |
r, |
|
|
const MomentumHandlerType & |
m_handler, |
|
|
double |
epsilon, |
|
|
bool |
reuse_adj |
|
) |
| |
|
inline |
Leapfrog algorithm. Expects theta, theta_adj, tp_adj, r to be Eigen column-like objects.
Updates theta, theta_adj, tp_adj, and r to contain the new leaped values and adjoints.
- Parameters
-
ad_expr | AD expression representing L(theta) It must be built such that values are read from theta and adjoints are placed into theta_adj. |
theta | theta at which we want to start leaping |
theta_adj | adjoint for theta. If not reusing, resets adjoints first. |
tp_adj | adjoints for transformed parameters. If not reusing, resets adjoints first. |
r | momentum vector to start leaping |
m_handler | momentum handler to compute correct dkinetic/dr |
epsilon | step size |
reuse_adj | flag to not compute gradient of L(theta) if user can guarantee that theta_adj currently has it. |
- Returns
- new potential energy (L(theta'))