Relativistic MCMC¶
- class pints.RelativisticMCMC(x0, sigma0=None)[source]¶
Implements Relativistic Monte Carlo as described in [1].
Uses a physical analogy of a particle moving across a landscape under Hamiltonian dynamics to aid efficient exploration of parameter space. Introduces an auxilary variable – the momentum (
p_i
) of a particle moving in dimensioni
of negative log posterior space – which supplements the position (q_i
) of the particle in parameter space. The particle’s motion is dictated by solutions to Hamilton’s equations,\[\begin{split}dq_i/dt &= \partial H/\partial p_i\\ dp_i/dt &= - \partial H/\partial q_i.\end{split}\]The Hamiltonian is given by,
\[\begin{split}H(q,p) &= U(q) + KE(p)\\ &= -\text{log}(p(q|X)p(q)) + mc^2 (\Sigma_{i=1}^{d} p_i^2 / (m^2 c^2) + 1)^{0.5}\end{split}\]where
d
is the dimensionality of model,m
is the scalar ‘mass’ given to each particle (chosen to be 1 as default) andc
is the speed of light (chosen to be 10 by default).To numerically integrate Hamilton’s equations, it is essential to use a sympletic discretisation routine, of which the most typical approach is the leapfrog method,
\[\begin{split}p_i(t + \epsilon/2) &= p_i(t) - (\epsilon/2) d U(q_i(t))/dq_i\\ q_i(t + \epsilon) &= q_i(t) + \epsilon M^{-1}(p_i(t + \epsilon/2)) p_i(t + \epsilon/2)\\ p_i(t + \epsilon) &= p_i(t + \epsilon/2) - (\epsilon/2) d U(q_i(t + \epsilon))/dq_i\end{split}\]where relativistic mass (a scalar) is,
\[M(p) = m (\Sigma_{i=1}^{d} p_i^2 / (m^2 c^2) + 1)^{0.5}\]In particular, the algorithm we implement follows eqs. in section 2.1 of [1].
Extends
SingleChainMCMC
.References
- hamiltonian_threshold()[source]¶
Returns threshold difference in Hamiltonian value from one iteration to next which determines whether an iteration is divergent.
- in_initial_phase()¶
For methods that need an initial phase (see
needs_initial_phase()
), this method returnsTrue
if the method is currently configured to be in its initial phase. For other methods aNotImplementedError
is returned.
- needs_initial_phase()¶
Returns
True
if this method needs an initial phase, for example an adaptation-free period for adaptive covariance methods, or a warm-up phase for DREAM.
- replace(current, current_log_pdf, proposed=None)¶
Replaces the internal current position, current LogPDF, and proposed point (if any) by the user-specified values.
This method can only be used once the initial position and LogPDF have been set (so after at least 1 round of ask-and-tell).
This is an optional method, and some samplers may not support it.
- set_hamiltonian_threshold(hamiltonian_threshold)[source]¶
Sets threshold difference in Hamiltonian value from one iteration to next which determines whether an iteration is divergent.
- set_hyper_parameters(x)[source]¶
The hyper-parameter vector is
[leapfrog_steps, leapfrog_step_size, mass, c]
.
- set_initial_phase(in_initial_phase)¶
For methods that need an initial phase (see
needs_initial_phase()
), this method toggles the initial phase algorithm. For other methods aNotImplementedError
is returned.
- set_leapfrog_steps(steps)[source]¶
Sets the number of leapfrog steps to carry out for each iteration.