buildAuxiliaryFilter {nimbleSMC} | R Documentation |
Create an auxiliary particle filter algorithm to estimate log-likelihood.
Description
Create an auxiliary particle filter algorithm for a given NIMBLE state space model.
Usage
buildAuxiliaryFilter(model, nodes, control = list())
Arguments
model |
A NIMBLE model object, typically representing a state space model or a hidden Markov model. |
nodes |
A character vector specifying the stochastic latent model nodes over which the particle filter will stochastically integrate to estimate the log-likelihood function. All provided nodes must be stochastic. Can be one of three forms: a variable name, in which case all elements in the variable are taken to be latent (e.g., 'x'); an indexed variable, in which case all indexed elements are taken to be latent (e.g., 'x[1:100]' or 'x[1:100, 1:2]'); or a vector of multiple nodes, one per time point, in increasing time order (e.g., c("x[1:2, 1]", "x[1:2, 2]", "x[1:2, 3]", "x[1:2, 4]")). |
control |
A list specifying different control options for the particle filter. Options are described in the details section below. |
Details
Each of the control()
list options are described in detail here:
- lookahead
The lookahead function used to calculate auxiliary weights. Can choose between
'mean'
and'simulate'
. Defaults to'simulate'
.- resamplingMethod
The type of resampling algorithm to be used within the particle filter. Can choose between
'default'
(which uses NIMBLE'srankSample()
function),'systematic'
,'stratified'
,'residual'
, and'multinomial'
. Defaults to'default'
. Resampling methods other than'default'
are currently experimental.- saveAll
Indicates whether to save state samples for all time points (
TRUE
), or only for the most recent time point (FALSE
)- smoothing
Decides whether to save smoothed estimates of latent states, i.e., samples from f(x[1:t]|y[1:t]) if
smoothing = TRUE
, or instead to save filtered samples from f(x[t]|y[1:t]) ifsmoothing = FALSE
.smoothing = TRUE
only works ifsaveAll = TRUE
.- timeIndex
An integer used to manually specify which dimension of the latent state variable indexes time. This need only be set if the number of time points is less than or equal to the size of the latent state at each time point.
- initModel
A logical value indicating whether to initialize the model before running the filtering algorithm. Defaults to TRUE.
The auxiliary particle filter modifies the bootstrap filter (buildBootstrapFilter
)
by adding a lookahead step to the algorithm: before propagating particles from one time
point to the next via the transition equation, the auxiliary filter calculates a weight
for each pre-propogated particle by predicting how well the particle will agree with the
next data point. These pre-weights are used to conduct an initial resampling step before
propagation.
The resulting specialized particle filter algorthm will accept a
single integer argument (m
, default 10,000), which specifies the number
of random \'particles\' to use for estimating the log-likelihood. The algorithm
returns the estimated log-likelihood value, and saves
unequally weighted samples from the posterior distribution of the latent
states in the mvWSamples
modelValues object, with corresponding logged weights in mvWSamples['wts',]
.
An equally weighted sample from the posterior can be found in the mvEWsamp
modelValues object.
The auxiliary particle filter uses a lookahead function to select promising particles before propagation. This function can eithre be the expected
value of the latent state at the next time point (lookahead = 'mean'
) or a simulation from the distribution of the latent state at the next time point (lookahead = 'simulate'
), conditioned on the current particle.
@section returnESS()
Method:
Calling the returnESS()
method of an auxiliary particle filter after that filter has been run()
for a given model will return a vector of ESS (effective
sample size) values, one value for each time point.
Author(s)
Nicholas Michaud
References
Pitt, M.K., and Shephard, N. (1999). Filtering via simulation: Auxiliary particle filters. Journal of the American Statistical Association 94(446): 590-599.
See Also
Other particle filtering methods:
buildBootstrapFilter
,
buildEnsembleKF
,
buildIteratedFilter2()
,
buildLiuWestFilter
Examples
## For illustration only.
exampleCode <- nimbleCode({
x0 ~ dnorm(0, var = 1)
x[1] ~ dnorm(.8 * x0, var = 1)
y[1] ~ dnorm(x[1], var = .5)
for(t in 2:10){
x[t] ~ dnorm(.8 * x[t-1], var = 1)
y[t] ~ dnorm(x[t], var = .5)
}
})
model <- nimbleModel(code = exampleCode, data = list(y = rnorm(10)),
inits = list(x0 = 0, x = rnorm(10)))
my_AuxF <- buildAuxiliaryFilter(model, 'x',
control = list(saveAll = TRUE, lookahead = 'mean'))
## Now compile and run, e.g.,
## Cmodel <- compileNimble(model)
## Cmy_AuxF <- compileNimble(my_AuxF, project = model)
## logLik <- Cmy_AuxF$run(m = 1000)
## ESS <- Cmy_AuxF$returnESS()
## aux_X <- as.matrix(Cmy_AuxF$mvEWSamples, 'x')