samplers: Particle Filtering MCMC Sampling Algorithms

SMCsamplersR Documentation

Particle Filtering MCMC Sampling Algorithms

Description

Details of the particle filtering MCMC sampling algorithms provided in nimbleSMC.

Usage

sampler_RW_PF(model, mvSaved, target, control)

sampler_RW_PF_block(model, mvSaved, target, control)

Arguments

model

(uncompiled) model on which the MCMC is to be run

mvSaved

modelValues object to be used to store MCMC samples

target

node(s) on which the sampler will be used

control

named list that controls the precise behavior of the sampler, with elements specific to samplertype. The default values for control list are specified in the setup code of each sampling algorithm. Descriptions of each sampling algorithm, and the possible customizations for each sampler (using the control argument) appear below.

RW_PF sampler

The particle filter sampler allows the user to perform particle MCMC (PMCMC) (Andrieu et al., 2010), primarily for state-space or hidden Markov models of time-series data. This method uses Metropolis-Hastings samplers for top-level parameters but uses the likelihood approximation of a particle filter (sequential Monte Carlo) to integrate over latent nodes in the time-series. The RW_PF sampler uses an adaptive Metropolis-Hastings algorithm with a univariate normal proposal distribution for a scalar parameter. Note that samples of the latent states can be retained as well, but the top-level parameter being sampled must be a scalar. A bootstrap, auxiliary, or user defined particle filter can be used to integrate over latent states.

For more information about user-defined samplers within a PMCMC sampler, see the NIMBLE User Manual.

The RW_PF sampler accepts the following control list elements:

  • adaptive. A logical argument, specifying whether the sampler should adapt the scale (proposal standard deviation) throughout the course of MCMC execution to achieve a theoretically desirable acceptance rate. (default = TRUE)

  • adaptInterval. The interval on which to perform adaptation. Every adaptInterval MCMC iterations (prior to thinning), the RW sampler will perform its adaptation procedure. This updates the scale variable, based upon the sampler's achieved acceptance rate over the past adaptInterval iterations. (default = 200)

  • scale. The initial value of the normal proposal standard deviation. If adaptive = FALSE, scale will never change. (default = 1)

  • pfNparticles. The number of particles to use in the approximation to the log likelihood of the data (default = 1000).

  • latents. Character vector specifying the nodes that are latent states over which the particle filter will operate to approximate the log-likelihood function.

  • pfType. Character argument specifying the type of particle filter that should be used for likelihood approximation. Choose from "bootstrap" and "auxiliary". Defaults to "bootstrap".

  • pfControl. A control list that is passed to the particle filter function. For details on control lists for bootstrap or auxiliary particle filters, see buildBootstrapFilter or buildAuxiliaryFilter respectively. Additionally, this can be used to pass custom arguments into a user-defined particle filter.

  • pfOptimizeNparticles. A logical argument, specifying whether to use an experimental procedure to automatically determine the optimal number of particles to use, based on Pitt and Shephard (2011). This will override any value of pfNparticles specified above.

  • pf. A user-defined particle filter object, if a bootstrap or auxiliary particle filter is not adequate.

RW_PF_block sampler

The particle filter block sampler allows the user to perform particle MCMC (PMCMC) (Andrieu et al., 2010) for multiple parameters jointly, primarily for state-space or hidden Markov models of time-series data. This method uses Metropolis-Hastings block samplers for top-level parameters but uses the likelihood approximation of a particle filter (sequential Monte Carlo) to integrate over latent nodes in the time-series. The RW_PF sampler uses an adaptive Metropolis-Hastings algorithm with a multivariate normal proposal distribution. Note that samples of the latent states can be retained as well, but the top-level parameter being sampled must be a scalar. A bootstrap, auxiliary, or user defined particle filter can be used to integrate over latent states.

For more information about user-defined samplers within a PMCMC sampler, see the NIMBLE User Manual.

The RW_PF_block sampler accepts the following control list elements:

  • adaptive. A logical argument, specifying whether the sampler should adapt the proposal covariance throughout the course of MCMC execution. (default = TRUE)

  • adaptScaleOnly. A logical argument, specifying whether adaptation should be done only for scale (TRUE) or also for provCov (FALSE). This argument is only relevant when adaptive = TRUE. When adaptScaleOnly = FALSE, both scale and propCov undergo adaptation; the sampler tunes the scaling to achieve a theoretically good acceptance rate, and the proposal covariance to mimic that of the empirical samples. When adaptScaleOnly = TRUE, only the proposal scale is adapted. (default = FALSE)

  • adaptInterval. The interval on which to perform adaptation. (default = 200)

  • scale. The initial value of the scalar multiplier for propCov. If adaptive = FALSE, scale will never change. (default = 1)

  • adaptFactorExponent. Exponent controling the rate of decay of the scale adaptation factor. See Shaby and Wells, 2011, for details. (default = 0.8)

  • propCov. The initial covariance matrix for the multivariate normal proposal distribution. This element may be equal to the 'identity', in which case the identity matrix of the appropriate dimension will be used for the initial proposal covariance matrix. (default is 'identity')

  • pfNparticles. The number of particles to use in the approximation to the log likelihood of the data (default = 1000).

  • latents. Character vector specifying the nodes that are latent states over which the particle filter will operate to approximate the log-likelihood function.

  • pfType. Character argument specifying the type of particle filter that should be used for likelihood approximation. Choose from "bootstrap" and "auxiliary". Defaults to "bootstrap".

  • pfControl. A control list that is passed to the particle filter function. For details on control lists for bootstrap or auxiliary particle filters, see buildBootstrapFilter or buildAuxiliaryFilter respectively. Additionally, this can be used to pass custom arguments into a user defined particle filter.

  • pfOptimizeNparticles. A logical argument, specifying whether to automatically determine the optimal number of particles to use, based on Pitt and Shephard (2011). This will override any value of pfNparticles specified above.

  • pf. A user-defined particle filter object, if a bootstrap or auxiliary particle filter is not adequate.


nimbleSMC documentation built on Sept. 11, 2023, 1:07 a.m.