sampler_NUTS | R Documentation |
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.
sampler_NUTS(model, mvSaved, target, control)
model |
An uncompiled nimble model object on which the MCMC will operate. |
mvSaved |
A nimble |
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. |
The NUTS sampler accepts the following control list elements:
messages. A logical argument, specifying whether to print informative messages (default = TRUE)
numWarnings. A numeric argument, specifying how many warnings messages to emit (for example, when NaN
values are encountered). See additional details below. (default = 0)
epsilon. A positive numeric argument, specifying the initial step-size value. If not provided, an appropriate initial value is selected.
gamma. A positive numeric argument, specifying the degree of shrinkage used during the initial period of step-size adaptation. (default = 0.05)
t0. A non-negative numeric argument, where larger values stabilize (attenuate) the initial period of step-size adaptation. (default = 10)
kappa. A numeric argument between zero and one, where smaller values give a higher weighting to more recent iterations during the initial period of step-size adaptation. (default = 0.75)
delta. A numeric argument, specifying the target acceptance probability used during the initial period of step-size adaptation. (default = 0.8)
deltaMax. A positive numeric argument, specifying the maximum allowable divergence from the Hamiltonian value. Paths which exceed this value are considered divergent, and will not proceed further. (default = 1000)
M. A vector of positive real numbers, with length equal to the number of dimensions being sampled. Elements of M
specify the diagonal elements of the diagonal mass matrix (or the metric) used for the auxiliary momentum variables in sampling. Sampling may be improved if the elements of M
approximate the marginal inverse variance (precision) of the (potentially transformed) parameters. (default: a vector of ones).
warmupMode. A character string, specifying the behavior for choosing the number of warmup iterations. Four values are possible. The value 'default' (the default) sets the number of warmup iterations as the number of burnin iterations (if a positive value for nburnin
is used) or half the number of MCMC iterations in each chain (if nburnin = 0
). The value 'burnin' sets the number of warmup iterations as the number of burnin iterations regardless of the length of the burnin period. The value 'fraction' sets the number of warmup iterations as fraction*niter
, where fraction
is the value of the warmup
control argument, and niter
is the number of MCMC iterations in each chain; in this case, the value of the warmup
control argument must be between 0 and 1. The value 'iterations' sets the number of warmup iterations as the value of the warmup
control argumnet, regardless of the length of the burnin period or the number of MCMC iterations; in this case the value of warmup
must be a non-negative integer. In all cases, the number of (pre-thinning) samples discarded equals nburnin
, as is always the case for MCMC in NIMBLE.
warmup. Numeric value used in determining the number of warmup iterations. This control argument is only used when warmupMode
is 'fraction' or 'iterations'.
maxTreeDepth. The maximum allowable depth of the binary leapfrog search tree for generating candidate transitions. (default = 10)
adaptWindow. Number of iterations in the first adaptation window used for adapting the mass matrix (M). Subsequent adaptation windows double in length, so long as enough warmup iterations are available. (default = 25)
initBuffer. Number of iterations in the initial warmup window, which occurs prior to the first adaptation of the metric M. (default = 75)
termBuffer. Number of iterations in the final (terminal) warmup window, before which the metric M is not adjusted(default = 50)
adaptive. A logical argument, specifying whether to do any adaptation whatsoever. When TRUE
, specific adaptation routines are controlled by the adaptEpsilon
and adaptM
control list elements. (default = TRUE)
adaptEpsilon. A logical argument, specifying whether to perform stepsize adaptation. Only used when adaptive = TRUE
. (default = TRUE)
adaptM. A logical argument, specifying whether to perform adaptation of the mass matrix (metric) M. Only used when adaptive = TRUE
. (default = TRUE)
initializeEpsilon. A logical argument, specifying whether to perform the epsilon (stepsize) initialization routine at the onset of each adaptation window. (default = TRUE)
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.
A object of class 'sampler_NUTS'.
Perry de Valpine and Daniel Turek
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.
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)
Add the following code to your website.
For more information on customizing the embed code, read Embedding Snippets.