abc {pomp}R Documentation

Approximate Bayesian computation

Description

The approximate Bayesian computation (ABC) algorithm for estimating the parameters of a partially-observed Markov process.

Usage

## S4 method for signature 'data.frame'
abc(
  data,
  Nabc = 1,
  proposal,
  scale,
  epsilon,
  probes,
  params,
  rinit,
  rprocess,
  rmeasure,
  dprior,
  ...,
  verbose = getOption("verbose", FALSE)
)

## S4 method for signature 'pomp'
abc(
  data,
  Nabc = 1,
  proposal,
  scale,
  epsilon,
  probes,
  ...,
  verbose = getOption("verbose", FALSE)
)

## S4 method for signature 'probed_pomp'
abc(data, probes, ..., verbose = getOption("verbose", FALSE))

## S4 method for signature 'abcd_pomp'
abc(
  data,
  Nabc,
  proposal,
  scale,
  epsilon,
  probes,
  ...,
  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.

Nabc

the number of ABC 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.

scale

named numeric vector of scales.

epsilon

ABC tolerance.

probes

a single probe or a list of one or more probes. A probe is simply a scalar- or vector-valued function of one argument that can be applied to the data array of a ‘pomp’. A vector-valued probe must always return a vector of the same size. A number of useful probes are provided with the package: see basic probes.

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.

rmeasure

simulator of the measurement model, 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 rmeasure=NULL removes the measurement model simulator. For more information, see rmeasure 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.

Running ABC

abc returns an object of class ‘abcd_pomp’. One or more ‘abcd_pomp’ objects can be joined to form an ‘abcList’ object.

Re-running ABC iterations

To re-run a sequence of ABC iterations, one can use the abc method on a ‘abcd_pomp’ object. By default, the same parameters used for the original ABC 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.

Continuing ABC iterations

One can continue a series of ABC iterations from where one left off using the continue method. A call to abc to perform Nabc=m iterations followed by a call to continue to perform Nabc=n iterations will produce precisely the same effect as a single call to abc to perform Nabc=m+n iterations. By default, all the algorithmic parameters are the same as used in the original call to abc. Additional arguments will override the defaults.

Methods

The following can be applied to the output of an abc operation:

abc

repeats the calculation, beginning with the last state

continue

continues the abc calculation

plot

produces a series of diagnostic plots

traces

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

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

References

J.-M. Marin, P. Pudlo, C. P. Robert, and R. J. Ryder. Approximate Bayesian computational methods. Statistics and Computing 22, 1167–1180, 2012. doi:10.1007/s11222-011-9288-2.

T. Toni and M. P. H. Stumpf. Simulation-based model selection for dynamical systems in systems and population biology. Bioinformatics 26, 104–110, 2010. doi:10.1093/bioinformatics/btp619.

T. Toni, D. Welch, N. Strelkowa, A. Ipsen, and M. P. H. Stumpf. Approximate Bayesian computation scheme for parameter inference and model selection in dynamical systems. Journal of the Royal Society Interface 6, 187–202, 2009. doi:10.1098/rsif.2008.0172.

See Also

More on methods based on summary statistics: basic_probes, nlf, probe(), probe_match, spect(), spect_match

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

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

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


[Package pomp version 5.11.0.0 Index]