Interfaces

Particle Filter

class pyparticleest.interfaces.SIR
copy_ind(particles, new_ind=None)

Copy select particles, can be overriden for models that require special handling of the particle representations when copying them

Args:

  • particles (array-like): Model specific representation of all particles, with first dimension = N (number of particles)
  • new_ind (array-like): Array of ints, specifying indices to copy
Returns:
(array-like) with first dimension = len(new_ind)
create_initial_estimate(N)

Sample particles from initial distribution

Args:
  • N (int): Number of particles to sample
Returns:
(array-like) with first dimension = N, model specific representation of all particles
class pyparticleest.interfaces.ParticleFiltering

Base class for particles to be used with particle filtering. particles are a model specific array where the first dimension indexes the different particles.

measure(particles, y, t)

Return the log-pdf value of the measurement

Args:

  • particles (array-like): Model specific representation of all particles, with first dimension = N (number of particles)
  • y (array-like): measurement
  • t (float): time-stamp
Returns:
(array-like) with first dimension = N, logp(y|x^i)
update(particles, u, t, noise)

Propagate estimate forward in time

Args:

  • particles (array-like): Model specific representation of all particles, with first dimension = N (number of particles)
  • u (array-like): input signal
  • t (float): time-stamp
  • noise (array-like): noise realization used for the calucations , with first dimension = N (number of particles)
Returns:
(array-like) with first dimension = N, particle estimate at time t+1
class pyparticleest.interfaces.FFProposeFromMeasure
propose_from_y(N, y, t)

Create N particles from p(x_t|y_t)

Particle Filter Non-Markov

class pyparticleest.interfaces.ParticleFilteringNonMarkov
cond_predict_single_step(part, past_trajs, pind, future_parts, find, ut, yt, tt, cur_ind)

Propagate states in ‘part’ conditioned on that the future state is ‘future_parts’. This is used for e.g. Rao-Blackwellized MHIPS, where we need to propagate forward in time conditioned on the nonlinear state, but we want to recompute the additional data stored, e.g to exclude measurements present in the sufficient statistics for future_parts.

Args:
  • part (array-like): Model specific representation of all particles, with first dimension = N (number of particles)
  • ptraj: array of trajectory step objects from previous time-steps, last index is step just before the current
  • anc (array-like): index of the ancestor of each particle in part
  • future_trajs (array-like): particle estimate for {t+1:T}
  • find (array-like): index in future_trajs corresponding to each particle in part
  • ut (array-like): input signals for {0:T}
  • yt (array-like): measurements for {0:T}
  • tt (array-like): time stamps for {0:T}
  • cur_ind (int): index of current timestep (in ut, yt and tt)
cond_sampled_initial(part, t)

Sample from initial distribution conditioned on the states being ‘part’ This is used for e.g. Rao-Blackwellized MHIPS, where we need to recompute the sufficient statistics without being affected by the intial measurement

Args: part: particles t: time-step

copy_ind(particles, new_ind=None)

Copy select particles, can be overriden for models that require special handling of the particle representations when copying them

Args:

  • particles (array-like): Model specific representation of all particles, with first dimension = N (number of particles)
  • new_ind (array-like): Array of ints, specifying indices to copy
Returns:
(array-like) with first dimension = len(new_ind)
create_initial_estimate(N)

Sample particles from initial distribution

Args:
  • N (int): Number of particles to sample
Returns:
(array-like) with first dimension = N, model specific representation of all particles
measure_full(particles, traj, uvec, yvec, tvec, ancestors)

Return the log-pdf value of the measurement

Args:
  • particles (array-like): Model specific representation of all particles, with first dimension = N (number of particles)
  • traj: array of trajectory step objects from previous time-steps, last index is step just before the current
  • ancestors (array-like): index of the ancestor of each particle in part
  • uvec (array-like): input signals for {0:t}
  • yvec (array-like): measurements for {0:t}
  • tvec (array-like): time stamps for {0:t}
Returns:
(array-like) with first dimension = N
sample_process_noise_full(ptraj, ancestors, ut, tt)

Sample process noise

Args:
  • ptraj: array of trajectory step objects from previous time-steps, last index is step just before the current
  • ancestors (array-like): index of the ancestor of each particle in part
  • ut (array-like): input signals for {0:T}
  • tt (array-like): time stamps for {0:T}
Returns:
(array-like) with first dimension = N
sample_smooth(part, ptraj, anc, future_trajs, find, ut, yt, tt, cur_ind)

Create sampled estimates for the smoothed trajectory. Allows the update representation of the particles used in the forward step to include additional data in the backward step, can also for certain models be used to update the points estimates based on the future information.

Default implementation uses the same format as forward in time it ss part of the ParticleFiltering interface since it is used also when calculating “ancestor” trajectories

