Sequential Monte Carlo

SMCUpdater - SMC-Based Particle Updater

Class Reference

class qinfer.SMCUpdater(model, n_particles, prior, resample_a=None, resampler=None, resample_thresh=0.5, debug_resampling=False, track_resampling_divergence=False, zero_weight_policy=u'error', zero_weight_thresh=None, canonicalize=True)[source]

Bases: qinfer.distributions.ParticleDistribution

Creates a new Sequential Monte carlo updater, using the algorithm of [GFWC12].

Parameters:
  • model (Model) – Model whose parameters are to be inferred.
  • n_particles (int) – The number of particles to be used in the particle approximation.
  • prior (Distribution) – A representation of the prior distribution.
  • resampler (callable) – Specifies the resampling algorithm to be used. See Resampling Algorithms for more details.
  • resample_thresh (float) – Specifies the threshold for \(N_{\text{ess}}\) to decide when to resample.
  • debug_resampling (bool) – If True, debug information will be generated on resampling performance, and will be written to the standard Python logger.
  • track_resampling_divergence (bool) – If true, then the divergences between the pre- and post-resampling distributions are tracked and recorded in the resampling_divergences attribute.
  • zero_weight_policy (str) – Specifies the action to be taken when the particle weights would all be set to zero by an update. One of ["ignore", "skip", "warn", "error", "reset"].
  • zero_weight_thresh (float) – Value to be used when testing for the zero-weight condition.
  • canonicalize (bool) – If True, particle locations will be updated to canonical locations as described by the model class after each prior sampling and resampling.
resample_count

Returns the number of times that the updater has resampled the particle approximation.

Type:int
just_resampled

True if and only if there has been no data added since the last resampling, or if there has not yet been a resampling step.

Type:bool
normalization_record

Returns the normalization record.

Type:float
log_total_likelihood

Returns the log-likelihood of all the data collected so far.

Equivalent to:

np.sum(np.log(updater.normalization_record))
Type:float
min_n_ess

Returns the smallest effective sample size (ESS) observed in the history of this updater.

Type:float
Returns:The minimum of observed effective sample sizes as reported by n_ess.
data_record

List of outcomes given to update().

Type:list of int
resampling_divergences

List of KL divergences between the pre- and post-resampling distributions, if that is being tracked. Otherwise, None.

Type:list of float or None
reset(n_particles=None, only_params=None, reset_weights=True)[source]

Causes all particle locations and weights to be drawn fresh from the initial prior.

Parameters:
  • n_particles (int) – Forces the size of the new particle set. If None, the size of the particle set is not changed.
  • only_params (slice) – Resets only some of the parameters. Cannot be set if n_particles is also given.
  • reset_weights (bool) – Resets the weights as well as the particles.
hypothetical_update(outcomes, expparams, return_likelihood=False, return_normalization=False)[source]

Produces the particle weights for the posterior of a hypothetical experiment.

Parameters:
  • outcomes (int or an ndarray of dtype int.) – Integer index of the outcome of the hypothetical experiment.
  • expparams (numpy.ndarray) – Experiments to be used for the hypothetical updates.
  • weights (ndarray, shape (n_outcomes, n_expparams, n_particles)) – Weights assigned to each particle in the posterior distribution \(\Pr(\omega | d)\).
update(outcome, expparams, check_for_resample=True)[source]

Given an experiment and an outcome of that experiment, updates the posterior distribution to reflect knowledge of that experiment.

After updating, resamples the posterior distribution if necessary.

Parameters:
  • outcome (int) – Label for the outcome that was observed, as defined by the Model instance under study.
  • expparams (ndarray of dtype given by the expparams_dtype property of the underlying model) – Parameters describing the experiment that was performed.
  • check_for_resample (bool) – If True, after performing the update, the effective sample size condition will be checked and a resampling step may be performed.
batch_update(outcomes, expparams, resample_interval=5)[source]

Updates based on a batch of outcomes and experiments, rather than just one.

Parameters:
  • outcomes (numpy.ndarray) – An array of outcomes of the experiments that were performed.
  • expparams (numpy.ndarray) – Either a scalar or record single-index array of experiments that were performed.
  • resample_interval (int) – Controls how often to check whether \(N_{\text{ess}}\) falls below the resample threshold.
