Metropolis-Adjusted Langevin Algorithm (MALA)

Implementation of the Metropolis-Adjusted Langevin Algorithm of Roberts and Tweedie [80] and Roberts and Stramer [79]. The sampler simulates autocorrelated draws from a distribution that can be specified up to a constant of proportionality. MALA is related to Hamiltonian Monte Carlo as described thoroughly by Girolami and Calderhead [39].

Model-Based Constructors

MALA(params::ElementOrVector{Symbol}, epsilon::Real; dtype::Symbol=:forward)
MALA(params::ElementOrVector{Symbol}, epsilon::Real, Sigma::Matrix{T<:Real}; dtype::Symbol=:forward)

Construct a Sampler object for MALA sampling. Parameters are assumed to be continuous, but may be constrained or unconstrained.

Arguments

  • params : stochastic node(s) to be updated with the sampler. Constrained parameters are mapped to unconstrained space according to transformations defined by the Stochastic unlist() function.
  • epsilon : factor by which the drift and covariance matrix of the proposal distribution are scaled.
  • Sigma : covariance matrix for the multivariate normal proposal distribution. The covariance matrix is relative to the unconstrained parameter space, where candidate draws are generated. If omitted, the identity matrix is assumed.
  • dtype : type of differentiation for gradient calculations. Options are
    • :central : central differencing.
    • :forward : forward differencing.

Value

Returns a Sampler{MALATune} type object.

Example

See the Dyes and other Examples.

Stand-Alone Function

sample!(v::MALAVariate)

Draw one sample from a target distribution using the MALA sampler. Parameters are assumed to be continuous and unconstrained.

Arguments

  • v : current state of parameters to be simulated.

Value

Returns v updated with simulated values and associated tuning parameters.

Example

The following example samples parameters in a simple linear regression model. Details of the model specification and posterior distribution can be found in the Supplement.

################################################################################
## Linear Regression
##   y ~ N(b0 + b1 * x, s2)
##   b0, b1 ~ N(0, 1000)
##   s2 ~ invgamma(0.001, 0.001)
################################################################################

using Mamba

## Data
data = Dict(
  :x => [1, 2, 3, 4, 5],
  :y => [1, 3, 3, 3, 5]
)

## Log-transformed Posterior(b0, b1, log(s2)) + Constant and Gradient Vector
logfgrad = function(x::DenseVector)
  b0 = x[1]
  b1 = x[2]
  logs2 = x[3]
  r = data[:y] - b0 - b1 * data[:x]
  logf = (-0.5 * length(data[:y]) - 0.001) * logs2 -
           (0.5 * dot(r, r) + 0.001) / exp(logs2) -
           0.5 * b0^2 / 1000 - 0.5 * b1^2 / 1000
  grad = [
    sum(r) / exp(logs2) - b0 / 1000,
    sum(data[:x] .* r) / exp(logs2) - b1 / 1000,
    -0.5 * length(data[:y]) - 0.001 + (0.5 * dot(r, r) + 0.001) / exp(logs2)
  ]
  logf, grad
end

## MCMC Simulation with Metropolis-Adjusted Langevin Algorithm
## Without (1) and with (2) a user-specified proposal covariance matrix
n = 5000
sim1 = Chains(n, 3, names = ["b0", "b1", "s2"])
sim2 = Chains(n, 3, names = ["b0", "b1", "s2"])
epsilon = 0.1
Sigma = eye(3)
theta1 = MALAVariate([0.0, 0.0, 0.0], epsilon, logfgrad)
theta2 = MALAVariate([0.0, 0.0, 0.0], epsilon, Sigma, logfgrad)
for i in 1:n
  sample!(theta1)
  sample!(theta2)
  sim1[i, :, 1] = [theta1[1:2]; exp(theta1[3])]
  sim2[i, :, 1] = [theta2[1:2]; exp(theta2[3])]
end
describe(sim1)
describe(sim2)

MALAVariate Type

Declaration

typealias MALAVariate SamplerVariate{MALATune}

Fields

  • value::Vector{Float64} : simulated values.
  • tune::MALATune : tuning parameters for the sampling algorithm.

Constructors

MALAVariate(x::AbstractVector{T<:Real}, epsilon::Real, logfgrad::Function)
MALAVariate(x::AbstractVector{T<:Real}, epsilon::Real, Sigma::Matrix{U<:Real}, logfgrad::Function)

Construct a MALAVariate object that stores simulated values and tuning parameters for MALA sampling.

Arguments

  • x : initial values.
  • epsilon : factor by which the drift and covariance matrix of the proposal distribution are scaled.
  • Sigma : covariance matrix for the multivariate normal proposal distribution. The covariance matrix is relative to the unconstrained parameter space, where candidate draws are generated. If omitted, the identity matrix is assumed.
  • logfgrad : function that takes a single DenseVector argument of parameter values at which to compute the log-transformed density (up to a normalizing constant) and gradient vector, and returns the respective results as a tuple.

Value

Returns a MALAVariate type object with fields set to the supplied x and tuning parameter values.

MALATune Type

Declaration

type MALATune <: SamplerTune

Fields

  • logfgrad::Nullable{Function} : function supplied to the constructor to compute the log-transformed density and gradient vector, or null if not supplied.
  • epsilon::Float64 : factor by which the drift and covariance matrix of the proposal distribution are scaled.
  • SigmaL::Union{UniformScaling{Int}, LowerTriangular{Float64}} : Cholesky factorization of the covariance matrix for the multivariate normal proposal distribution.