Args:
  • part (array-like): Model specific representation of all particles, with first dimension = N (number of particles)
  • ptraj: array of trajectory step objects from previous time-steps, last index is step just before the current
  • anc (array-like): index of the ancestor of each particle in part
  • future_trajs (array-like): particle estimate for {t+1:T}
  • find (array-like): index in future_trajs corresponding to each particle in part
  • ut (array-like): input signals for {0:T}
  • yt (array-like): measurements for {0:T}
  • tt (array-like): time stamps for {0:T}
  • cur_ind (int): index of current timestep (in ut, yt and tt)
Returns:
(array-like) with first dimension = N
update_full(particles, traj, uvec, yvec, tvec, ancestors, noise)

Propagate estimate forward in time

Args:
  • particles (array-like): Model specific representation of all particles, with first dimension = N (number of particles)
  • traj: array of trajectory step objects from previous time-steps, last index is step just before the current
  • ancestors (array-like): index of the ancestor of each particle in part
  • uvec (array-like): input signals for {0:t}
  • yvec (array-like): measurements for {0:t}
  • tvec (array-like): time stamps for {0:t}
  • noise (array-like): samples noise for time t
Returns:
(array-like) with first dimension = N

Auxiliary Particle Filter

class pyparticleest.interfaces.AuxiliaryParticleFiltering

Base class for particles to be used with auxiliary particle filtering

eval_1st_stage_weights(particles, u, y, t)

Evaluate “first stage weights” for the auxiliary particle filter. (log-probability of measurement using some propagated statistic, such as the mean, for the future state)

Args:

  • particles (array-like): Model specific representation of all particles, with first dimension = N (number of particles)
  • u (array-like): input signal
  • y (array-like): measurement
  • t (float): time-stamp
Returns:
(array-like) with first dimension = N, logp(y_{t+1}|hat{x}_{t+1|t}^i)

Forward-Filter Backward Simulator

class pyparticleest.interfaces.FFBSi

Base class for particles to be used with particle smoothing (Backward Simulation)

logp_xnext(particles, next_part, u, t)

Return the log-pdf value for the possible future state ‘next’ given input u

Args:

  • particles (array-like): Model specific representation of all particles, with first dimension = N (number of particles)
  • next_part (array-like): particle estimate for t+1
  • u (array-like): input signal
  • t (float): time stamps
Returns:
(array-like) with first dimension = N, logp(x_{t+1}|x_t^i)
logp_xnext_full(part, past_trajs, pind, future_trajs, find, ut, yt, tt, cur_ind)

Return the log-pdf value for the entire future trajectory. Useful for non-markovian modeles, that result from e.g marginalized state-space models.

Default implemention just calls logp_xnext which is enough for Markovian models

Args:
  • part (array-like): Model specific representation of all particles, with first dimension = N (number of particles)
  • past_trajs: array of trajectory step objects from previous time-steps, last index is step just before the current
  • pind (array-like): index of the ancestor of each particle in part
  • future_trajs (array-like): particle estimate for {t+1:T}
  • find (array-like): index in future_trajs corresponding to each particle in part
  • ut (array-like): input signals for {0:T}
  • yt (array-like): measurements for {0:T}
  • tt (array-like): time stamps for {0:T}
  • cur_ind (int): index of current timestep (in ut, yt and tt)
Returns:
(array-like) with first dimension = N, logp(x_{t+1:T}|x_t^i)
logp_xnext_singlestep(part, past_trajs, pind, future_parts, find, ut, yt, tt, cur_ind)

Return the log-pdf value for the first step of the future trajectory. Needed in e.g MHIPS

Args:

  • part (array-like): Model specific representation of all particles, with first dimension = N (number of particles)
  • past_trajs: Trajectory leading up to current time
  • pind: Indices relating part to past_trajs
  • future_parts (array-like): particle estimate for {t+1}, stored using ‘filtered’ particle representation, ie. sample_smooth has not been performed on them
  • find: Indices relatin part and future_parts
  • ut (array-like): input signals for {1:T}
  • yt (array-like): measurements for {1:T}
  • tt (array-like): time stamps for {1:T}
  • cur_ind: index for current time
Returns:
(array-like) with first dimension = N, logp(x_{t+1:T}|x_t^i)

FFBSi Non-Markov

class pyparticleest.interfaces.FFBSiNonMarkov

FFBSi with Rejection Sampling

class pyparticleest.interfaces.FFBSiRS

Base class for models to be used with rejection sampling methods

logp_xnext_max(particles, u, t)

Return the max log-pdf value for all possible future states’ given input u

Args:

  • particles (array-like): Model specific representation of all particles, with first dimension = N (number of particles)
  • next_part (array-like): particle estimate for t+1
  • u (array-like): input signal
  • t (float): time stamps
Returns:
(array-like) with first dimension = N, argmax_{x_{t+1}} logp(x_{t+1}|x_t)
class pyparticleest.interfaces.FFBSiRSNonMarkov

Base class for models to be used with rejection sampling methods

logp_xnext_max_full(part, past_trajs, pind, uvec, yvec, tvec, cur_ind)

Return the max log-pdf value for all possible future states’ given input u

