Particle Filtering
Gen.initialize_particle_filter
— Functionstate = initialize_particle_filter(model::GenerativeFunction, model_args::Tuple,
observations::ChoiceMap, proposal::GenerativeFunction, proposal_args::Tuple,
num_particles::Int)
Initialize the state of a particle filter using a custom proposal for the initial latent state.
state = initialize_particle_filter(model::GenerativeFunction, model_args::Tuple,
observations::ChoiceMap, num_particles::Int)
Initialize the state of a particle filter, using the default proposal for the initial latent state.
Gen.particle_filter_step!
— Function(log_incremental_weights,) = particle_filter_step!(
state::ParticleFilterState, new_args::Tuple, argdiffs,
observations::ChoiceMap, proposal::GenerativeFunction, proposal_args::Tuple)
Perform a particle filter update, where the model arguments are adjusted, new observations are added, and some combination of a custom proposal and the model's internal proposal is used for proposing new latent state. That is, for each particle,
- The proposal function
proposal
is evaluated with argumentsTuple(t_old, proposal_args...)
(wheret_old
is the old model trace), and produces its own trace (call itproposal_trace
); and - The old model trace is replaced by a new model trace (call it
t_new
).
The choicemap of t_new
satisfies the following conditions:
get_choices(t_old)
is a subset ofget_choices(t_new)
;observations
is a subset ofget_choices(t_new)
;get_choices(proposal_trace)
is a subset ofget_choices(t_new)
.
Here, when we say one choicemap a
is a "subset" of another choicemap b
, we mean that all keys that occur in a
also occur in b
, and the values at those addresses are equal.
It is an error if no trace t_new
satisfying the above conditions exists in the support of the model (with the new arguments). If such a trace exists, then the random choices not determined by the above requirements are sampled using the internal proposal.
(log_incremental_weights,) = particle_filter_step!(
state::ParticleFilterState, new_args::Tuple, argdiffs,
observations::ChoiceMap)
Perform a particle filter update, where the model arguments are adjusted, new observations are added, and the default proposal is used for new latent state.
Gen.maybe_resample!
— Functiondid_resample::Bool = maybe_resample!(state::ParticleFilterState;
ess_threshold::Float64=length(state.traces)/2, verbose=false)
Do a resampling step if the effective sample size is below the given threshold. Return true
if a resample thus occurred, false
otherwise.
Gen.log_ml_estimate
— Functionestimate = log_ml_estimate(state::ParticleFilterState)
Return the particle filter's current estimate of the log marginal likelihood.
Gen.get_traces
— Functiontraces = get_traces(state::ParticleFilterState)
Return the vector of traces in the current state, one for each particle.
Gen.get_log_weights
— Functionlog_weights = get_log_weights(state::ParticleFilterState)
Return the vector of log weights for the current state, one for each particle.
The weights are not normalized, and are in log-space.
Gen.sample_unweighted_traces
— Functiontraces::Vector = sample_unweighted_traces(state::ParticleFilterState, num_samples::Int)
Sample a vector of num_samples
traces from the weighted collection of traces in the given particle filter state.