AdaptiveMCMC.AdaptiveMetropolisMethod
r = AdaptiveMetropolis(x0; kwargs)
r = AdaptiveMetropolis(x0, [scale, [step]])

Constructor for AdaptiveMetropolis state.

Arguments

  • x0: The initial state vector

Keyword arguments

  • scale: Scaling parameter; default 2.38/sqrt(d) where d is dimension.
  • step: Step size object; default PolynomialStepSize(1.0)
  • S_init: Initial proposal shape Cholesky factor; default identity_cholesky(x0)

If s is RWMState, then proposal samples may be drawn calling draw!(s, r) and adaptation is performed with adapt!(r, s) or adapt_rb!(r, s, α).

AdaptiveMCMC.AdaptiveScalingMetropolisMethod
r = AdaptiveScalingMetropolis(x0; kwargs)
r = AdaptiveScalingMetropolis(x0, [acc_target, [scale, [step]]])
r = AdaptiveScalingMetropolis(acc_target, scale, step)

Constructor for AdaptiveScalingMetropolis state.

Arguments

  • x0: The initial state vector (not used).

Keyword arguments

  • acc_target: Desired mean accept rate; default 0.44 for univariate and 0.234 for multivariate.
  • scale: Initial scaling; default 1.0.
  • step: Step size object; default PolynomialStepSize(0.66).

If s is RWMState, then proposal samples may be drawn calling draw!(s, r) and adaptation is performed with adapt!(r, s, α).

AdaptiveMCMC.AdaptiveScalingWithinAdaptiveMetropolisMethod
r = AdaptiveScalingWithinAdaptiveMetropolis(x0; kwargs)
r = AdaptiveScalingWithinAdaptiveMetropolis(x0, [acc_target, 
      [scale, [stepAM, [stepASM]]]])

Constructor for AdaptiveScalingWithinAdaptiveMetropolis state.

Arguments

  • x0: The initial state vector.

Keyword arguments

  • acc_target: Desired mean accept rate; default 0.234.
  • scale: Initial scaling; default 2.38/sqrt(d) where d is the dimension.
  • stepAM: Step size object for covariance adaptation; default PolynomialStepSize(0.66).
  • stepASM: Step size object for scaling adaptation; default PolynomialStepSize(0.66).

If s is RWMState, then proposal samples may be drawn calling draw!(s, r) and adaptation is performed with adapt!(r, s, α) or adapt_rb!(r, s, α).

AdaptiveMCMC.PolynomialStepSizeMethod
gamma = PolynomialStepSize(eta::AbstractFloat, [c::AbstractFloat=1.0]))

Constructor for PolynomialStepSize.

Arguments

  • eta: The step size exponent, should be within (1/2,1].
  • c: Scaling factor; default 1.0.
AdaptiveMCMC.RAMStepSizeMethod
gamma = RAMStepSize(eta::AbstractFloat, d::Int)

Constructor for RAM step size.

Arguments

  • eta: The step size exponent, should be within (1/2,1].
  • d: State dimension.
AdaptiveMCMC.RWMStateMethod
s = RWMState(x0, [rng, [q!]])

Constructor for RWMState (random walk Metropolis state).

Arguments

  • x0: The initial state vector
  • rng: Random number generator; default Random.GLOBAL_RNG.
  • q!: Symmetric proposal distribution; default randn!. Called by q!(rng, u) which puts a draw to vector u

If s is RWMState, then s.x is the current state. New proposal may be drawn to s.y by calling draw!, and the proposal is accepted by calling accept!(s).

AdaptiveMCMC.RobustAdaptiveMetropolisMethod
r = RobustAdaptiveMetropolis(x0; kwargs)
r = RobustAdaptiveMetropolis(x0, [acc_target, [step]])

Constructor for RobustAdaptiveMetropolis state.

Arguments

  • x0: The initial state vector

