AdvancedHMC.jl
Documentation for AdvancedHMC.jl
Types
AdvancedHMC.ClassicNoUTurn
— Typestruct ClassicNoUTurn{F<:AbstractFloat} <: AdvancedHMC.DynamicTerminationCriterion
Classic No-U-Turn criterion as described in Eq. (9) in [1].
Informally, this will terminate the trajectory expansion if continuing the simulation either forwards or backwards in time will decrease the distance between the left-most and right-most positions.
Fields
max_depth::Int64
Δ_max::AbstractFloat
References
- Hoffman, M. D., & Gelman, A. (2014). The No-U-Turn Sampler: adaptively setting path lengths in Hamiltonian Monte Carlo. Journal of Machine Learning Research, 15(1), 1593-1623. (arXiv)
AdvancedHMC.HMCSampler
— TypeHMCSampler
An AbstractMCMC.AbstractSampler
for kernels in AdvancedHMC.jl.
Fields
metric
: Choice of initial metricAbstractMetric
. The metric type will be preserved during adaption.adaptor
:AdvancedHMC.Adaptation.AbstractAdaptor
.
Notes
Note that all the fields have the prefix initial_
to indicate that these will not necessarily correspond to the kernel
, metric
, and adaptor
after sampling.
To access the updated fields, use the resulting HMCState
.
AdvancedHMC.HMC
— TypeHMC(ϵ::Real, n_leapfrog::Int)
Hamiltonian Monte Carlo sampler with static trajectory.
Fields
n_leapfrog
: Number of leapfrog steps.integrator
: Choice of integrator, specified either using aSymbol
orAbstractIntegrator
metric
: Choice of initial metric;Symbol
means it is automatically initialised. The metric type will be preserved during automatic initialisation and adaption.
AdvancedHMC.NUTS
— TypeNUTS(δ::Real; max_depth::Int=10, Δ_max::Real=1000, integrator = :leapfrog, metric = :diagonal)
No-U-Turn Sampler (NUTS) sampler.
Fields
δ
: Target acceptance rate for dual averaging.max_depth
: Maximum doubling tree depth.Δ_max
: Maximum divergence during doubling tree.integrator
: Choice of integrator, specified either using aSymbol
orAbstractIntegrator
metric
: Choice of initial metric;Symbol
means it is automatically initialised. The metric type will be preserved during automatic initialisation and adaption.
AdvancedHMC.HMCDA
— TypeHMCDA(δ::Real, λ::Real, integrator = :leapfrog, metric = :diagonal)
Hamiltonian Monte Carlo sampler with Dual Averaging algorithm.
Fields
δ
: Target acceptance rate for dual averaging.λ
: Target leapfrog length.integrator
: Choice of integrator, specified either using aSymbol
orAbstractIntegrator
metric
: Choice of initial metric;Symbol
means it is automatically initialised. The metric type will be preserved during automatic initialisation and adaption.
Notes
For more information, please view the following paper (arXiv link):
- Hoffman, Matthew D., and Andrew Gelman. "The No-U-turn sampler: adaptively setting path lengths in Hamiltonian Monte Carlo." Journal of Machine Learning Research 15, no. 1 (2014): 1593-1623.
Functions
StatsBase.sample
— Functionsample([rng], a, [wv::AbstractWeights])
Select a single random element of a
. Sampling probabilities are proportional to the weights given in wv
, if provided.
Optionally specify a random number generator rng
as the first argument (defaults to Random.default_rng()
).
sample([rng], a, [wv::AbstractWeights], n::Integer; replace=true, ordered=false)
Select a random, optionally weighted sample of size n
from an array a
using a polyalgorithm. Sampling probabilities are proportional to the weights given in wv
, if provided. replace
dictates whether sampling is performed with replacement. ordered
dictates whether an ordered sample (also called a sequential sample, i.e. a sample where items appear in the same order as in a
) should be taken.
Optionally specify a random number generator rng
as the first argument (defaults to Random.default_rng()
).
sample([rng], a, [wv::AbstractWeights], dims::Dims; replace=true, ordered=false)
Select a random, optionally weighted sample from an array a
specifying the dimensions dims
of the output array. Sampling probabilities are proportional to the weights given in wv
, if provided. replace
dictates whether sampling is performed with replacement. ordered
dictates whether an ordered sample (also called a sequential sample, i.e. a sample where items appear in the same order as in a
) should be taken.
Optionally specify a random number generator rng
as the first argument (defaults to Random.default_rng()
).
sample([rng], wv::AbstractWeights)
Select a single random integer in 1:length(wv)
with probabilities proportional to the weights given in wv
.
Optionally specify a random number generator rng
as the first argument (defaults to Random.default_rng()
).
sample(
rng::Random.AbatractRNG=Random.default_rng(),
model::AbstractModel,
sampler::AbstractSampler,
N_or_isdone;
kwargs...,
)
Sample from the model
with the Markov chain Monte Carlo sampler
and return the samples.
If N_or_isdone
is an Integer
, exactly N_or_isdone
samples are returned.
Otherwise, sampling is performed until a convergence criterion N_or_isdone
returns true
. The convergence criterion has to be a function with the signature
isdone(rng, model, sampler, samples, state, iteration; kwargs...)
where state
and iteration
are the current state and iteration of the sampler, respectively. It should return true
when sampling should end, and false
otherwise.
Keyword arguments
See https://turinglang.org/AbstractMCMC.jl/dev/api/#Common-keyword-arguments for common keyword arguments.
sample(
rng::Random.AbstractRNG=Random.default_rng(),
model::AbstractModel,
sampler::AbstractSampler,
parallel::AbstractMCMCEnsemble,
N::Integer,
nchains::Integer;
kwargs...,
)
Sample nchains
Monte Carlo Markov chains from the model
with the sampler
in parallel using the parallel
algorithm, and combine them into a single chain.
Keyword arguments
See https://turinglang.org/AbstractMCMC.jl/dev/api/#Common-keyword-arguments for common keyword arguments.
sample(
rng::Random.AbstractRNG=Random.default_rng(),
logdensity,
sampler::AbstractSampler,
N_or_isdone;
kwargs...,
)
Wrap the logdensity
function in a LogDensityModel
, and call sample
with the resulting model instead of logdensity
.
The logdensity
function has to support the LogDensityProblems.jl interface.
sample(
rng::Random.AbstractRNG=Random.default_rng(),
logdensity,
sampler::AbstractSampler,
parallel::AbstractMCMCEnsemble,
N::Integer,
nchains::Integer;
kwargs...,
)
Wrap the logdensity
function in a LogDensityModel
, and call sample
with the resulting model instead of logdensity
.
The logdensity
function has to support the LogDensityProblems.jl interface.
sample(
rng::AbstractRNG,
h::Hamiltonian,
κ::AbstractMCMCKernel,
θ::AbstractVecOrMat{T},
n_samples::Int,
adaptor::AbstractAdaptor=NoAdaptation(),
n_adapts::Int=min(div(n_samples, 10), 1_000);
drop_warmup::Bool=false,
verbose::Bool=true,
progress::Bool=false
)
Sample n_samples
samples using the proposal κ
under Hamiltonian h
.
- The randomness is controlled by
rng
.- If
rng
is not provided,GLOBAL_RNG
will be used.
- If
- The initial point is given by
θ
. - The adaptor is set by
adaptor
, for which the default is no adaptation.- It will perform
n_adapts
steps of adaptation, for which the default is the minimum of1_000
and 10% ofn_samples
- It will perform
drop_warmup
controls to drop the samples during adaptation phase or notverbose
controls the verbosityprogress
controls whether to show the progress meter or not
More types
AdvancedHMC.AbstractIntegrator
— Typeabstract type AbstractIntegrator
Represents an integrator used to simulate the Hamiltonian system.
Implementation
A AbstractIntegrator
is expected to have the following implementations:
stat
(@ref)nom_step_size
(@ref)step_size
(@ref)
AdvancedHMC.AbstractMCMCKernel
— Typeabstract type AbstractMCMCKernel
Abstract type for HMC kernels.
AdvancedHMC.AbstractMetric
— Typeabstract type AbstractMetric
Abstract type for preconditioning metrics.
AdvancedHMC.AbstractTerminationCriterion
— Typeabstract type AbstractTerminationCriterion
Abstract type for termination criteria for Hamiltonian trajectories, e.g. no-U-turn and fixed number of leapfrog integration steps.
AdvancedHMC.AbstractTrajectorySampler
— TypeHow to sample a phase-point from the simulated trajectory.
AdvancedHMC.BinaryTree
— TypeA full binary tree trajectory with only necessary leaves and information stored.
AdvancedHMC.ClassicNoUTurn
— Typestruct ClassicNoUTurn{F<:AbstractFloat} <: AdvancedHMC.DynamicTerminationCriterion
Classic No-U-Turn criterion as described in Eq. (9) in [1].
Informally, this will terminate the trajectory expansion if continuing the simulation either forwards or backwards in time will decrease the distance between the left-most and right-most positions.
Fields
max_depth::Int64
Δ_max::AbstractFloat
References
- Hoffman, M. D., & Gelman, A. (2014). The No-U-Turn Sampler: adaptively setting path lengths in Hamiltonian Monte Carlo. Journal of Machine Learning Research, 15(1), 1593-1623. (arXiv)
AdvancedHMC.DynamicTerminationCriterion
— Typeabstract type DynamicTerminationCriterion <: AdvancedHMC.AbstractTerminationCriterion
Abstract type for dynamic Hamiltonian trajectory termination criteria.
AdvancedHMC.EndPointTS
— TypeSamples the end-point of the trajectory.
AdvancedHMC.FixedIntegrationTime
— Typestruct FixedIntegrationTime{F<:AbstractFloat} <: AdvancedHMC.StaticTerminationCriterion
Standard HMC implementation with a fixed integration time.
Fields
λ::AbstractFloat
: Total length of the trajectory, i.e. takefloor(λ / integrator_step_size)
number of leapfrog steps.
References
- Neal, R. M. (2011). MCMC using Hamiltonian dynamics. Handbook of Markov chain Monte Carlo, 2(11), 2. (arXiv)
AdvancedHMC.FixedNSteps
— Typestruct FixedNSteps <: AdvancedHMC.StaticTerminationCriterion
Static HMC with a fixed number of leapfrog steps.
Fields
L::Int64
: Number of steps to simulate, i.e. length of trajectory will beL + 1
.
References
- Neal, R. M. (2011). MCMC using Hamiltonian dynamics. Handbook of Markov chain Monte Carlo, 2(11), 2. (arXiv)
AdvancedHMC.FullMomentumRefreshment
— TypeCompletly resample new momentum.
AdvancedHMC.GeneralisedNoUTurn
— Typestruct GeneralisedNoUTurn{F<:AbstractFloat} <: AdvancedHMC.DynamicTerminationCriterion
Generalised No-U-Turn criterion as described in Section A.4.2 in [1].
Fields
max_depth::Int64
Δ_max::AbstractFloat
References
- Betancourt, M. (2017). A Conceptual Introduction to Hamiltonian Monte Carlo. arXiv preprint arXiv:1701.02434.
AdvancedHMC.GeneralizedLeapfrog
— Typestruct GeneralizedLeapfrog{T<:(Union{AbstractVector{var"#s17"}, var"#s17"} where var"#s17"<:AbstractFloat)} <: AdvancedHMC.AbstractLeapfrog{T<:(Union{AbstractVector{var"#s17"}, var"#s17"} where var"#s17"<:AbstractFloat)}
Generalized leapfrog integrator with fixed step size ϵ
.
Fields
ϵ::Union{AbstractVector{var"#s17"}, var"#s17"} where var"#s17"<:AbstractFloat
: Step size.n::Int64
References
- Girolami, Mark, and Ben Calderhead. "Riemann manifold Langevin and Hamiltonian Monte Carlo methods." Journal of the Royal Statistical Society Series B: Statistical Methodology 73, no. 2 (2011): 123-214.
AdvancedHMC.HMC
— TypeHMC(ϵ::Real, n_leapfrog::Int)
Hamiltonian Monte Carlo sampler with static trajectory.
Fields
n_leapfrog
: Number of leapfrog steps.integrator
: Choice of integrator, specified either using aSymbol
orAbstractIntegrator
metric
: Choice of initial metric;Symbol
means it is automatically initialised. The metric type will be preserved during automatic initialisation and adaption.
AdvancedHMC.HMCDA
— TypeHMCDA(δ::Real, λ::Real, integrator = :leapfrog, metric = :diagonal)
Hamiltonian Monte Carlo sampler with Dual Averaging algorithm.
Fields
δ
: Target acceptance rate for dual averaging.λ
: Target leapfrog length.integrator
: Choice of integrator, specified either using aSymbol
orAbstractIntegrator
metric
: Choice of initial metric;Symbol
means it is automatically initialised. The metric type will be preserved during automatic initialisation and adaption.
Notes
For more information, please view the following paper (arXiv link):
- Hoffman, Matthew D., and Andrew Gelman. "The No-U-turn sampler: adaptively setting path lengths in Hamiltonian Monte Carlo." Journal of Machine Learning Research 15, no. 1 (2014): 1593-1623.
AdvancedHMC.HMCProgressCallback
— TypeHMCProgressCallback
A callback to be used with AbstractMCMC.jl's interface, replicating the logging behavior of the non-AbstractMCMC sample
.
Fields
pm
:Progress
meter from ProgressMeters.jl.progress
: Specifies whether or not to use display a progress bar.verbose
: Ifprogress
is not specified and this istrue
some information will be logged upon completion of adaptation.num_divergent_transitions
: Number of divergent transitions fo far.num_divergent_transitions_during_adaption
AdvancedHMC.HMCSampler
— TypeHMCSampler
An AbstractMCMC.AbstractSampler
for kernels in AdvancedHMC.jl.
Fields
metric
: Choice of initial metricAbstractMetric
. The metric type will be preserved during adaption.adaptor
:AdvancedHMC.Adaptation.AbstractAdaptor
.
Notes
Note that all the fields have the prefix initial_
to indicate that these will not necessarily correspond to the kernel
, metric
, and adaptor
after sampling.
To access the updated fields, use the resulting HMCState
.
AdvancedHMC.HMCState
— TypeHMCState
Represents the state of a HMCSampler
.
Fields
i
: Index of current iteration.transition
: CurrentTransition
.metric
: CurrentAbstractMetric
, possibly adapted.κ
: CurrentAbstractMCMCKernel
.adaptor
: CurrentAbstractAdaptor
.
AdvancedHMC.JitteredLeapfrog
— Typestruct JitteredLeapfrog{FT<:AbstractFloat, T<:Union{AbstractArray{FT<:AbstractFloat, 1}, FT<:AbstractFloat}} <: AdvancedHMC.AbstractLeapfrog{T<:Union{AbstractArray{FT<:AbstractFloat, 1}, FT<:AbstractFloat}}
Leapfrog integrator with randomly "jittered" step size ϵ
for every trajectory.
Fields
ϵ0::Union{AbstractVector{FT}, FT} where FT<:AbstractFloat
: Nominal (non-jittered) step size.jitter::AbstractFloat
: The proportion of the nominal step sizeϵ0
that may be added or subtracted.ϵ::Union{AbstractVector{FT}, FT} where FT<:AbstractFloat
: Current (jittered) step size.
Description
This is the same as LeapFrog
(@ref) but with a "jittered" step size. This means that at the beginning of each trajectory we sample a step size ϵ
by adding or subtracting from the nominal/base step size ϵ0
some random proportion of ϵ0
, with the proportion specified by jitter
, i.e. ϵ = ϵ0 - jitter * ϵ0 * rand()
. p Jittering might help alleviate issues related to poor interactions with a fixed step size:
- In regions with high "curvature" the current choice of step size might mean over-shoot leading to almost all steps being rejected. Randomly sampling the step size at the beginning of the trajectories can therefore increase the probability of escaping such high-curvature regions.
- Exact periodicity of the simulated trajectories might occur, i.e. you might be so unlucky as to simulate the trajectory forwards in time
L ϵ
and ending up at the same point (which results in non-ergodicity; see Section 3.2 in [1]). If momentum is refreshed before each trajectory, then this should not happen exactly but it can still be an issue in practice. Randomly choosing the step-sizeϵ
might help alleviate such problems.
References
- Neal, R. M. (2011). MCMC using Hamiltonian dynamics. Handbook of Markov chain Monte Carlo, 2(11), 2. (arXiv)
AdvancedHMC.Leapfrog
— Typestruct Leapfrog{T<:(Union{AbstractVector{var"#s26"}, var"#s26"} where var"#s26"<:AbstractFloat)} <: AdvancedHMC.AbstractLeapfrog{T<:(Union{AbstractVector{var"#s26"}, var"#s26"} where var"#s26"<:AbstractFloat)}
Leapfrog integrator with fixed step size ϵ
.
Fields
ϵ::Union{AbstractVector{var"#s26"}, var"#s26"} where var"#s26"<:AbstractFloat
: Step size.
AdvancedHMC.MultinomialTS
— Methodstruct MultinomialTS{F<:AbstractFloat} <: AdvancedHMC.AbstractTrajectorySampler
Multinomial sampler for a trajectory consisting only a leaf node.
- tree weight is the (unnormalised) energy of the leaf.
AdvancedHMC.MultinomialTS
— Methodstruct MultinomialTS{F<:AbstractFloat} <: AdvancedHMC.AbstractTrajectorySampler
Multinomial sampler for the starting single leaf tree. (Log) weights for leaf nodes are their (unnormalised) Hamiltonian energies.
Ref: https://github.com/stan-dev/stan/blob/develop/src/stan/mcmc/hmc/nuts/base_nuts.hpp#L226
AdvancedHMC.MultinomialTS
— Typestruct MultinomialTS{F<:AbstractFloat} <: AdvancedHMC.AbstractTrajectorySampler
Multinomial trajectory sampler carried during the building of the tree. It contains the weight of the tree, defined as the total probabilities of the leaves.
Fields
zcand::AdvancedHMC.PhasePoint
: Sampled candidatePhasePoint
.ℓw::AbstractFloat
: Total energy for the given tree, i.e. the sum of energies of all leaves.
AdvancedHMC.NUTS
— TypeNUTS(δ::Real; max_depth::Int=10, Δ_max::Real=1000, integrator = :leapfrog, metric = :diagonal)
No-U-Turn Sampler (NUTS) sampler.
Fields
δ
: Target acceptance rate for dual averaging.max_depth
: Maximum doubling tree depth.Δ_max
: Maximum divergence during doubling tree.integrator
: Choice of integrator, specified either using aSymbol
orAbstractIntegrator
metric
: Choice of initial metric;Symbol
means it is automatically initialised. The metric type will be preserved during automatic initialisation and adaption.
AdvancedHMC.PartialMomentumRefreshment
— Typestruct PartialMomentumRefreshment{F<:AbstractFloat} <: AdvancedHMC.AbstractMomentumRefreshment
Partial momentum refreshment with refresh rate α
.
Fields
α::AbstractFloat
See equation (5.19) [1]
r' = α⋅r + sqrt(1-α²)⋅G
where r is the momentum and G is a Gaussian random variable.
References
- Neal, Radford M. "MCMC using Hamiltonian dynamics." Handbook of markov chain monte carlo 2.11 (2011): 2.
AdvancedHMC.SliceTS
— Methodstruct SliceTS{F<:AbstractFloat} <: AdvancedHMC.AbstractTrajectorySampler
Slice sampler for the starting single leaf tree. Slice variable is initialized.
AdvancedHMC.SliceTS
— Methodstruct SliceTS{F<:AbstractFloat} <: AdvancedHMC.AbstractTrajectorySampler
Create a slice sampler for a single leaf tree:
- the slice variable is copied from the passed-in sampler
s
and - the number of acceptable candicates is computed by comparing the slice variable against the current energy.
AdvancedHMC.SliceTS
— Typestruct SliceTS{F<:AbstractFloat} <: AdvancedHMC.AbstractTrajectorySampler
Trajectory slice sampler carried during the building of the tree. It contains the slice variable and the number of acceptable condidates in the tree.
Fields
zcand::AdvancedHMC.PhasePoint
: Sampled candidatePhasePoint
.ℓu::AbstractFloat
: Slice variable in log-space.n::Int64
: Number of acceptable candidates, i.e. those with probability larger than slice variableu
.
AdvancedHMC.StaticTerminationCriterion
— Typeabstract type StaticTerminationCriterion <: AdvancedHMC.AbstractTerminationCriterion
Abstract type for a fixed number of leapfrog integration steps.
AdvancedHMC.StrictGeneralisedNoUTurn
— Typestruct StrictGeneralisedNoUTurn{F<:AbstractFloat} <: AdvancedHMC.DynamicTerminationCriterion
Generalised No-U-Turn criterion as described in Section A.4.2 in [1] with added U-turn check as described in [2].
Fields
max_depth::Int64
Δ_max::AbstractFloat
References
- Betancourt, M. (2017). A Conceptual Introduction to Hamiltonian Monte Carlo. arXiv preprint arXiv:1701.02434.
- https://github.com/stan-dev/stan/pull/2800
AdvancedHMC.TemperedLeapfrog
— Typestruct TemperedLeapfrog{FT<:AbstractFloat, T<:Union{AbstractArray{FT<:AbstractFloat, 1}, FT<:AbstractFloat}} <: AdvancedHMC.AbstractLeapfrog{T<:Union{AbstractArray{FT<:AbstractFloat, 1}, FT<:AbstractFloat}}
Tempered leapfrog integrator with fixed step size ϵ
and "temperature" α
.
Fields
ϵ::Union{AbstractVector{FT}, FT} where FT<:AbstractFloat
: Step size.α::AbstractFloat
: Temperature parameter.
Description
Tempering can potentially allow greater exploration of the posterior, e.g. in a multi-modal posterior jumps between the modes can be more likely to occur.
AdvancedHMC.Termination
— MethodTermination(s, nt, H0, H′)
Check termination of a Hamiltonian trajectory.
AdvancedHMC.Termination
— MethodTermination(s, nt, H0, H′)
Check termination of a Hamiltonian trajectory.
AdvancedHMC.Termination
— TypeTermination
Termination reasons
dynamic
: due to stoping criterianumerical
: due to large energy deviation from starting (possibly numerical errors)
AdvancedHMC.Trajectory
— Typestruct Trajectory{TS<:AdvancedHMC.AbstractTrajectorySampler, I<:AdvancedHMC.AbstractIntegrator, TC<:AdvancedHMC.AbstractTerminationCriterion}
Numerically simulated Hamiltonian trajectories.
AdvancedHMC.Transition
— Typestruct Transition{P<:AdvancedHMC.PhasePoint, NT<:NamedTuple}
A transition that contains the phase point and other statistics of the transition.
Fields
z::AdvancedHMC.PhasePoint
: Phase-point for the transition.stat::NamedTuple
: Statistics related to the transition, e.g. energy.
AdvancedHMC.Adaptation.AbstractAdaptor
— Typeabstract type AbstractAdaptor
Abstract type for HMC adaptors.
AdvancedHMC.Adaptation.NesterovDualAveraging
— TypeAn implementation of the Nesterov dual averaging algorithm to tune step size.
References
Hoffman, M. D., & Gelman, A. (2014). The No-U-Turn Sampler: adaptively setting path lengths in Hamiltonian Monte Carlo. Journal of Machine Learning Research, 15(1), 1593-1623. Nesterov, Y. (2009). Primal-dual subgradient methods for convex problems. Mathematical programming, 120(1), 221-259.