resample()[source]

Forces the updater to perform a resampling step immediately.

bayes_risk(expparams)[source]

Calculates the Bayes risk for hypothetical experiments, assuming the quadratic loss function defined by the current model’s scale matrix (see qinfer.abstract_model.Simulatable.Q).

Parameters:expparams (ndarray of dtype given by the current model’s expparams_dtype property, and of shape (1,)) – The experiments at which to compute the risk.
Return np.ndarray:
 The Bayes risk for the current posterior distribution at each hypothetical experiment in expparams, therefore has shape (expparams.size,)
expected_information_gain(expparams)[source]

Calculates the expected information gain for each hypothetical experiment.

Parameters:expparams (ndarray of dtype given by the current model’s expparams_dtype property, and of shape (n,)) – The experiments at which to compute expected information gain.
Return float:The expected information gain for each hypothetical experiment in expparams.
risk(x0)[source]
posterior_marginal(idx_param=0, res=100, smoothing=0, range_min=None, range_max=None)[source]

Returns an estimate of the marginal distribution of a given model parameter, based on taking the derivative of the interpolated cdf.

Parameters:
  • idx_param (int) – Index of parameter to be marginalized.
  • res1 (int) – Resolution of of the axis.
  • smoothing (float) – Standard deviation of the Gaussian kernel used to smooth; same units as parameter.
  • range_min (float) – Minimum range of the output axis.
  • range_max (float) – Maximum range of the output axis.
plot_posterior_marginal(idx_param=0, res=100, smoothing=0, range_min=None, range_max=None, label_xaxis=True, other_plot_args={}, true_model=None)[source]

Plots a marginal of the requested parameter.

Parameters:
  • idx_param (int) – Index of parameter to be marginalized.
  • res1 (int) – Resolution of of the axis.
  • smoothing (float) – Standard deviation of the Gaussian kernel used to smooth; same units as parameter.
  • range_min (float) – Minimum range of the output axis.
  • range_max (float) – Maximum range of the output axis.
  • label_xaxis (bool) – Labels the \(x\)-axis with the model parameter name given by this updater’s model.
  • other_plot_args (dict) – Keyword arguments to be passed to matplotlib’s plot function.
  • true_model (np.ndarray) – Plots a given model parameter vector as the “true” model for comparison.
plot_covariance(corr=False, param_slice=None, tick_labels=None, tick_params=None)[source]

Plots the covariance matrix of the posterior as a Hinton diagram.

Note

This function requires that mpltools is installed.

Parameters:
  • corr (bool) – If True, the covariance matrix is first normalized by the outer product of the square root diagonal of the covariance matrix such that the correlation matrix is plotted instead.
  • param_slice (slice) – Slice of the modelparameters to be plotted.
  • tick_labels (list) – List of tick labels for each component; by default, these are drawn from the model itself.
posterior_mesh(idx_param1=0, idx_param2=1, res1=100, res2=100, smoothing=0.01)[source]

Returns a mesh, useful for plotting, of kernel density estimation of a 2D projection of the current posterior distribution.

Parameters:
  • idx_param1 (int) – Parameter to be treated as \(x\) when plotting.
  • idx_param2 (int) – Parameter to be treated as \(y\) when plotting.
  • res1 (int) – Resolution along the \(x\) direction.
  • res2 (int) – Resolution along the \(y\) direction.
  • smoothing (float) – Standard deviation of the Gaussian kernel used to smooth the particle approximation to the current posterior.
plot_posterior_contour(idx_param1=0, idx_param2=1, res1=100, res2=100, smoothing=0.01)[source]

Plots a contour of the kernel density estimation of a 2D projection of the current posterior distribution.

Parameters:
  • idx_param1 (int) – Parameter to be treated as \(x\) when plotting.
  • idx_param2 (int) – Parameter to be treated as \(y\) when plotting.
  • res1 (int) – Resolution along the \(x\) direction.
  • res2 (int) – Resolution along the \(y\) direction.
  • smoothing (float) – Standard deviation of the Gaussian kernel used to smooth the particle approximation to the current posterior.