pfilter function

Particle filter

Particle filter

A plain vanilla sequential Monte Carlo (particle filter) algorithm. Resampling is performed at each observation.

## S4 method for signature 'data.frame' pfilter( data, ..., Np, params, rinit, rprocess, dmeasure, pred.mean = FALSE, pred.var = FALSE, filter.mean = FALSE, filter.traj = FALSE, save.states = c("no", "filter", "prediction", "weighted", "unweighted", "FALSE", "TRUE"), verbose = getOption("verbose", FALSE) ) ## S4 method for signature 'pomp' pfilter( data, ..., Np, pred.mean = FALSE, pred.var = FALSE, filter.mean = FALSE, filter.traj = FALSE, save.states = c("no", "filter", "prediction", "weighted", "unweighted", "FALSE", "TRUE"), verbose = getOption("verbose", FALSE) ) ## S4 method for signature 'pfilterd_pomp' pfilter(data, ..., Np, verbose = getOption("verbose", FALSE)) ## S4 method for signature 'objfun' pfilter(data, ...)

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.

  • ...: additional arguments are passed to pomp. This allows one to set, unset, or modify basic model components within a call to this function.

  • 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 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 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 .

  • pred.mean: logical; if TRUE, the prediction means are calculated for the state variables and parameters.

  • pred.var: logical; if TRUE, the prediction variances are calculated for the state variables and parameters.

  • filter.mean: logical; if TRUE, the filtering means are calculated for the state variables and parameters.

  • filter.traj: logical; if TRUE, a filtered trajectory is returned for the state variables and parameters. See filter_traj for more information.

  • save.states: character; If save.states="no" (the default), information on the latent states is not saved. If save.states="filter", the state-vector for each filtered particle Xn,jFX_{n,j}^F at each time nn is saved. If save.states="prediction", the state-vector for each prediction particle Xn,jPX_{n,j}^P at each time nn is saved, along with the corresponding weight wn,j=fYnXn(yXn,jP;θ)w_{n,j} = f_{Y_n|X_n}(y^*|X_{n, j}^P;\theta). The options "unweighted", "weighted", TRUE, and FALSE are deprecated and will issue a warning if used, mapping to the new values for backward compatibility. The options "unweighted" and TRUE are synonymous with "filter"; the option "weighted" is synonymous with "prediction"; the option FALSE is synonymous with "no". To retrieve the saved states, apply saved_states to the result of the pfilter computation.

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

Returns

An object of class pfilterd_pomp , which extends class pomp . Information can be extracted from this object using the methods documented below.

Methods

  • logLik: the estimated log likelihood
  • cond_logLik: the estimated conditional log likelihood
  • eff_sample_size: the (time-dependent) estimated effective sample size
  • pred_mean, pred_var: the mean and variance of the approximate prediction distribution
  • filter_mean: the mean of the filtering distribution
  • filter_traj: retrieve one particle trajectory. Useful for building up the smoothing distribution.
  • saved_states: retrieve saved states
  • as.data.frame: coerce to a data frame
  • plot: diagnostic plots

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.

Examples

pf <- pfilter(gompertz(),Np=1000) ## use 1000 particles plot(pf) logLik(pf) cond_logLik(pf) ## conditional log-likelihoods eff_sample_size(pf) ## effective sample size logLik(pfilter(pf)) ## run it again with 1000 particles ## run it again with 2000 particles pf <- pfilter(pf,Np=2000,filter.mean=TRUE,filter.traj=TRUE,save.states="filter") fm <- filter_mean(pf) ## extract the filtering means ft <- filter_traj(pf) ## one draw from the smoothing distribution ss <- saved_states(pf,format="d") ## the latent-state portion of each particle as(pf,"data.frame") |> head()

References

M.S. Arulampalam, S. Maskell, N. Gordon, and T. Clapp. A tutorial on particle filters for online nonlinear, non-Gaussian Bayesian tracking. IEEE Transactions on Signal Processing 50 , 174--188, 2002. tools:::Rd_expr_doi("10.1109/78.978374") .

A. Bhadra and E.L. Ionides. Adaptive particle allocation in iterated sequential Monte Carlo via approximating meta-models. Statistics and Computing 26 , 393--407, 2016. tools:::Rd_expr_doi("10.1007/s11222-014-9513-x") .

See Also

More on pomp elementary algorithms: elementary_algorithms, kalman, pomp-package, probe(), simulate(), spect(), trajectory(), wpfilter()

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

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

Author(s)

Aaron A. King

  • Maintainer: Aaron A. King
  • License: GPL-3
  • Last published: 2025-01-08