bsmc2 {pomp}R Documentation

The Liu and West Bayesian particle filter

Description

Modified version of the Liu & West (2001) algorithm.

Usage

## S4 method for signature 'data.frame'
bsmc2(
  data,
  Np,
  smooth = 0.1,
  params,
  rprior,
  rinit,
  rprocess,
  dmeasure,
  partrans,
  ...,
  verbose = getOption("verbose", FALSE)
)

## S4 method for signature 'pomp'
bsmc2(data, Np, smooth = 0.1, ..., 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.

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.

smooth

Kernel density smoothing parameter. The compensating shrinkage factor will be sqrt(1-smooth^2). Thus, smooth=0 means that no noise will be added to parameters. The general recommendation is that the value of smooth should be chosen close to 0 (e.g., shrink ~ 0.1).

params

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

rprior

optional; prior distribution sampler, 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 rprior=NULL removes the prior distribution sampler.

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.

partrans

optional parameter transformations, constructed using parameter_trans.

Many algorithms for parameter estimation search an unconstrained space of parameters. When working with such an algorithm and a model for which the parameters are constrained, it can be useful to transform parameters. One should supply the partrans argument via a call to parameter_trans. For more information, see parameter_trans. Setting partrans=NULL removes the parameter transformations, i.e., sets them to the identity transformation.

...

additional arguments are passed to pomp.

verbose

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

Details

bsmc2 uses a version of the original algorithm (Liu & West 2001), but discards the auxiliary particle filter. The modification appears to give superior performance for the same amount of effort.

Samples from the prior distribution are drawn using the rprior component. This is allowed to depend on elements of params, i.e., some of the elements of params can be treated as “hyperparameters”. Np draws are made from the prior distribution.

Value

An object of class ‘bsmcd_pomp’. The following methods are avaiable:

plot

produces diagnostic plots

as.data.frame

puts the prior and posterior samples into a data frame

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)

Michael Lavine, Matthew Ferrari, Aaron A. King, Edward L. Ionides

References

J. Liu and M. West. Combining parameter and state estimation in simulation-based filtering. In A. Doucet, N. de Freitas, and N. J. Gordon, (eds.), Sequential Monte Carlo Methods in Practice, pp. 197–224. Springer, New York, 2001. doi:10.1007/978-1-4757-3437-9_10.

See Also

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

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

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

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


[Package pomp version 5.11.0.0 Index]