pmcmc: The particle Markov chain Metropolis-Hastings algorithm

pmcmcR Documentation

The particle Markov chain Metropolis-Hastings algorithm

Description

The Particle MCMC algorithm for estimating the parameters of a partially-observed Markov process. Running pmcmc causes a particle random-walk Metropolis-Hastings Markov chain algorithm to run for the specified number of proposals.

Usage

## S4 method for signature 'data.frame'
pmcmc(
  data,
  Nmcmc = 1,
  proposal,
  Np,
  params,
  rinit,
  rprocess,
  dmeasure,
  dprior,
  ...,
  verbose = getOption("verbose", FALSE)
)

## S4 method for signature 'pomp'
pmcmc(
  data,
  Nmcmc = 1,
  proposal,
  Np,
  ...,
  verbose = getOption("verbose", FALSE)
)

## S4 method for signature 'pfilterd_pomp'
pmcmc(
  data,
  Nmcmc = 1,
  proposal,
  Np,
  ...,
  verbose = getOption("verbose", FALSE)
)

## S4 method for signature 'pmcmcd_pomp'
pmcmc(data, Nmcmc, proposal, ..., verbose = getOption("verbose", FALSE))

Arguments

data

either a data frame holding the time series data, or an object of class ‘pomp’, i.e., the output of another pomp calculation. Internally, data will be coerced to an array with storage-mode double.

Nmcmc

The number of PMCMC iterations to perform.

proposal

optional function that draws from the proposal distribution. Currently, the proposal distribution must be symmetric for proper inference: it is the user's responsibility to ensure that it is. Several functions that construct appropriate proposal function are provided: see MCMC proposals for more information.

Np

the number of particles to use. This may be specified as a single positive integer, in which case the same number of particles will be used at each timestep. Alternatively, if one wishes the number of particles to vary across timesteps, one may specify Np either as a vector of positive integers of length

length(time(object,t0=TRUE))

or as a function taking a positive integer argument. In the latter case, Np(k) must be a single positive integer, representing the number of particles to be used at the k-th timestep: Np(0) is the number of particles to use going from timezero(object) to time(object)[1], Np(1), from timezero(object) to time(object)[1], and so on, while when T=length(time(object)), Np(T) is the number of particles to sample at the end of the time-series.

params

optional; named numeric vector of parameters. This will be coerced internally to storage mode double.

rinit

simulator of the initial-state distribution. This can be furnished either as a C snippet, an R function, or the name of a pre-compiled native routine available in a dynamically loaded library. Setting rinit=NULL sets the initial-state simulator to its default. For more information, see rinit specification.

rprocess

simulator of the latent state process, specified using one of the rprocess plugins. Setting rprocess=NULL removes the latent-state simulator. For more information, see rprocess specification for the documentation on these plugins.

dmeasure

evaluator of the measurement model density, specified either as a C snippet, an R function, or the name of a pre-compiled native routine available in a dynamically loaded library. Setting dmeasure=NULL removes the measurement density evaluator. For more information, see dmeasure specification.

dprior

optional; prior distribution density evaluator, specified either as a C snippet, an R function, or the name of a pre-compiled native routine available in a dynamically loaded library. For more information, see prior specification. Setting dprior=NULL resets the prior distribution to its default, which is a flat improper prior.

...

additional arguments are passed to pomp.

verbose

logical; if TRUE, diagnostic messages will be printed to the console.

Value

An object of class ‘pmcmcd_pomp’.

Methods

The following can be applied to the output of a pmcmc operation:

pmcmc

repeats the calculation, beginning with the last state

continue

continues the pmcmc calculation

plot

produces a series of diagnostic plots

filter_traj

extracts a random sample from the smoothing distribution

traces

produces an mcmc object, to which the various coda convergence diagnostics can be applied

Re-running PMCMC Iterations

To re-run a sequence of PMCMC iterations, one can use the pmcmc method on a ‘pmcmc’ object. By default, the same parameters used for the original PMCMC run are re-used (except for verbose, the default of which is shown above). If one does specify additional arguments, these will override the defaults.

Note for Windows users

Some Windows users report problems when using C snippets in parallel computations. These appear to arise when the temporary files created during the C snippet compilation process are not handled properly by the operating system. To circumvent this problem, use the cdir and cfile options to cause the C snippets to be written to a file of your choice, thus avoiding the use of temporary files altogether.

Author(s)

Edward L. Ionides, Aaron A. King, Sebastian Funk

References

\Andrieu

2010

See Also

More on pomp estimation algorithms: abc(), bsmc2(), estimation_algorithms, mif2(), nlf, pomp-package, probe_match, spect_match

More on sequential Monte Carlo methods: bsmc2(), cond_logLik(), eff_sample_size(), filter_mean(), filter_traj(), kalman, mif2(), pfilter(), pred_mean(), pred_var(), saved_states(), wpfilter()

More on full-information (i.e., likelihood-based) methods: bsmc2(), mif2(), pfilter(), wpfilter()

More on Markov chain Monte Carlo methods: abc(), proposals

More on Bayesian methods: abc(), bsmc2(), dprior(), prior_spec, rprior()


pomp documentation built on Sept. 13, 2024, 1:08 a.m.