Keyword arguments

  • acc_target: Desired mean accept rate; default 0.234.
  • step: Step size object; default RAMStepSize(0.66,d) where d is state dimension.
  • S_init: Initial proposal shape Cholesky factor; default identity_cholesky(x0)

If s is RWMState, then proposal samples may be drawn calling draw!(s, r) and adaptation is performed with adapt!(r, s, α).

AdaptiveMCMC.accept!Method

accept!(s::RWMState)

Accept the proposal, that is, set copy s.y to s.x.

AdaptiveMCMC.adapt_rb!Method

adapt_rb!(s, r, α)

Rao-Blackwellised adaptation step of Adaptive Metropolis as suggested by Andrieu & Thoms (Statist. Comput. 2008).

Arguments:

  • s: AdaptiveMetropolis object
  • r: RWMState object
  • α: Acceptance rate

NB: This function should be called before calling accept!(r).

AdaptiveMCMC.adaptive_rwmMethod
out = adaptive_rwm(x0, log_p, n; kwargs)

Generic adaptive random walk Metropolis algorithm from initial state vector x0 targetting log probability density log_p run for n iterations, including adaptive parallel tempering.

Arguments

  • x0::Vector{<:AbstractFloat}: The initial state vector
  • log_p::Function: Function that returns log probability density values (up to an additive constant) for any state vector.
  • n::Int: Total number of iterations

Keyword arguments

  • algorithm::Symbol: The random walk adaptation algorithm; current choices are :ram (default), :am, :asm, :aswam and :rwm. (Alternatively, if algorithm is a vector of AdaptState, then this will be used as an initial state for adaptation.)
  • b::Int: Burn-in length: b:th sample is the first saved sample. Default ⌊n/5⌋
  • thin::Int: Thinning factor; only every thin:th sample is stored; default 1
  • fulladapt::Bool: Whether to adapt after burn-in; default true
  • Sp: Saved adaptive state from output to restart MCMC; default nothing
  • Rp: Saved rng state from output to restart MCMC; default nothing
  • indp::Int: Index of saved adaptive state to restart MCMC; default 0
  • rng::AbstractRNG: Random number generator; default Random.GLOBAL_RNG
  • q::Function: Zero-mean symmetric proposal generator (with arguments x and rng); default q=randn!(x, rng)
  • L::Int: Number of parallel tempering levels
  • acc_sw::AbstractFloat: Desired acceptance rate between level swaps; default 0.234
  • all_levels::Bool: Whether to store output of all levels; default false
  • log_pr::Function: Log-prior density function; default log_pr(x) = 0.0.
  • swaps::Symbol: Swap strategy, one of: :single (default, single randomly picked swap) :randperm (swap in random order) :sweep (up- or downward sweep, picked at random) :nonrev (alternate even/odd sites as in Syed, Bouchard-Côté, Deligiannidis, Doucet, arXiv:1905.02939)
  • progress::Union{Bool,Progress}: Whether a progress meter is shown; default false

Note that if log_pr is supplied, then log_p(x) is regarded as the log-likelihood (or, equivalently, log-target is log_p(x) + log_pr(x)). Tempering is only applied to log_p, not to log_pr.

The output out.X contains the simulated samples (column vectors).out.allX[k]fork>=2` contain higher temperature auxiliary chains (if requested)

Examples

log_p(x) = -.5*sum(x.^2)
o = adaptive_rwm(zeros(2), log_p, 10_000; algorithm=:am)
using MCMCChains, StatsPlots # Assuming MCMCChains & StatsPlots are installed...
c = Chains(o.X[1]', start=o.params.b, thin=o.params.thin); plot(c)
AdaptiveMCMC.draw!Method

draw!(s::RWMState, sc::AbstractFloat)

Draw proposal from s.x to s.y using scaling sc.

AdaptiveMCMC.draw!Method

draw!(s::RWMState, L::Cholesky, sc::AbstractFloat)

Draw proposal from s.x to s.y using scaling sc*L.L, with default sc=1.0.