simdkalman primitives module

Low-level Kalman filter computation steps with multi-dimensional input arrays. Unlike with the KalmanFilter class, all inputs must be numpy arrays. However, their dimensions can flexibly vary form 1 to 3 as long as they are reasonable from the point of view of matrix multiplication and numpy broadcasting rules. Matrix operations are applied on the last two axes of the arrays.

For example, to apply a prediction step to 3 different states \(x\) of dimension 2 with 3 different state transition matrices \(A\), one could do

from simdkalman.primitives import predict
import numpy

# different states
m1, m2, m3 = [numpy.array([[0],[i]]) for i in range(3)]
# with the same covariance (initially)
P = numpy.eye(2)

# different transition matrices
A1, A2, A3 = [numpy.eye(2)*i for i in range(3)]
# same process noise
Q = numpy.eye(2)*0.01

# stack correctly
m = numpy.vstack([
    m1[numpy.newaxis, ...],
    m2[numpy.newaxis, ...],
    m3[numpy.newaxis, ...]
])
A = numpy.vstack([
    A1[numpy.newaxis, ...],
    A2[numpy.newaxis, ...],
    A3[numpy.newaxis, ...]
])

# predict
m, P = predict(m, P, A, Q)

print(m.shape)
print(P.shape)
(3, 2, 1)
(3, 2, 2)
simdkalman.primitives.predict(mean, covariance, state_transition, process_noise)

Kalman filter prediction step

Parameters:
  • mean\({\mathbb E}[x_{j-1}]\), the filtered mean form the previous step

  • covariance\({\rm Cov}[x_{j-1}]\), the filtered covariance form the previous step

  • state_transition – matrix \(A\)

  • process_noise – matrix \(Q\)

Return type:

(prior_mean, prior_cov) predicted mean and covariance \({\mathbb E}[x_j]\), \({\rm Cov}[x_j]\)

simdkalman.primitives.update(prior_mean, prior_covariance, observation_model, observation_noise, measurement)

Kalman filter update step

Parameters:
  • prior_mean\({\mathbb E}[x_j|y_1,\ldots,y_{j-1}]\), the prior mean of \(x_j\)

  • prior_covariance\({\rm Cov}[x_j|y_1,\ldots,y_{j-1}]\), the prior covariance of \(x_j\)

  • observation_model – matrix \(H\)

  • observation_noise – matrix \(R\)

  • measurement – observation \(y_j\)

Return type:

(posterior_mean, posterior_covariance) posterior mean and covariance \({\mathbb E}[x_j|y_1,\ldots,y_j]\), \({\rm Cov}[x_j|y_1,\ldots,y_j]\) after observing \(y_j\)

simdkalman.primitives.update_with_nan_check(prior_mean, prior_covariance, observation_model, observation_noise, measurement)

Kalman filter update with a check for NaN observations. Like update but returns (prior_mean, prior_covariance) if measurement is NaN

simdkalman.primitives.smooth(posterior_mean, posterior_covariance, state_transition, process_noise, next_smooth_mean, next_smooth_covariance)

Kalman smoother backwards step

Parameters:
  • posterior_mean\({\mathbb E}[x_j|y_1,\ldots,y_j]\), the filtered mean of \(x_j\)

  • posterior_covariance\({\rm Cov}[x_j|y_1,\ldots,y_j]\), the filtered covariance of \(x_j\)

  • state_transition – matrix \(A\)

  • process_noise – matrix \(Q\)

  • next_smooth_mean\({\mathbb E}[x_{j+1}|y_1,\ldots,y_T]\)

  • next_smooth_covariance\({\rm Cov}[x_{j+1}|y_1,\ldots,y_T]\)

Return type:

(smooth_mean, smooth_covariance, smoothing_gain) smoothed mean \({\mathbb E}[x_j|y_1,\ldots,y_T]\), and covariance \({\rm Cov}[x_j|y_1,\ldots,y_T]\)

simdkalman.primitives.predict_observation(mean, covariance, observation_model, observation_noise)

Compute probability distribution of the observation \(y\), given the distribution of \(x\).

Parameters:
  • mean\({\mathbb E}[x]\)

  • covariance\({\rm Cov}[x]\)

  • observation_model – matrix \(H\)

  • observation_noise – matrix \(R\)

Return type:

mean \({\mathbb E}[y]\) and covariance \({\rm Cov}[y]\)