Base classes¶
Mixed Linear/Nonlinear Gaussian¶
Model definition for base class for Mixed Linear/Nonlinear Gaussian systems
@author: Jerker Nordh
-
class
pyparticleest.models.mlnlg.
MixedNLGaussianMarginalized
(lxi, lz, Az=None, C=None, Qz=None, R=None, fz=None, Axi=None, Qxi=None, Qxiz=None, fxi=None, h=None, params=None, **kwargs)¶ This class implements a fully marginalized smoother for mixed linear/nonlinear models, in contrast to the MixedNLGaussian class it never samples the linear states.
This is somewhat slower, and doesn’t readily admit using rejection sampling, it is up to the end user which method is best for their particular problem
-
calc_prop1
(particles, next_part, u, t)¶ internal helper function
-
calc_prop3
(particles, Omega, Lambda, u, t)¶ internal helper function
-
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)
-
sample_smooth
(part, ptraj, anc, future_trajs, find, ut, yt, tt, cur_ind)¶ Calculate statistics needed when evaluating the logp_xnext_full for the marginalized trajectory
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
-
-
class
pyparticleest.models.mlnlg.
MixedNLGaussianSampled
(lxi, lz, Az=None, C=None, Qz=None, R=None, fz=None, Axi=None, Qxi=None, Qxiz=None, fxi=None, h=None, params=None, **kwargs)¶ Base class for particles of the type mixed linear/non-linear with additive gaussian noise.
Implement this type of system by extending this class and provide the methods for returning the system matrices at each time instant.
xi_{t+1} = f_xi + A_xi*z_t + v_xi, z_{t+1} = f_z + A_z*z_t + v_z, y_t = h + C*z_t + e, (v_xi, v_z)^T ~ N(0, (Q_xi, Qxiz Qxiz^T, Qz)) e ~ N (0, R)- Args:
- lxi (int): number of nonlinear states
- lz (int): number of linear states
- Az (arraylike): Az (if constant)
- C (arraylike): C (if constant)
- Qz (arraylike): Qz (if constant)
- R (arraylike): R (if constant)
- fz (arraylike): fz (if constant)
- Axi (arraylike): Axi (if constant)
- Qxi (arraylike): Qxi (if constant)
- Qxiz (arraylike): Qxiz (if constant)
- fxi (arraylike): fxi (if constant)
- h (arraylike): h (if constant)
- params (array-like): model parameters (if any)
-
calc_A_f_Q
(particles, u, t)¶ Calculate the A, f and Q matrices for the particles. Where A, f and Q are the stacked matrices of (A_xi, A_z) and so on
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:
- (A, f, Q, A_identical f_identical, Q_identical). The last elements indicate if the matrices are identical for all particles
-
calc_cond_dynamics
(particles, xi_next, u, t)¶ Calculates the linear dynamics for each particle
- Args:
- particles (array-like): Model specific representation of all particles, with first dimension = N (number of particles)
- xi_next (array-like): next non linear state
- u (array-like): input signal
- t (float): time stamp
- Returns:
- (Az, fz, Qz):
- Az (array-like): Az matrix for each particle
- fz (array-like): fz vector for each particle
- Qz (array-lie): Noise covariance for each particle
-
calc_l2
(xin, zn, Pn, zl, Pl, A, f, M)¶ Internal helper function
-
calc_l2_grad
(xin, zn, Pn, zl, Pl, A, f, M, f_grad, A_grad)¶ Internal helper function
-
calc_l3
(y, zl, Pl, Cl, hl)¶ internal helper function
-
calc_l3_grad
(y, zl, Pl, Cl, hl, C_grad, h_grad)¶ internal helper function
-
calc_xi_next
(particles, noise, u, t)¶ Calculate the next nonlinear state given the input and noise realization
- 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 for each particle
- Returns:
- (array-like): xi values for future particles
-
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)
-
eval_logp_x0
(particles, t)¶ Evaluate 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
-
eval_logp_x0_val_grad
(particles, t)¶ Evaluate gradient of 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
-
eval_logp_xi0
(xil)¶ Evaluate logprob of the initial non-linear state eta, default implementation assumes all are equal, override this if another behavior is desired
- Args:
- xil (list): Initial xi states
-
eval_logp_xi0_grad
(xil)¶ Evaluate logprob of the initial non-linear state eta, default implementation assumes all are equal, override this if another behavior is desired
- Args:
- xil (list): Initial xi states
-
eval_logp_xnext
(particles, x_next, u, t)¶ Evaluate 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: (float)
-
eval_logp_xnext_val_grad
(particles, x_next, u, t)¶ Evaluate value and gradient 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: ((float), (array-like))
-
eval_logp_y
(particles, y, t)¶ Evaluate value of 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: (float)
-
eval_logp_y_val_grad
(particles, y, t)¶ Evaluate value and gradient of 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: ((float), (array-like))
-
get_cross_covariance
(particles, u, t)¶ Return cross-covariance between noise for nonlinear and linear states
- 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
- Returns:
- (array-like): Qxiz(xi_t, u_t, t) for each particle
-
get_meas_dynamics_grad
(particles, y, t)¶ Override this method if (C, h, R) depends on the parameters
- Args:
- particles (array-like): Model specific representation of all particles, with first dimension = N (number of particles)
- y (array-like): measurment
- t (float): time stamps
- Returns:
- (C_grad, h_grad, R_grad): Element-wise gradients with respect to all the parameters for the system matrices
-
get_pred_dynamics_grad
(particles, u, t)¶ Override this method if (A, f, Q) depends on the parameters
- 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 stamps
- Returns:
- (A_grad, f_grad, Q_grad): Element-wise gradients with respect to all the parameters for the system matrices
-
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)
-
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_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)
-
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)
-
meas_xi_next
(particles, xi_next, u, t)¶ Update estimate using observation of next state
- Args:
- particles (array-like): Model specific representation of all particles, with first dimension = N (number of particles)
- xi_next (array-like): future nonlinear states
- u (array-like): input signal
- t (float): time stamp
-
measure
(particles, y, t)¶ Return the log-pdf value of the measurement and update the statistics for the linear states
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)
-
pred_xi
(particles, u, t)¶ Predict the next nonlinear state given the input
- 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
- Returns:
- (array-like): xi values for future particles
-
propose_smooth
(ptraj, anc, future_trajs, find, yt, ut, tt, cur_ind)¶ Sample from a distribution q(x_t | x_{0:t-1}, x_{t+1:T}, y_0: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)
-
sample_process_noise
(particles, u, t)¶ Return sampled process noise for the non-linear states
- 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
- 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
-
set_dynamics
(Az=None, fz=None, Qz=None, R=None, Axi=None, fxi=None, Qxi=None, Qxiz=None, C=None, h=None)¶ Update dynamics, typically used when changing the system dynamics due to a parameter change
- Args:
- lxi (int): number of nonlinear states
- lz (int): number of linear states
- Az (arraylike): Az (if constant)
- C (arraylike): C (if constant)
- Qz (arraylike): Qz (if constant)
- R (arraylike): R (if constant)
- fz (arraylike): fz (if constant)
- Axi (arraylike): Axi (if constant)
- Qxi (arraylike): Qxi (if constant)
- Qxiz (arraylike): Qxiz (if constant)
- fxi (arraylike): fxi (if constant)
- h (arraylike): h (if constant)
-
set_params
(params)¶ This methods should be overriden if the system dynamics depends on any parameters, this method should however be called to store the new parameter values correctly
- Args:
- params (array-like): new parameter values
-
pyparticleest.models.mlnlg.
factor_psd
(A)¶ internal helper function
Nonlinear Gaussian¶
Model definition for base class for Nonlinear Gaussian systems
@author: Jerker Nordh
-
class
pyparticleest.models.nlg.
NonlinearGaussian
(lxi, f=None, g=None, Q=None, R=None)¶ Base class for particles of the type mixed linear/non-linear with additive gaussian noise.
Implement this type of system by extending this class and provide the methods for returning the system matrices at each time instant.
x_{t+1} = f(x_t, u_t) + v, v ~ N(0, Q(x_t, u_t)) y_t = g(x_t) + e, e ~ N(=, R(x_t))
This class currently doesn’t support analytic gradients when performing parameter estimation, however using numerical gradients is typically fine
- Args:
- lxi (int): number of states in model
- f (array-like): f (if constaint)
- g (array-like): g (if constaint)
- Q (array-like): Q (if constaint)
- R (array-like): R (if constaint)
-
calc_Q
(particles, u, t)¶ Calucate Q
- 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
- Returns:
- (array-like): Q for all particles
-
calc_R
(particles, t)¶ Calucate R
- Args:
- particles (array-like): Model specific representation of all particles, with first dimension = N (number of particles)
- t (float): time stamp
- Returns:
- (array-like): R for all particles
-
calc_f
(particles, u, t)¶ Calucate f
- 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
- Returns:
- (array-like): f for all particles
-
calc_g
(particles, t)¶ Calucate g
- Args:
- particles (array-like): Model specific representation of all particles, with first dimension = N (number of particles)
- t (float): time stamp
- Returns:
- (array-like): g for all particles
-
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)
-
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)
-
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_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)
-
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)
-
propose_smooth
(ptraj, anc, future_trajs, find, yt, ut, tt, cur_ind)¶ Sample from a distribution q(x_t | x_{0: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)
-
sample_process_noise
(particles, u, t)¶ Sample process noise
- 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
- Returns:
- (array-like) with first dimension = N
-
set_params
(params)¶ This methods should be overriden if the system dynamics depends on any parameters, this method should however be called to store the new parameter values correctly
- Args:
- params (array-like): new parameter values
-
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.models.nlg.
NonlinearGaussianInitialGaussian
(x0=None, Px0=None, lxi=None, **kwargs)¶ Nonlinear gaussian system with initial Gaussian distribution.
- Args:
- x0 (array-like): mean value of initial state, defaults to 0
- Px0 (array-like): covariance of initial state, defaults to 0
- lxi (int): number of states, only needed if neither x0 or Px0 specified
-
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
-
eval_logp_x0
(particles, t)¶ Evaluate 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
Linear Time-Varying¶
Model definition for base class for Linear Time-varying systems @author: Jerker Nordh
-
class
pyparticleest.models.ltv.
LTV
(z0, P0, A=None, C=None, Q=None, R=None, f=None, h=None, params=None, **kwargs)¶ Base class for particles of the type linear time varying with additive gaussian noise.
Implement this type of system by extending this class and provide the methods for returning the system matrices at each time instant
z_{t+1} = A*z_t + f + v, v ~ N(0, Q) y_t = C*z_t + h + e, e ~ N(0,R)
- Args:
- z0: Initial mean value of the state estimate
- P0: Coviariance of initial z estimate
- A (array-like): A matrix (if constant)
- C (array-like): C matrix (if constant)
- Q (array-like): Q matrix (if constant)
- R (array-like): R matrix (if constant)
- f (array-like): f vector (if constant)
- h (array-like): h vector (if constant)
- params (array-like): model parameters (if any)
-
calc_l1
(z, P, z0, P0)¶ internal helper function
-
calc_l1_grad
(z, P, z0, P0, z0_grad)¶ internal helper function
-
calc_l2
(zn, Pn, z, P, A, f, M)¶ internal helper function
-
calc_l2_grad
(zn, Pn, z, P, A, f, M, A_grad, f_grad)¶ internal helper function
-
calc_l3
(y, z, P)¶ internal helper function
-
calc_l3_grad
(y, z, P, C_grad, h_grad)¶ internal helper function
-
create_initial_estimate
(N)¶ Sample particles from initial distribution
- Args:
- N (int): Number of particles to sample, since the estimate is deterministic there is no reason for N > 1
- Returns:
- (array-like) with first dimension = N, model specific representation of all particles
-
eval_logp_x0
(particles, t)¶ Evaluate 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
-
eval_logp_x0_val_grad
(particles, t)¶ Evaluate gradient of 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
-
eval_logp_xnext
(particles, x_next, u, t)¶ Evaluate 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)
-
eval_logp_xnext_val_grad
(particles, x_next, u, t)¶ Evaluate value and gradient of 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), (array-like))
-
eval_logp_y
(particles, y, t)¶ Evaluate value of 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)
-
eval_logp_y_val_grad
(particles, y, t)¶ Evaluate value and gradient of 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), (array-like))
-
fwd_peak_density
(u, t)¶ No need for rejections sampling for this type of model, always returns 0.0 since all particles are equivalent
- Args:
- u: Unused
- t: Unused
- Returns
- (float) 0.0
-
get_initial_grad
()¶ Default implementation has no dependence on xi, override if needed
Calculate gradient estimate of initial state for linear state condition on the nonlinear estimate
- Args:
- xi0 (array-like): Initial xi states
- Returns:
- (z,P): z is a list of element-wise gradients for the inital mean values, P is a list of element-wise gradients for the covariance matrices
-
get_meas_dynamics
(y, t)¶ Return matrices describing affine relation of measurement and current state estimates
y_t = C*z_t + h + e, e ~ N(0,R)
- 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:
- (y, C, h, R): y is a preprocessed measurement, the rest are lists with the corresponding matrix for each particle. None indicates that the matrix is identical for all particles and the value stored in this class should be used instead
-
get_meas_dynamics_grad
(y, t)¶ Override this method if (C, h, R) depends on the parameters
- Args:
- particles (array-like): Model specific representation of all particles, with first dimension = N (number of particles)
- y (array-like): measurment
- t (float): time stamps
- Returns:
- (C_grad, h_grad, R_grad): Element-wise gradients with respect to all the parameters for the system matrices
-
get_pred_dynamics
(u, t)¶ Return matrices describing affine relation of next nonlinear state conditioned on the current time and input signal
z_{t+1} = A*z_t + f + v, v ~ N(0, Q)
- 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
- Returns:
- (A, f, Q) where each element is a list with the corresponding matrix for each particle. None indicates that the matrix is identical for all particles and the value stored in this class should be used instead
-
get_pred_dynamics_grad
(u, t)¶ Override this method if (A, f, Q) depends on the parameters
- 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 stamps
- Returns:
- (A_grad, f_grad, Q_grad): Element-wise gradients with respect to all the parameters for the system matrices
-
get_states
(particles)¶ Return the estimates contained in the particles array
- Args:
- particles (array-like): Model specific representation of all particles, with first dimension = N (number of particles)
- Returns
- (zl, Pl):
- zl: list of mean values for z
- Pl: list of covariance matrices for z
-
logp_xnext
(particles, next_part, u, t)¶ Return the log-pdf value for the possible future state ‘next’ given input u.
Always returns zeros since all particles are always equivalent for this type of model
Args:
- particles (array-like): Model specific representation of all particles, with first dimension = N (number of particles)
- next_part: Unused
- u: Unused
- t: Unused
- Returns:
- (array-like) with first dimension = N, numpu.zeros((N,))
-
measure
(particles, y, t)¶ Return the log-pdf value of the measurement and update the statistics for the states
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)
-
sample_process_noise
(particles, u, t)¶ There is no need to sample noise for this type of model
Args:
- particles: Unused
- next_part: Unused
- u: Unused
- t: Unused
- Returns:
- None
-
sample_smooth
(part, ptraj, anc, future_trajs, find, ut, yt, tt, cur_ind)¶ Update sufficient statistics based on the future states
- 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
-
set_states
(particles, z_list, P_list)¶ Set the estimate of the states
- Args:
- particles (array-like): Model specific representation of all particles, with first dimension = N (number of particles)
- z_list (list): list of mean values for z for each particle
- P_list (list): list of covariance matrices for z for each particle
-
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: Unused for this type of model
- Returns:
- (array-like) with first dimension = N, particle estimate at time t+1
Hierarchial¶
Model definition for base class for hierarchical systems
@author: Jerker Nordh
-
class
pyparticleest.models.hierarchial.
HierarchicalBase
(len_xi, len_z, **kwargs)¶ Base class for Rao-Blackwellization of hierarchical models
- Args:
- len_xi (int): number of nonlinear states
- len_z (int): number of linear states
-
calc_cond_dynamics
(particles, xi_next, u, t)¶ Calculates the linear dynamics for each particle
- Args:
- particles (array-like): Model specific representation of all particles, with first dimension = N (number of particles)
- xi_next (array-like): next non linear state
- u (array-like): input signal
- t (float): time stamp
- Returns:
- (Az, fz, Qz):
- Az (array-like): Az matrix for each particle
- fz (array-like): fz vector for each particle
- Qz (array-lie): Noise covariance for each particle
-
calc_xi_next
(particles, u, t, noise)¶ Calculate the next nonlinear state given the input and noise realization
- 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 for each particle
- Returns:
- (array-like): xi values for future particles
-
logp_xnext
(particles, next_part, u, t)¶ Return the log-pdf value for the possible future state ‘next_part’ given input u
If Nn = 1 all particle are evaluated against the same future state, otherwise N must equal Nn and each particle is only evaluated against the future state with the same index.
- Args:
- particles (array-like): Model specific representation of all particles, with first dimension = N (number of particles)
- next_part (array-like): future states, with first dimension = Nn
- u (array-like): input signal
- t (float): time stamp
- Returns:
- (array-like):
- log-probability of the future state for each particle
-
logp_xnext_xi
(particles, next_xi, u, t)¶ Evaluate the log-probability of the next nonlinear state
- Args:
- particles (array-like): Model specific representation of all particles, with first dimension = N (number of particles)
- next_xi (array-like): future nonlinear state
- u (array-like): input signal
- t (float): time stamp
- Returns:
- (array-like):
- log-probability of the future nonlinear state for each particle
-
measure
(particles, y, t)¶ Return the log-pdf value of the measurement and update the statistics for the linear states
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)
-
measure_nonlin
(particles, y, t)¶ Measurement probability for the nonlinear parts of the measurement equations
- 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):
- log-probability of the measurement for each particle
-
sample_smooth
(part, ptraj, anc, future_trajs, find, ut, yt, tt, cur_ind)¶ Sampled linear state conditioned on future_trajs
- 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