_kalman_smooth

This module implements Kalman filters

pynumdiff.kalman_smooth._kalman_smooth.constant_acceleration(x, dt, params, options=None)

Run a forward-backward constant acceleration RTS Kalman smoother to estimate the derivative.

Parameters
  • x (np.array (float)) – array of time series to differentiate

  • dt (float) – time step size

  • params (dict {'forwardbackward': boolean}, optional) –

    a list of two elements:

    • r: covariance of the x noise

    • q: covariance of the constant velocity model

  • options – a dictionary indicating whether to run smoother forwards and backwards (usually achieves better estimate at end points)

Returns

a tuple consisting of:

  • x_hat: estimated (smoothed) x

  • dxdt_hat: estimated derivative of x

Return type

tuple -> (np.array, np.array)

pynumdiff.kalman_smooth._kalman_smooth.constant_jerk(x, dt, params, options=None)

Run a forward-backward constant jerk RTS Kalman smoother to estimate the derivative.

Parameters
  • x (np.array (float)) – array of time series to differentiate

  • dt (float) – time step size

  • params (dict {'forwardbackward': boolean}, optional) –

    a list of two elements:

    • r: covariance of the x noise

    • q: covariance of the constant velocity model

  • options – a dictionary indicating whether to run smoother forwards and backwards (usually achieves better estimate at end points)

Returns

a tuple consisting of:

  • x_hat: estimated (smoothed) x

  • dxdt_hat: estimated derivative of x

Return type

tuple -> (np.array, np.array)

pynumdiff.kalman_smooth._kalman_smooth.constant_velocity(x, dt, params, options=None)

Run a forward-backward constant velocity RTS Kalman smoother to estimate the derivative.

Parameters
  • x (np.array (float)) – array of time series to differentiate

  • dt (float) – time step size

  • params (dict {'forwardbackward': boolean}, optional) –

    a list of two elements:

    • r: covariance of the x noise

    • q: covariance of the constant velocity model

  • options – a dictionary indicating whether to run smoother forwards and backwards (usually achieves better estimate at end points)

Returns

a tuple consisting of:

  • x_hat: estimated (smoothed) x

  • dxdt_hat: estimated derivative of x

Return type

tuple -> (np.array, np.array)

pynumdiff.kalman_smooth._kalman_smooth.known_dynamics(x, params, u=None, options=None)

Run a forward RTS Kalman smoother given known dynamics to estimate the derivative.

Parameters
  • x (np.array (float)) – matrix of time series of (noisy) measurements

  • params (dict {'smooth': boolean}, optional) – a list of: - x0: inital condition, matrix of Nx1, N = number of states - P0: initial covariance matrix of NxN - A: dynamics matrix, NxN - B: control input matrix, NxM, M = number of measurements - C: measurement dynamics, MxN - R: covariance matrix for the measurements, MxM - Q: covariance matrix for the model, NxN

  • u (np.array (float)) – matrix of time series of control inputs

  • options – a dictionary indicating whether to run smoother

Returns

matrix: - xhat_smooth: smoothed estimates of the full state x

Return type

tuple -> (np.array, np.array)

pynumdiff.kalman_smooth._kalman_smooth.savgol_const_accel(x, dt, params, options=None)

Run a forward-backward constant acceleration RTS Kalman smoother to estimate the derivative, where initial estimates of the velocity are first estimated using the savitzky-golay filter.

Parameters
  • x (np.array (float)) – array of time series to differentiate

  • dt (float) – time step size

  • params (dict {'forwardbackward': boolean}, optional) – a list of six elements: - N: for savgoldiff, order of the polynomial - window_size: for savgoldiff, size of the sliding window, must be odd (if not, 1 is added) - smoothing_win: for savgoldiff, size of the window used for gaussian smoothing, a good default is window_size, but smaller for high frequnecy data - r1: covariance of the x noise - r2: covariance of the vel noise - q: covariance of the constant velocity model

  • options – a dictionary indicating whether to run smoother forwards and backwards (usually achieves better estimate at end points)

Returns

a tuple consisting of:

  • x_hat: estimated (smoothed) x

  • dxdt_hat: estimated derivative of x

Return type

tuple -> (np.array, np.array)