Args:

  • part (array-like): Model specific representation of all particles, with first dimension = N (number of particles)
  • past_trajs: Trajectory leading up to current time
  • pind: Indices relating part to past_trajs
  • uvec (array-like): input signals for {1:T}
  • yvec (array-like): measurements for {1:T}
  • tvec (array-like): time stamps for {1:T}
  • cur_ind: index for current time
Returns:
(array-like) with first dimension = N, argmax_{x_{t+1}} logp(x_{t+1}|x_t)

Proposal methods, e.g. MHIPS and MHBP

class pyparticleest.interfaces.SampleProposer

Base class for models to be used with methods that require drawing of new samples. Here ‘q’ is the name we give to the proposal distribtion.

logp_proposal(prop_part, ptraj, anc, future_trajs, find, yt, ut, tt, cur_ind)

Eval the log-propability of the proposal distribution

Args:
  • prop_part (array-like): Proposed particle estimate, first dimension has length = N
  • ptraj: array of trajectory step objects from previous time-steps, last index is step just before the current
  • anc (array-like): index of the ancestor of each particle in part
  • future_trajs (array-like): particle estimate for {t+1:T}
  • find (array-like): index in future_trajs corresponding to each generated sample
  • ut (array-like): input signals for {0:T}
  • yt (array-like): measurements for {0:T}
  • tt (array-like): time stamps for {0:T}
  • cur_ind (int): index of current timestep (in ut, yt and tt)
Returns
(array-like) with first dimension = N, log q(x_t | x_{t-1}, x_{t+1:T}, y_t:T)
propose_smooth(ptraj, anc, future_trajs, find, yt, ut, tt, cur_ind)

Sample from a distribution q(x_t | x_{t-1}, x_{t+1:T}, y_t:T)

Args:
  • ptraj: array of trajectory step objects from previous time-steps, last index is step just before the current
  • anc (array-like): index of the ancestor of each particle in part
  • future_trajs (array-like): particle estimate for {t+1:T}
  • find (array-like): index in future_trajs corresponding to each generated sample
  • ut (array-like): input signals for {0:T}
  • yt (array-like): measurements for {0:T}
  • tt (array-like): time stamps for {0:T}
  • cur_ind (int): index of current timestep (in ut, yt and tt)
Returns:
(array-like) of dimension N, wher N is the dimension of partp and/or future_trajs (one of which may be ‘None’ at the start/end of the dataset)

Parameter estimation

Numerical gradient

class pyparticleest.paramest.interfaces.ParamEstInterface

Interface s for particles to be used with the parameter estimation algorithm presented in [1] [1] - ‘System identification of nonlinear state-space models’ by Schon, Wills and Ninness

eval_logp_xnext(particles, particles_next, u, t)

Calculate gradient of a term of the I2 integral approximation as specified in [1].

Eg. Evaluate log p(x_{t+1}|x_t) or sum log p(x_{t+1}|x_t)

Args:
  • particles (array-like): Model specific representation of all particles, with first dimension = N (number of particles)
  • x_next (array-like): future states
  • t (float): time stamp

Returns: (array-like) or (float)

eval_logp_y(particles, y, t)

Calculate gradient of a term of the I3 integral approximation as specified in [1].

Eg. Evaluate log p(y_t|x_t) or sum log p(y_t|x_t)

Args:
  • particles (array-like): Model specific representation of all particles, with first dimension = N (number of particles)
  • y (array-like): measurement
  • t (float): time stamp

Returns: (array-like) or (float)

Analytic gradient

class pyparticleest.paramest.interfaces.ParamEstInterface_GradientSearch

Interface s for particles to be used with the parameter estimation algorithm presented in [1] using analytic gradients

eval_logp_x0_val_grad(particles, t)

Calculate term of the I1 integral approximation as specified in [1]. Eg. Evaluate log p(x_0) or sum log p(x_0)

Args:
  • particles (array-like): Model specific representation of all particles, with first dimension = N (number of particles)
  • t (float): time stamp

The gradient is an array where each element is the derivative with respect to the corresponding parameter

Returns: ((array-like) or (float), array-like) (value, gradient)

eval_logp_xnext_val_grad(particles, particles_next, u, t)

Calculate gradient of a term of the I2 integral approximation as specified in [1].

Eg. Evaluate log p(x_{t+1}|x_t) or sum log p(x_{t+1}|x_t)

Args:
  • particles (array-like): Model specific representation of all particles, with first dimension = N (number of particles)
  • particles_next (array-like): future states
  • u (array-like): input signal
  • t (float): time stamp

The gradient is an array where each element is the derivative with respect to the corresponding parameter

Returns: ((array-like) or (float), array-like) (value, gradient)

eval_logp_y_val_grad(particles, y, t)

Calculate gradient of a term of the I3 integral approximation as specified in [1].

Eg. Evaluate log p(y_t|x_t) or sum log p(y_t|x_t)

Args:
  • particles (array-like): Model specific representation of all particles, with first dimension = N (number of particles)
  • y (array-like): measurement
  • t (float): time stamp

The gradient is an array where each element is the derivative with respect to the corresponding parameter

Returns: ((array-like) or (float), array-like) (value, gradient)