Binary Metropolised Gibbs (BMG)¶
Implementation of the binary-state Metropolised Gibbs sampler described by Schafer [83][84] in which components are drawn sequentially from full conditional marginal distributions and accepted together in a single Metropolis-Hastings step. The sampler simulates autocorrelated draws from a distribution that can be specified up to a constant of proportionality.
Model-Based Constructor¶
-
BMG
(params::ElementOrVector{Symbol}; k::BMGForm=1)¶ Construct a
Sampler
object for BMG sampling. Parameters are assumed to have binary numerical values (0 or 1).Arguments
params
: stochastic node(s) to be updated with the sampler.k
: number of parameters or vector of parameter indices to select at random for simultaneous updating in each call of the sampler.
Value
Returns aSampler{BMGTune{typeof(k)}}
type object.Example
Stand-Alone Function¶
-
sample!
(v::SamplerVariate{BMGTune{F<:BMGForm}})¶ Draw one sample from a target distribution using the BMG sampler. Parameters are assumed to have binary numerical values (0 or 1).
Arguments
v
: current state of parameters to be simulated.
Value
Returnsv
updated with simulated values and associated tuning parameters.Example
################################################################################ ## Linear Regression ## y ~ MvNormal(X * (beta0 .* gamma), 1) ## gamma ~ DiscreteUniform(0, 1) ################################################################################ using Mamba ## Data n, p = 25, 10 X = randn(n, p) beta0 = randn(p) gamma0 = rand(0:1, p) y = X * (beta0 .* gamma0) + randn(n) ## Log-transformed Posterior(gamma) + Constant logf = function(gamma::DenseVector) logpdf(MvNormal(X * (beta0 .* gamma), 1.0), y) end ## MCMC Simulation with Binary Metropolised Gibbs t = 10000 sim1 = Chains(t, p, names = map(i -> "gamma[$i]", 1:p)) sim2 = Chains(t, p, names = map(i -> "gamma[$i]", 1:p)) gamma1 = BMGVariate(zeros(p), logf) gamma2 = BMGVariate(zeros(p), logf, k=Vector{Int}[[i] for i in 1:p]) for i in 1:t sample!(gamma1) sample!(gamma2) sim1[i, :, 1] = gamma1 sim2[i, :, 1] = gamma2 end describe(sim1) describe(sim2)
BMG Variate Type¶
Declaration¶
SamplerVariate{BMGTune{F<:BMGForm}}
Fields¶
value::Vector{Float64}
: simulated values.tune::BMGTune{F}
: tuning parameters for the sampling algorithm.
Constructor¶
-
BMGVariate
(x::AbstractVector{T<:Real}, logf::Function; k::BMGForm=1)¶ Construct a
SamplerVariate
object that stores simulated values and tuning parameters for BMG sampling.Arguments
x
: initial values.logf
: function that takes a singleDenseVector
argument of parameter values at which to compute the log-transformed density (up to a normalizing constant).k
: number of parameters or vector of parameter indices to select at random for simultaneous updating in each call of the sampler.
Value
Returns aSamplerVariate{BMGTune{typeof(k)}}
type object with fields set to the suppliedx
and tuning parameter values.
BMGTune Type¶
Declaration¶
typealias BMGForm Union{Int, Vector{Vector{Int}}}
type BMGTune{F<:BMGForm} <: SamplerTune
Fields¶
logf::Nullable{Function}
: function supplied to the constructor to compute the log-transformed density, or null if not supplied.k::F
: number of parameters or vector of parameter indices to select at random for simultaneous updating in each call of the sampler.