zigzagHMC: Sample from a truncated Gaussian distribution

View source: R/zigzagHMC.R

zigzagHMCR Documentation

Sample from a truncated Gaussian distribution

Description

Generate MCMC samples from a d-dimensional truncated Gaussian distribution with element-wise truncations using the Zigzag Hamiltonian Monte Carlo sampler (Zigzag-HMC).

Usage

zigzagHMC(
  nSample,
  burnin = 0,
  mean,
  prec,
  lowerBounds,
  upperBounds,
  init = NULL,
  stepsize = NULL,
  nutsFlg = FALSE,
  precondition = FALSE,
  seed = NULL,
  diagnosticMode = FALSE
)

Arguments

nSample

number of samples after burn-in.

burnin

number of burn-in samples (default = 0).

mean

a d-dimensional mean vector.

prec

a d-by-d precision matrix of the Gaussian distribution.

lowerBounds

a d-dimensional vector specifying the lower bounds. -Inf is accepted.

upperBounds

a d-dimensional vector specifying the upper bounds. Inf is accepted.

init

a d-dimensional vector of the initial value. init must satisfy all constraints. If init = NULL, a random initial value will be used.

stepsize

step size for Zigzag-HMC or Zigzag-NUTS (if nutsFlg = TRUE). Default value is the empirically optimal choice: sqrt(2)(lambda)^(-1/2) for Zigzag-HMC and 0.1(lambda)^(-1/2) for Zigzag-NUTS, where lambda is the minimal eigenvalue of the precision matrix.

nutsFlg

logical. If TRUE the No-U-Turn sampler will be used (Zigzag-NUTS).

precondition

logical. If TRUE, the precision matrix will be preconditioned so that its diagonals (i.e. conditional variances) are all 1.

seed

random seed (default = 1).

diagnosticMode

logical. TRUE for also returning diagnostic information such as the stepsize used.

Value

an nSample-by-d matrix of samples. If diagnosticMode is TRUE, a list with additional diagnostic information is returned.

References

\insertRef

nishimura2024zigzaghdtg

\insertRef

nishimura2020discontinuoushdtg

Examples

set.seed(1)
d <- 10
A <- matrix(runif(d^2)*2-1, ncol=d)
covMat <- t(A) %*% A
precMat <- solve(covMat)
initial <- rep(1, d)
results <- zigzagHMC(nSample = 1000, burnin = 1000, mean = rep(0, d), prec = precMat,
lowerBounds = rep(0, d), upperBounds = rep(Inf, d))


hdtg documentation built on June 8, 2025, 1:07 p.m.

Related to zigzagHMC in hdtg...