sampler_NUTS {nimbleHMC}R Documentation

No-U-Turn (NUTS) Hamiltonian Monte Carlo (HMC) Sampler

Description

The NUTS sampler implements No-U-Turn (NUTS) Hamiltonian Monte Carlo (HMC) sampling following the algorithm of version 2.32.2 of Stan. Internally, any posterior dimensions with bounded support are transformed, so sampling takes place on an unconstrained space. In contrast to standard HMC (Neal, 2011), the NUTS algorithm removes the tuning parameters of the leapfrog step size and the number of leapfrog steps, thus providing a sampling algorithm that can be used without hand-tuning or trial runs.

Usage

sampler_NUTS(model, mvSaved, target, control)

Arguments

model

An uncompiled nimble model object on which the MCMC will operate.

mvSaved

A nimble modelValues object to be used to store MCMC samples.

target

A character vector of node names on which the sampler will operate.

control

A named list that controls the precise behavior of the sampler. The default values for control list elements are specified in the setup code of the sampler. A description of the possible control list elements appear in the details section.

Details

The NUTS sampler accepts the following control list elements:

NaN values may be encountered in the course of the leapfrog procedure. In particular, when the stepsize (epsilon) is too large, the leapfrog procedure can step too far and arrive at an invalid region of parameter space, thus generating a NaN value in the likelihood evaluation or in the gradient calculation. These situation are handled by the sampler by rejecting the NaN value, and reducing the stepsize.

Value

A object of class 'sampler_NUTS'.

Author(s)

Perry de Valpine and Daniel Turek

References

Hoffman, Matthew D., and Gelman, Andrew (2014). The No-U-Turn Sampler: Adaptively setting path lengths in Hamiltonian Monte Carlo. Journal of Machine Learning Research, 15(1): 1593-1623.

Stan Development Team. 2023. Stan Modeling Language Users Guide and Reference Manual, 2.32.2. https://mc-stan.org.

Examples

code <- nimbleCode({
    b0 ~ dnorm(0, 0.001)
    b1 ~ dnorm(0, 0.001)
    sigma ~ dunif(0, 10000)
    for(i in 1:N) {
        mu[i] <- b0 + b1 * x[i]
        y[i] ~ dnorm(mu[i], sd = sigma)
    }
})

set.seed(0)
N <- 100
x <- rnorm(N)
y <- 1 + 0.3*x + rnorm(N)
constants <- list(N = N, x = x)
data <- list(y = y)
inits <- list(b0 = 1, b1 = 0.1, sigma = 1)

Rmodel <- nimbleModel(code, constants, data, inits, buildDerivs = TRUE)

conf <- configureMCMC(Rmodel, nodes = NULL)

conf$addSampler(target = c('b0', 'b1', 'sigma'), type = 'NUTS')

Rmcmc <- buildMCMC(conf)


[Package nimbleHMC version 0.2.2 Index]