Skip to content

Filters

Filtering routines for ODE models.

ode_filters.filters.PController dataclass

Proportional (single-error) step-size controller.

The proposed step is

h_new = h * safety * err^(-alpha)

clipped to [min_factor, max_factor] * h. The previous error is ignored, so this controller has no memory.

Default alpha = 1 / order matches the standard "I-controller" of Hairer-Wanner-Norsett (1993) and the elementary controller in scipy's solve_ivp. For an EKF1 with IWP(q) prior, order = q.

Parameters:

Name Type Description Default
order int

Convergence order of the local error estimate.

required
safety float

Safety factor on the proposed step.

0.9
alpha float | None

Gain on the current error; defaults to 1.0 / order.

None
min_factor float

Lower clip on h_new / h.

0.2
max_factor float

Upper clip on h_new / h.

5.0

ode_filters.filters.PController.propose(h: float, err: float, err_prev: float | None = None) -> float

Return the next proposed step size.

Parameters:

Name Type Description Default
h float

Current step size.

required
err float

Normalised error of the current step (1.0 = at tolerance).

required
err_prev float | None

Ignored; accepted for protocol compatibility with :class:PIController.

None

ode_filters.filters.PIController dataclass

Gustafsson-style proportional-integral controller.

The proposed step is

h_new = h * safety * err^(-alpha) * (err_prev / err)^(beta)

clipped to [min_factor, max_factor] * h. When err_prev is unknown (e.g. on the first step or right after a reject) the I-term is dropped and the update reduces to the :class:PController form h_new = h * safety * err^(-alpha).

Defaults follow Gustafsson (1991): alpha = 0.7 / order, beta = 0.4 / order, which are also the values used in Bosch et al. (2021). For an EKF1 with IWP(q) prior, order = q.

Parameters:

Name Type Description Default
order int

Convergence order of the local error estimate.

required
safety float

Safety factor on the proposed step.

0.9
alpha float | None

Proportional gain; defaults to 0.7 / order.

None
beta float | None

Integral gain; defaults to 0.4 / order. Set beta=0 for a pure proportional controller -- :class:PController is the standalone equivalent with the more natural default alpha.

None
min_factor float

Lower clip on h_new / h.

0.2
max_factor float

Upper clip on h_new / h.

5.0

ode_filters.filters.PIController.propose(h: float, err: float, err_prev: float | None) -> float

Return the next proposed step size.

Parameters:

Name Type Description Default
h float

Current step size.

required
err float

Normalised error of the current step (1.0 = at tolerance).

required
err_prev float | None

Normalised error of the previously accepted step, or None if unavailable.

required

ode_filters.filters.StepSizeController

Bases: Protocol

Duck-typed interface for step-size controllers.

err_prev may be None to indicate the absence of memory (first step, or right after a reject).

ode_filters.filters.StepSizeController.propose(h: float, err: float, err_prev: float | None) -> float

Return the proposed next step size.

ode_filters.filters.AdaptiveLoopResult

Bases: NamedTuple

Output of :func:ekf1_sqr_adaptive_loop.

All sequences are aligned to accepted steps. t_seq has length N_accepted + 1 (initial time plus one entry per accepted step); the per-step sequences (m_pred_seq, Pz_seq_sqr, ...) and the smoothing quantities (G_back_seq, d_back_seq, P_back_seq_sqr) have length N_accepted. Filtered means/covariances (m_seq, P_seq_sqr) have length N_accepted + 1 and include the initial state.

ode_filters.filters.ekf1_sqr_adaptive_loop(mu_0: Array, Sigma_0_sqr: Array, prior: BasePrior, measure: BaseODEInformation, tspan: tuple[float, float], *, atol: float = 0.0001, rtol: float = 0.01, h_init: float | None = None, h_min: float = 1e-10, h_max: float | None = None, controller: StepSizeController | None = None, calibration: CalibrationMode = 'dynamic', sigma_in_error: SigmaInError = 'running_mean', min_sigma_sqr: float = 0.0, calibrate: bool | None = None, max_steps: int = 100000) -> AdaptiveLoopResult

Adaptive-step square-root EKF with per-step diffusion calibration.

Python while driver around a jitted per-step body. Every accepted step contributes to the returned sequences and the result carries the smoothing-relevant outputs, so it can be fed directly to :func:rts_sqr_smoother_loop.

Parameters:

Name Type Description Default
mu_0 Array

Initial state mean.

required
Sigma_0_sqr Array

Initial state covariance (square-root form).

required
prior BasePrior

Prior model (e.g. :class:IWP). The Kronecker structure Q(h) = kron(_Q(h), xi) is used to build the per-step process noise. prior.q sets the default controller order.

required
measure BaseODEInformation

Measurement model (e.g. :class:ODEInformation).

required
tspan tuple[float, float]

Time interval (t_start, t_end). Must be a tuple (hashable for jit-compatible inner step).

required
atol float

Absolute tolerance on the function-value local error.

0.0001
rtol float

Relative tolerance.

0.01
h_init float | None

Initial step. Defaults to (t_end - t_start) / 100.

None
h_min float

Minimum step. Steps below this raise RuntimeError.

1e-10
h_max float | None

Maximum step. Defaults to t_end - t_start.

None
controller StepSizeController | None

Step controller implementing :class:~ode_filters.filters.adaptive_controller.StepSizeController. Defaults to PIController(order=prior.q); pass :class:~ode_filters.filters.adaptive_controller.PController for a memoryless proportional-only controller.

None
calibration CalibrationMode

How the per-step sigma_hat^2 enters the stored posterior. See the module docstring for the four available modes:

  • "dynamic" (default) -- Bosch et al. 2021 Eq. 32. Scalar sigma_hat^2 from H Q(h) H.T baked into Q_h. Honest per-step but scalar: fails on multi-scale ODE systems (one component's residual dominates sigma_hat^2 and starves the others).
  • "diagonal_ekf0" -- per-component sigma_hat^2_i = m_z[i]^2 / (E_1 Q E_1.T)_ii using the EK0 observation matrix H_0 = E_1. The denominator is exactly diagonal when xi is diagonal, so the per-component estimator is the genuine MLE in a block-diagonal sub-model. The recommended choice for multi-component problems. Requires prior.xi diagonal.
  • "diagonal" -- same formula, but with the EK1 Jacobian H_t = E_1 - J_f E_0 in the denominator. Heuristic: taking the diagonal of a generally-dense H_t Q(h) H_t.T ignores cross-component coupling from J_f. Equivalent to "diagonal_ekf0" when J_f is block-diagonal (e.g. fully decoupled RHS); biased in proportion to the off-diagonal entries of J_f otherwise. Requires prior.xi diagonal.
  • "cumulative" -- legacy multiplicative scheme; scalar sigma from the full H P_pred H.T post-multiplies the whole step. Non-Markovian; robust on multi-scale by accident (inflated P inherited from earlier transitions). Prefer "diagonal_ekf0".
  • "none" -- propagate sigma=1, still report sigma_hat^2 for post-hoc rescaling / diagnostics.

All four modes work with :class:JointPrior / :class:PrecondJointPrior: the diagonal modes scale only the state block per-component and require prior._prior_x.xi to be diagonal.

For the diagonal modes the per-step entries of result.sigma_sqr_seq are length-d arrays; for the scalar modes they are Python floats.

'dynamic'
sigma_in_error SigmaInError

How sigma_hat^2 enters the local-error estimate.

  • "running_mean" (default) -- use the cumulative running mean of past accepted-step sigma_hat^2 values in the error formula (sigma_hat^2 is still per-step in the Q calibration). The running mean has variance 2/(nd), so estimator noise vanishes after a few accepted steps -- h becomes much smoother and reject counts typically drop by several times. For "diagonal" modes the running mean is per-component.
  • "per_step" -- use the just-computed sigma_hat^2 (Bosch et al. 2021 original recipe). Honest but inherits the chi-squared estimator noise (relative std sqrt(2/d)), which propagates into step-size decisions and causes visible h oscillations in flat regions.
'running_mean'
min_sigma_sqr float

Lower bound applied to the per-step sigma_hat^2 (and to each component of the per-component vector in diagonal modes) before it is baked into Q_step_sqr. Default 0.0 preserves the unclamped behavior. Use a small positive value (e.g. 1e-30) on problems where the trivial zero-residual fixed point would otherwise collapse the state-block diffusion to zero and propagate NaN.

0.0
calibrate bool | None

Deprecated. True maps to calibration="dynamic", False maps to calibration="none". Pass calibration directly instead.

None
max_steps int

Hard cap on iterations (rejected + accepted) as a safety valve against infinite loops.

100000

Returns:

Type Description
AdaptiveLoopResult

class:AdaptiveLoopResult with the accepted trajectory.

Raises:

Type Description
RuntimeError

If the controller proposes a step below h_min or the iteration cap is exceeded.

ValueError

If tspan is non-increasing or calibration / sigma_in_error is unknown.

ode_filters.filters.ekf1_sqr_loop(mu_0: Array, Sigma_0_sqr: Array, prior: BasePrior, measure: BaseODEInformation, tspan: tuple[float, float], N: int) -> LoopResult

Run a square-root EKF over N observation steps.

Parameters:

Name Type Description Default
mu_0 Array

Initial state mean estimate.

required
Sigma_0_sqr Array

Initial state covariance (square-root form).

required
prior BasePrior

Prior model (e.g., IWP or PrecondIWP).

required
measure BaseODEInformation

Measurement model (e.g., ODEInformation or subclass).

required
tspan tuple[float, float]

Time interval (t_start, t_end).

required
N int

Number of filter steps.

required

Returns:

Type Description
LoopResult

Tuple of 9 arrays and a scalar log marginal likelihood.

ode_filters.filters.ekf1_sqr_loop_dynamic(mu_0: Array, Sigma_0_sqr: Array, prior: BasePrior, measure: BaseODEInformation, tspan: tuple[float, float], N: int, *, calibration: str = 'dynamic', min_sigma_sqr: float = 0.0, calibrate: bool | None = None) -> DynamicLoopResult

Fixed-step square-root EKF with online per-step diffusion calibration.

Same fixed grid as :func:ekf1_sqr_loop, but at every step the per-step quasi-MLE diffusion is estimated and -- depending on calibration -- baked into the current step's Q_h before propagation. See :func:~ode_filters.filters.ekf1_sqr_adaptive_loop for the full description of the modes; this loop supports the same names:

  • "dynamic" (default, scalar): Bosch et al. 2021 Eq. 32. P_pred = A P_prev A.T + sigma_hat^2_n * Q_h.
  • "diagonal_ekf0": per-component sigma_hat^2_i from (E_1 Q E_1.T)_ii (exact-diagonal denominator). Genuine MLE in a block-diagonal sub-model for component i. Recommended default for multi-component problems.
  • "diagonal": same per-component formula but with the EK1 Jacobian H_t = E_1 - J_f E_0 in the denominator. Heuristic -- equivalent to "diagonal_ekf0" when J_f is block-diagonal in components, biased in proportion to J_f's off-diagonal entries otherwise.
  • "none": propagate sigma=1, report sigma_hat^2 for post-hoc rescaling (combine with :func:posthoc_mle_sigma_sqr / :func:rescale_sqr_seq).

"cumulative" is intentionally not exposed on the fixed-step loop: the cumulative scheme's appeal is robust adaptive-step control, and at fixed step there's nothing to gain over post-hoc rescaling.

The returned log_likelihood reflects the calibrated posterior -- each step's contribution uses the calibrated Pz_sqr. With calibration="none" it matches the log-likelihood from :func:ekf1_sqr_loop; in the calibrated modes it is the log-likelihood under the corresponding generative model and is not directly comparable across modes.

Parameters:

Name Type Description Default
mu_0 Array

Initial state mean.

required
Sigma_0_sqr Array

Initial state covariance (square-root form).

required
prior BasePrior

Prior model (e.g. IWP).

required
measure BaseODEInformation

Measurement model (e.g. ODEInformation).

required
tspan tuple[float, float]

Time interval (t_start, t_end).

required
N int

Number of filter steps.

required
calibration str

"dynamic" (default), "diagonal", or "none". All modes work with :class:JointPrior / :class:PrecondJointPrior; diagonal modes scale the state block only and require prior._prior_x.xi to be diagonal.

'dynamic'
min_sigma_sqr float

Lower bound on the per-step sigma_hat^2 (and on each component in diagonal modes) before it is baked into Q_step_sqr. Default 0.0 preserves unclamped behavior.

0.0
calibrate bool | None

Deprecated; True -> "dynamic", False -> "none".

None

Returns:

Type Description
DynamicLoopResult

Tuple of 10 sequences plus the scalar log-marginal-likelihood:

DynamicLoopResult

``(m_seq, P_seq_sqr, m_pred_seq, P_pred_seq_sqr, G_back_seq,

DynamicLoopResult

d_back_seq, P_back_seq_sqr, mz_seq, Pz_seq_sqr, sigma_sqr_seq,

DynamicLoopResult

log_likelihood). In the diagonal modessigma_sqr_seq``

DynamicLoopResult

entries are length-d_ode arrays; otherwise Python floats.

ode_filters.filters.ekf1_sqr_loop_dynamic_scan(mu_0: Array, Sigma_0_sqr: Array, prior: BasePrior, measure: BaseODEInformation, tspan: tuple[float, float], N: int, *, calibration: str = 'dynamic', min_sigma_sqr: float = 0.0, calibrate: bool | None = None) -> DynamicScanLoopResult

jax.lax.scan variant of :func:ekf1_sqr_loop_dynamic.

Same per-step semantics, executed inside a single jax.lax.scan for full JIT/compile-once efficiency. Returns stacked Array sequences rather than lists. Supports calibration in {"dynamic", "diagonal", "diagonal_ekf0", "none"}; in the diagonal modes sigma_sqr_seq has shape (N, d) rather than (N,).

Parameters:

Name Type Description Default
mu_0 Array

Initial state mean.

required
Sigma_0_sqr Array

Initial state covariance (square-root form).

required
prior BasePrior

Prior (e.g. IWP).

required
measure BaseODEInformation

Measurement model.

required
tspan tuple[float, float]

Time interval (t_start, t_end).

required
N int

Number of filter steps.

required
calibration str

"dynamic" (default), "diagonal", "diagonal_ekf0", or "none". All modes work with :class:JointPrior / :class:PrecondJointPrior; diagonal modes require prior._prior_x.xi to be diagonal.

'dynamic'
min_sigma_sqr float

Lower bound applied to per-step sigma_hat^2 (or each component in diagonal modes) before baking into Q_step_sqr. Default 0.0 preserves unclamped behavior.

0.0
calibrate bool | None

Deprecated; True -> "dynamic", False -> "none".

None

Returns:

Type Description
DynamicScanLoopResult

class:DynamicScanLoopResult -- 10 arrays plus the scalar

DynamicScanLoopResult

log-marginal-likelihood. sigma_sqr_seq has shape (N, d_ode)

DynamicScanLoopResult

in the diagonal modes, (N,) otherwise.

ode_filters.filters.ekf1_sqr_loop_preconditioned(mu_0: Array, P_0_sqr: Array, prior: BasePrior, measure: BaseODEInformation, tspan: tuple[float, float], N: int) -> tuple[Array, Array, Array, Array, Array, Array, Array, Array, Array, Array, Array, Array, float]

Run a preconditioned square-root EKF over N observation steps.

Parameters:

Name Type Description Default
mu_0 Array

Initial state mean estimate.

required
P_0_sqr Array

Initial state covariance (square-root form).

required
prior BasePrior

Prior model (typically PrecondIWP).

required
measure BaseODEInformation

Measurement model (e.g., ODEInformation or subclass).

required
tspan tuple[float, float]

Time interval (t_start, t_end).

required
N int

Number of filter steps.

required

Returns:

Type Description
Array

Tuple of 11 arrays, preconditioning matrix, and scalar log marginal

Array

likelihood.

ode_filters.filters.ekf1_sqr_loop_preconditioned_dynamic(mu_0: Array, P_0_sqr: Array, prior: BasePrior, measure: BaseODEInformation, tspan: tuple[float, float], N: int, *, calibration: str = 'dynamic', min_sigma_sqr: float = 0.0, calibrate: bool | None = None) -> PrecondDynamicLoopResult

Preconditioned fixed-step square-root EKF with online calibration.

Same algorithm as :func:ekf1_sqr_loop_dynamic, executed in the preconditioned (bar) space. The calibration estimators are invariant under preconditioning: H Q H.T = (H T) Q_bar (H T).T, so sigma_hat^2 computed in bar space matches the original-space value exactly.

Supported calibration modes: "dynamic", "diagonal", "diagonal_ekf0", "none" -- see :func:ekf1_sqr_adaptive_loop for the full description. ("cumulative" is not exposed for fixed-step loops.)

For the diagonal modes the bar-space time-block Q_bar_time is obtained from prior._Q_bar if present (the PrecondIWP convention, where Q_bar_time is h-independent); otherwise this function raises AttributeError. Other preconditioned priors that expose a different bar-space convention are not supported in the diagonal modes.

Parameters:

Name Type Description Default
mu_0 Array

Initial state mean (original space).

required
P_0_sqr Array

Initial state covariance (square-root form, original space).

required
prior BasePrior

Preconditioned prior (typically PrecondIWP).

required
measure BaseODEInformation

Measurement model (e.g., ODEInformation).

required
tspan tuple[float, float]

Time interval (t_start, t_end).

required
N int

Number of filter steps.

required
calibration str

"dynamic" (default), "diagonal", "diagonal_ekf0", or "none". All modes work with :class:JointPrior / :class:PrecondJointPrior; diagonal modes scale the state block only and require prior._prior_x.xi to be diagonal.

'dynamic'
min_sigma_sqr float

Lower bound applied to per-step sigma_hat^2 before baking into Q_step_sqr_bar. Default 0.0 preserves unclamped behavior.

0.0
calibrate bool | None

Deprecated; True -> "dynamic", False -> "none".

None

Returns:

Type Description
PrecondDynamicLoopResult

Tuple of 12 sequences, the preconditioning matrix T_h, and the

PrecondDynamicLoopResult

scalar log-marginal-likelihood. sigma_sqr_seq entries are

PrecondDynamicLoopResult

length-d_ode arrays for the diagonal modes; floats otherwise.

ode_filters.filters.ekf1_sqr_loop_preconditioned_dynamic_scan(mu_0: Array, P_0_sqr: Array, prior: BasePrior, measure: BaseODEInformation, tspan: tuple[float, float], N: int, *, calibration: str = 'dynamic', min_sigma_sqr: float = 0.0, calibrate: bool | None = None) -> PrecondDynamicScanLoopResult

jax.lax.scan variant of :func:ekf1_sqr_loop_preconditioned_dynamic.

Operates in preconditioned (bar) space. Calibration is invariant under preconditioning: H Q H.T = (HT) Q_bar (HT).T, so all four modes behave identically to the original-space variants.

Parameters:

Name Type Description Default
mu_0 Array

Initial state mean (original space).

required
P_0_sqr Array

Initial state covariance (square-root form, original space).

required
prior BasePrior

Preconditioned prior (e.g. PrecondIWP).

required
measure BaseODEInformation

Measurement model.

required
tspan tuple[float, float]

Time interval (t_start, t_end).

required
N int

Number of filter steps.

required
calibration str

"dynamic" (default), "diagonal", "diagonal_ekf0", or "none". All modes work with :class:PrecondJointPrior; diagonal modes require prior._prior_x.xi to be diagonal.

'dynamic'
min_sigma_sqr float

Lower bound applied to per-step sigma_hat^2 before baking into Q_step_sqr_bar. Default 0.0.

0.0
calibrate bool | None

Deprecated; True -> "dynamic", False -> "none".

None

Returns:

Type Description
PrecondDynamicScanLoopResult

class:PrecondDynamicScanLoopResult -- 12 arrays, T_h, and the

PrecondDynamicScanLoopResult

scalar log-marginal-likelihood. sigma_sqr_seq shape is

PrecondDynamicScanLoopResult

(N, d_ode) in diagonal modes, (N,) otherwise.

ode_filters.filters.ekf1_sqr_loop_preconditioned_sequential(mu_0: Array, P_0_sqr: Array, prior: BasePrior, measure: BaseODEInformation, tspan: tuple[float, float], N: int, *, observations: list[Measurement] | None = None) -> tuple

Run a preconditioned sequential square-root EKF over N steps.

Parameters:

Name Type Description Default
mu_0 Array

Initial state mean estimate.

required
P_0_sqr Array

Initial state covariance (square-root form).

required
prior BasePrior

Prior model (typically PrecondIWP).

required
measure BaseODEInformation

Measurement model (e.g., ODEInformation or subclass).

required
tspan tuple[float, float]

Time interval (t_start, t_end).

required
N int

Number of filter steps.

required
observations list[Measurement] | None

Optional list of :class:Measurement constraints.

None

Returns:

Type Description
tuple

Tuple of arrays, preconditioning matrix, and two scalars

tuple

(log_likelihood_ode, log_likelihood_obs).

ode_filters.filters.ekf1_sqr_loop_preconditioned_sequential_scan(mu_0: Array, P_0_sqr: Array, prior: BasePrior, measure: BaseODEInformation, tspan: tuple[float, float], N: int, obs_model: ObsModel | None = None) -> SeqPrecondScanLoopResult

Run a preconditioned sequential square-root EKF using jax.lax.scan.

Parameters:

Name Type Description Default
mu_0 Array

Initial state mean estimate.

required
P_0_sqr Array

Initial state covariance (square-root form).

required
prior BasePrior

Prior model (typically PrecondIWP).

required
measure BaseODEInformation

Measurement model (ODE + Conservation only).

required
tspan tuple[float, float]

Time interval (t_start, t_end).

required
N int

Number of filter steps.

required
obs_model ObsModel | None

Pre-computed observation data from :func:prepare_observations, or None for ODE-only.

None

Returns:

Type Description
SeqPrecondScanLoopResult

Tuple of 13 arrays, preconditioning matrix, and two scalar

SeqPrecondScanLoopResult

log-likelihoods (log_likelihood_ode, log_likelihood_obs).

ode_filters.filters.ekf1_sqr_loop_sequential(mu_0: Array, Sigma_0_sqr: Array, prior: BasePrior, measure: BaseODEInformation, tspan: tuple[float, float], N: int, *, observations: list[Measurement] | None = None) -> SeqLoopResult

Run a sequential square-root EKF over N steps.

At each step the ODE + Conservation update is applied first. Then, if observations are provided, an observation update is applied using the ODE-updated state as the prior.

Parameters:

Name Type Description Default
mu_0 Array

Initial state mean estimate.

required
Sigma_0_sqr Array

Initial state covariance (square-root form).

required
prior BasePrior

Prior model (e.g., IWP or PrecondIWP).

required
measure BaseODEInformation

Measurement model (e.g., ODEInformation or subclass).

required
tspan tuple[float, float]

Time interval (t_start, t_end).

required
N int

Number of filter steps.

required
observations list[Measurement] | None

Optional list of :class:Measurement constraints. When provided, build_obs_at_time is called at each step to construct an (H, c, R_sqr) tuple for the observation update.

None

Returns:

Type Description
SeqLoopResult

Tuple of 11 lists and two scalars

SeqLoopResult

(log_likelihood_ode, log_likelihood_obs).

ode_filters.filters.ekf1_sqr_loop_sequential_scan(mu_0: Array, Sigma_0_sqr: Array, prior: BasePrior, measure: BaseODEInformation, tspan: tuple[float, float], N: int, obs_model: ObsModel | None = None) -> SeqScanLoopResult

Run a sequential square-root EKF using jax.lax.scan.

At each step the ODE + Conservation update is applied first. Then, if obs_model is provided, an observation update is applied using the ODE-updated state as the prior. Inactive observation steps are masked via jax.lax.select so that the state and log-likelihood are unaffected.

Parameters:

Name Type Description Default
mu_0 Array

Initial state mean estimate.

required
Sigma_0_sqr Array

Initial state covariance (square-root form).

required
prior BasePrior

Prior model (e.g., IWP).

required
measure BaseODEInformation

Measurement model (ODE + Conservation only, no :class:Measurement constraints bundled in).

required
tspan tuple[float, float]

Time interval (t_start, t_end).

required
N int

Number of filter steps.

required
obs_model ObsModel | None

Pre-computed observation data from :func:prepare_observations, or None for ODE-only.

None

Returns:

Type Description
SeqScanLoopResult

Tuple of 11 arrays and two scalar log-likelihoods

SeqScanLoopResult

(log_likelihood_ode, log_likelihood_obs).

ode_filters.filters.rts_sqr_smoother_loop(m_N: Array, P_N_sqr: Array, G_back_seq: Array, d_back_seq: Array, P_back_seq_sqr: Array, N: int) -> tuple[Array, Array]

Run a Rauch-Tung-Striebel smoother over N steps.

Parameters:

Name Type Description Default
m_N Array

Final filtered state mean.

required
P_N_sqr Array

Final filtered state covariance (square-root form).

required
G_back_seq Array

Backward pass gain sequence from filter.

required
d_back_seq Array

Backward pass offset sequence from filter.

required
P_back_seq_sqr Array

Backward pass covariance sequence (square-root form) from filter.

required
N int

Number of smoothing steps.

required

Returns:

Type Description
tuple[Array, Array]

Tuple of smoothed state means and covariances (square-root form).

ode_filters.filters.rts_sqr_smoother_loop_preconditioned(m_N: Array, P_N_sqr: Array, m_N_bar: Array, P_N_sqr_bar: Array, G_back_seq_bar: Array, d_back_seq_bar: Array, P_back_seq_sqr_bar: Array, N: int, T_h: Array) -> tuple[Array, Array]

Run a preconditioned Rauch-Tung-Striebel smoother over N steps.

Parameters:

Name Type Description Default
m_N Array

Final filtered state mean (original space).

required
P_N_sqr Array

Final filtered state covariance (square-root form, original space).

required
m_N_bar Array

Final filtered state mean (preconditioned space).

required
P_N_sqr_bar Array

Final filtered state covariance (square-root form, preconditioned space).

required
G_back_seq_bar Array

Backward pass gain sequence (preconditioned).

required
d_back_seq_bar Array

Backward pass offset sequence (preconditioned).

required
P_back_seq_sqr_bar Array

Backward pass covariance sequence (square-root form, preconditioned).

required
N int

Number of smoothing steps.

required
T_h Array

Preconditioning transformation matrix.

required

Returns:

Type Description
tuple[Array, Array]

Smoothed state means and covariances (square-root form, original space).

ode_filters.filters.ekf1_sqr_filter_step(A_t: Array, b_t: Array, Q_t_sqr: Array, m_prev: Array, P_prev_sqr: Array, measure: BaseODEInformation, t: float = 0.0) -> FilterStepResult

Perform a single square-root EKF prediction and update step.

Parameters:

Name Type Description Default
A_t Array

State transition matrix for current step.

required
b_t Array

Drift vector for current step.

required
Q_t_sqr Array

Square-root of process noise covariance.

required
m_prev Array

Previous state mean estimate.

required
P_prev_sqr Array

Previous state covariance (square-root form).

required
measure BaseODEInformation

Measurement model (e.g., ODEInformation or subclass).

required
t float

Current time (default 0.0).

0.0

Returns:

Type Description
FilterStepResult

Tuple of 4 tuples containing prediction, backward pass, and update results.

ode_filters.filters.ekf1_sqr_filter_step_preconditioned(A_bar: Array, b_bar: Array, Q_sqr_bar: Array, T_t: Array, m_prev_bar: Array, P_prev_sqr_bar: Array, measure: BaseODEInformation, t: float = 0.0) -> PreconditionedFilterStepResult

Perform a single preconditioned square-root EKF step.

Parameters:

Name Type Description Default
A_bar Array

Stepsize-independent state transition matrix.

required
b_bar Array

Stepsize-independent drift vector.

required
Q_sqr_bar Array

Square-root of stepsize-independent process noise covariance.

required
T_t Array

Preconditioning transformation matrix for current step.

required
m_prev_bar Array

Previous state mean estimate (preconditioned space).

required
P_prev_sqr_bar Array

Previous state covariance (square-root form, preconditioned space).

required
measure BaseODEInformation

Measurement model (e.g., ODEInformation or subclass).

required
t float

Current time (default 0.0).

0.0

Returns:

Type Description
PreconditionedFilterStepResult

Tuple of 5 tuples with preconditioned and original-space results.

ode_filters.filters.ekf1_sqr_filter_step_preconditioned_sequential(A_bar: Array, b_bar: Array, Q_sqr_bar: Array, T_t: Array, m_prev_bar: Array, P_prev_sqr_bar: Array, measure: BaseODEInformation, t: float = 0.0, *, obs: tuple[Array, Array, Array] | None = None) -> SeqPrecondFilterStepResult

Perform a sequential preconditioned square-root EKF step.

Parameters:

Name Type Description Default
A_bar Array

Stepsize-independent state transition matrix.

required
b_bar Array

Stepsize-independent drift vector.

required
Q_sqr_bar Array

Square-root of stepsize-independent process noise covariance.

required
T_t Array

Preconditioning transformation matrix for current step.

required
m_prev_bar Array

Previous state mean estimate (preconditioned space).

required
P_prev_sqr_bar Array

Previous state covariance (square-root form, preconditioned).

required
measure BaseODEInformation

Measurement model (e.g., ODEInformation or subclass).

required
t float

Current time (default 0.0).

0.0
obs tuple[Array, Array, Array] | None

Optional observation tuple (H_obs, c_obs, R_obs_sqr).

None

Returns:

Type Description
SeqPrecondFilterStepResult

Tuple of 6 tuples with preconditioned and original-space results.

ode_filters.filters.ekf1_sqr_filter_step_preconditioned_sequential_scan(A_bar: Array, b_bar: Array, Q_sqr_bar: Array, T_t: Array, m_prev_bar: Array, P_prev_sqr_bar: Array, measure: BaseODEInformation, t: float, H_obs: Array, c_obs: Array, R_obs_sqr: Array, obs_active: Array) -> tuple[tuple[Array, Array], tuple[Array, Array, Array], tuple[Array, Array], tuple[Array, Array], tuple[Array, Array], tuple[Array, Array]]

Scan-compatible preconditioned sequential square-root EKF step.

Parameters:

Name Type Description Default
A_bar Array

Stepsize-independent state transition matrix.

required
b_bar Array

Stepsize-independent drift vector.

required
Q_sqr_bar Array

Square-root stepsize-independent process noise.

required
T_t Array

Preconditioning transformation matrix.

required
m_prev_bar Array

Previous state mean (preconditioned space).

required
P_prev_sqr_bar Array

Previous state covariance (preconditioned, square-root).

required
measure BaseODEInformation

Measurement model (ODE + Conservation only).

required
t float

Current time.

required
H_obs Array

Observation Jacobian, shape [obs_dim, state_dim].

required
c_obs Array

Observation offset, shape [obs_dim].

required
R_obs_sqr Array

Square-root observation noise, shape [obs_dim, obs_dim].

required
obs_active Array

Scalar boolean — whether observations are active.

required

Returns:

Type Description
tuple[tuple[Array, Array], tuple[Array, Array, Array], tuple[Array, Array], tuple[Array, Array], tuple[Array, Array], tuple[Array, Array]]

Tuple of 6 tuples with preconditioned and original-space results.

ode_filters.filters.ekf1_sqr_filter_step_sequential(A_t: Array, b_t: Array, Q_t_sqr: Array, m_prev: Array, P_prev_sqr: Array, measure: BaseODEInformation, t: float = 0.0, *, obs: tuple[Array, Array, Array] | None = None) -> SeqFilterStepResult

Perform a sequential square-root EKF step: ODE update then observation.

Unlike the joint update in ekf1_sqr_filter_step, this function first incorporates ODE + Conservation information, then applies an optional observation update using the ODE-updated state as the prior.

Parameters:

Name Type Description Default
A_t Array

State transition matrix for current step.

required
b_t Array

Drift vector for current step.

required
Q_t_sqr Array

Square-root of process noise covariance.

required
m_prev Array

Previous state mean estimate.

required
P_prev_sqr Array

Previous state covariance (square-root form).

required
measure BaseODEInformation

Measurement model (e.g., ODEInformation or subclass).

required
t float

Current time (default 0.0).

0.0
obs tuple[Array, Array, Array] | None

Optional observation tuple (H_obs, c_obs, R_obs_sqr). When provided, the observation update uses the ODE-updated state (m_ode, P_ode_sqr) as its prior.

None

Returns:

Type Description
SeqFilterStepResult

Tuple of 5 tuples: prediction, backward pass, ODE observation

SeqFilterStepResult

marginal, observation marginal, and final updated state.

ode_filters.filters.ekf1_sqr_filter_step_sequential_scan(A_t: Array, b_t: Array, Q_t_sqr: Array, m_prev: Array, P_prev_sqr: Array, measure: BaseODEInformation, t: float, H_obs: Array, c_obs: Array, R_obs_sqr: Array, obs_active: Array) -> tuple[tuple[Array, Array], tuple[Array, Array, Array], tuple[Array, Array], tuple[Array, Array], tuple[Array, Array]]

Scan-compatible sequential square-root EKF step.

Always executes the observation update with fixed-shape arrays. When obs_active is False, the ODE-only state is selected via jnp.where so that observations have no effect.

Parameters:

Name Type Description Default
A_t Array

State transition matrix.

required
b_t Array

Drift vector.

required
Q_t_sqr Array

Square-root process noise covariance.

required
m_prev Array

Previous state mean.

required
P_prev_sqr Array

Previous state covariance (square-root form).

required
measure BaseODEInformation

Measurement model (ODE + Conservation only).

required
t float

Current time.

required
H_obs Array

Observation Jacobian, shape [obs_dim, state_dim].

required
c_obs Array

Observation offset, shape [obs_dim].

required
R_obs_sqr Array

Square-root observation noise, shape [obs_dim, obs_dim].

required
obs_active Array

Scalar boolean — whether observations are active.

required

Returns:

Type Description
tuple[Array, Array]

Tuple of 5 tuples: prediction, backward pass, ODE observation

tuple[Array, Array, Array]

marginal, observation marginal, and final updated state.

ode_filters.filters.rts_sqr_smoother_step(G_back: Array, d_back: Array, P_back_sqr: Array, m_s: Array, P_s_sqr: Array) -> tuple[Array, Array]

Perform a single Rauch-Tung-Striebel backward smoothing step.

Parameters:

Name Type Description Default
G_back Array

Backward pass gain matrix.

required
d_back Array

Backward pass offset vector.

required
P_back_sqr Array

Backward pass covariance (square-root form).

required
m_s Array

Smoothed state mean from next time step.

required
P_s_sqr Array

Smoothed state covariance (square-root form) from next time step.

required

Returns:

Type Description
tuple[Array, Array]

Tuple of previous smoothed mean and covariance (square-root form).

ode_filters.filters.rts_sqr_smoother_step_preconditioned(G_back_bar: Array, d_back_bar: Array, P_back_sqr_bar: Array, m_s_bar: Array, P_s_sqr_bar: Array, T_t: Array) -> tuple[tuple[Array, Array], tuple[Array, Array]]

Perform a single preconditioned Rauch-Tung-Striebel backward smoothing step.

Parameters:

Name Type Description Default
G_back_bar Array

Backward pass gain matrix (preconditioned space).

required
d_back_bar Array

Backward pass offset vector (preconditioned space).

required
P_back_sqr_bar Array

Backward pass covariance (square-root form, preconditioned space).

required
m_s_bar Array

Smoothed state mean from next step (preconditioned space).

required
P_s_sqr_bar Array

Smoothed state covariance (square-root form, preconditioned space).

required
T_t Array

Preconditioning transformation matrix.

required

Returns:

Type Description
tuple[Array, Array]

Tuple containing:

tuple[Array, Array]
  • (mean, covariance) in preconditioned space
tuple[tuple[Array, Array], tuple[Array, Array]]
  • (mean, covariance) in original space

ode_filters.filters.adaptive_controller

Step-size controllers for the adaptive EKF loop.

A controller's job is to propose a new step size h_new given the current step h and a normalised error estimate err (one means "at tolerance", greater than one means "rejected, shrink"). Accept/reject is decided by the loop driver, not by the controller.

Two controllers are provided:

  • :class:PController -- proportional-only on the current error. Sometimes called the I-controller or elementary controller in the adaptive-ODE-solver literature (Hairer-Wanner-Norsett, Soederlind) because h itself is the integral of the log-step adjustment.
  • :class:PIController -- adds a Gustafsson PI term using the previous accepted step's error to damp oscillations around the tolerance.

Both implement the :class:StepSizeController protocol.

ode_filters.filters.adaptive_controller.StepSizeController

Bases: Protocol

Duck-typed interface for step-size controllers.

err_prev may be None to indicate the absence of memory (first step, or right after a reject).

ode_filters.filters.adaptive_controller.StepSizeController.propose(h: float, err: float, err_prev: float | None) -> float

Return the proposed next step size.

ode_filters.filters.adaptive_controller.PController dataclass

Proportional (single-error) step-size controller.

The proposed step is

h_new = h * safety * err^(-alpha)

clipped to [min_factor, max_factor] * h. The previous error is ignored, so this controller has no memory.

Default alpha = 1 / order matches the standard "I-controller" of Hairer-Wanner-Norsett (1993) and the elementary controller in scipy's solve_ivp. For an EKF1 with IWP(q) prior, order = q.

Parameters:

Name Type Description Default
order int

Convergence order of the local error estimate.

required
safety float

Safety factor on the proposed step.

0.9
alpha float | None

Gain on the current error; defaults to 1.0 / order.

None
min_factor float

Lower clip on h_new / h.

0.2
max_factor float

Upper clip on h_new / h.

5.0

ode_filters.filters.adaptive_controller.PController.propose(h: float, err: float, err_prev: float | None = None) -> float

Return the next proposed step size.

Parameters:

Name Type Description Default
h float

Current step size.

required
err float

Normalised error of the current step (1.0 = at tolerance).

required
err_prev float | None

Ignored; accepted for protocol compatibility with :class:PIController.

None

ode_filters.filters.adaptive_controller.PIController dataclass

Gustafsson-style proportional-integral controller.

The proposed step is

h_new = h * safety * err^(-alpha) * (err_prev / err)^(beta)

clipped to [min_factor, max_factor] * h. When err_prev is unknown (e.g. on the first step or right after a reject) the I-term is dropped and the update reduces to the :class:PController form h_new = h * safety * err^(-alpha).

Defaults follow Gustafsson (1991): alpha = 0.7 / order, beta = 0.4 / order, which are also the values used in Bosch et al. (2021). For an EKF1 with IWP(q) prior, order = q.

Parameters:

Name Type Description Default
order int

Convergence order of the local error estimate.

required
safety float

Safety factor on the proposed step.

0.9
alpha float | None

Proportional gain; defaults to 0.7 / order.

None
beta float | None

Integral gain; defaults to 0.4 / order. Set beta=0 for a pure proportional controller -- :class:PController is the standalone equivalent with the more natural default alpha.

None
min_factor float

Lower clip on h_new / h.

0.2
max_factor float

Upper clip on h_new / h.

5.0

ode_filters.filters.adaptive_controller.PIController.propose(h: float, err: float, err_prev: float | None) -> float

Return the next proposed step size.

Parameters:

Name Type Description Default
h float

Current step size.

required
err float

Normalised error of the current step (1.0 = at tolerance).

required
err_prev float | None

Normalised error of the previously accepted step, or None if unavailable.

required

ode_filters.filters.ode_filter_adaptive

Adaptive-step EKF loop with online sigma calibration.

The trajectory length is data-dependent (steps may be rejected and retried at a smaller h), so this loop cannot use jax.lax.scan. Instead it is a Python while driver around a jitted per-step body. The per-step body does the prediction + update + per-step quasi-MLE sigma and returns the normalised local-error estimate; the Python layer makes the accept/reject decision.

Four calibration schemes are exposed via the calibration parameter:

  • "dynamic" (default, Bosch et al. 2021 Eq. 32): scalar sigma_hat^2_n = m_z.T (H Q(h) H.T)^{-1} m_z / d, baked into the current step's Q_h before propagation. Past steps keep their own sigma_hat^2 -- statistically honest per-step, Markov-preserving. Fails on multi-scale ODEs (the scalar average is dominated by the larger-residual component, which starves smaller-scale components).

  • "diagonal_ekf0": per-component sigma_hat^2_i = m_z[i]^2 / (E_1 Q E_1.T)_ii. Because E_1 selects only the first-derivative block, E_1 Q E_1.T is exactly diagonal when xi is diagonal -- so sigma_hat^2_i is the genuine MLE in a block-diagonal sub-model for component i. This is the standard per-component recipe (Bosch et al. 2021 §4.2, DynamicMVDiffusion) and is the recommended choice for multi-component problems. Requires diagonal xi.

  • "diagonal" (EK1-flavoured heuristic): same formula as "diagonal_ekf0" but with the EK1 Jacobian H_t = E_1 - J_f E_0 in the denominator. H_t Q H_t.T is then generally dense, and taking its diagonal ignores cross-component coupling introduced by J_f. When J_f is block-diagonal in components (e.g. fully decoupled RHS) this collapses to "diagonal_ekf0"; when J_f couples components strongly, the diagonal-of-dense recipe biases the per-component MLE. Use only when you understand the trade-off; prefer "diagonal_ekf0" by default. Requires diagonal xi.

  • "cumulative" (legacy multiplicative scheme): scalar sigma_hat^2_n from the full noise-free predicted-obs covariance H P_pred H.T; post-multiplies the whole step's covariance, so past contributions inherit every subsequent sigma_hat^2. Non-Markovian by design; an online empirical-Bayes update of a single global unknown sigma. Not per-step honest, but stays robust on multi-scale problems (the inflated carried P accidentally compensates). "diagonal_ekf0" is strictly better when applicable; keep "cumulative" for non-diagonal Xi or for backward compatibility with earlier runs.

  • "none": propagate uncalibrated, still report sigma_hat^2 for diagnostics / post-hoc rescaling.

A separate sigma_in_error parameter controls how sigma_hat^2 enters the local-error estimate used for step-size selection. The quasi-MLE estimator has chi-squared noise (relative std sqrt(2/d)) which propagates into h decisions, causing visible oscillations. The default sigma_in_error="running_mean" substitutes the cumulative running mean in the error formula -- noise vanishes after a few accepted steps, h becomes much smoother, reject counts drop. Per-step sigma_hat^2 is still used in the actual Q calibration. Pass sigma_in_error="per_step" to recover the original Bosch et al. 2021 recipe. Caveat: the running mean lags real changes in problem stiffness; in regimes where the local diffusion changes abruptly (entering a stiff transition), "per_step" is more responsive.

The accepted-step outputs are stored in lists matching the shape conventions of :func:ekf1_sqr_loop, so the result can be passed directly to :func:rts_sqr_smoother_loop.

The returned log_likelihood is the post-calibration marginal likelihood: each step's contribution uses the calibrated Pz_sqr (which includes sigma_hat^2 in the dynamic/diagonal modes and the post-multiplied scale in cumulative mode). It is therefore the right quantity for inference given the chosen calibration model, but is not directly comparable across calibration modes -- different calibration settings define different generative models.

ode_filters.filters.ode_filter_adaptive.AdaptiveLoopResult

Bases: NamedTuple

Output of :func:ekf1_sqr_adaptive_loop.

All sequences are aligned to accepted steps. t_seq has length N_accepted + 1 (initial time plus one entry per accepted step); the per-step sequences (m_pred_seq, Pz_seq_sqr, ...) and the smoothing quantities (G_back_seq, d_back_seq, P_back_seq_sqr) have length N_accepted. Filtered means/covariances (m_seq, P_seq_sqr) have length N_accepted + 1 and include the initial state.

ode_filters.filters.ode_filter_adaptive.ekf1_sqr_adaptive_loop(mu_0: Array, Sigma_0_sqr: Array, prior: BasePrior, measure: BaseODEInformation, tspan: tuple[float, float], *, atol: float = 0.0001, rtol: float = 0.01, h_init: float | None = None, h_min: float = 1e-10, h_max: float | None = None, controller: StepSizeController | None = None, calibration: CalibrationMode = 'dynamic', sigma_in_error: SigmaInError = 'running_mean', min_sigma_sqr: float = 0.0, calibrate: bool | None = None, max_steps: int = 100000) -> AdaptiveLoopResult

Adaptive-step square-root EKF with per-step diffusion calibration.

Python while driver around a jitted per-step body. Every accepted step contributes to the returned sequences and the result carries the smoothing-relevant outputs, so it can be fed directly to :func:rts_sqr_smoother_loop.

Parameters:

Name Type Description Default
mu_0 Array

Initial state mean.

required
Sigma_0_sqr Array

Initial state covariance (square-root form).

required
prior BasePrior

Prior model (e.g. :class:IWP). The Kronecker structure Q(h) = kron(_Q(h), xi) is used to build the per-step process noise. prior.q sets the default controller order.

required
measure BaseODEInformation

Measurement model (e.g. :class:ODEInformation).

required
tspan tuple[float, float]

Time interval (t_start, t_end). Must be a tuple (hashable for jit-compatible inner step).

required
atol float

Absolute tolerance on the function-value local error.

0.0001
rtol float

Relative tolerance.

0.01
h_init float | None

Initial step. Defaults to (t_end - t_start) / 100.

None
h_min float

Minimum step. Steps below this raise RuntimeError.

1e-10
h_max float | None

Maximum step. Defaults to t_end - t_start.

None
controller StepSizeController | None

Step controller implementing :class:~ode_filters.filters.adaptive_controller.StepSizeController. Defaults to PIController(order=prior.q); pass :class:~ode_filters.filters.adaptive_controller.PController for a memoryless proportional-only controller.

None
calibration CalibrationMode

How the per-step sigma_hat^2 enters the stored posterior. See the module docstring for the four available modes:

  • "dynamic" (default) -- Bosch et al. 2021 Eq. 32. Scalar sigma_hat^2 from H Q(h) H.T baked into Q_h. Honest per-step but scalar: fails on multi-scale ODE systems (one component's residual dominates sigma_hat^2 and starves the others).
  • "diagonal_ekf0" -- per-component sigma_hat^2_i = m_z[i]^2 / (E_1 Q E_1.T)_ii using the EK0 observation matrix H_0 = E_1. The denominator is exactly diagonal when xi is diagonal, so the per-component estimator is the genuine MLE in a block-diagonal sub-model. The recommended choice for multi-component problems. Requires prior.xi diagonal.
  • "diagonal" -- same formula, but with the EK1 Jacobian H_t = E_1 - J_f E_0 in the denominator. Heuristic: taking the diagonal of a generally-dense H_t Q(h) H_t.T ignores cross-component coupling from J_f. Equivalent to "diagonal_ekf0" when J_f is block-diagonal (e.g. fully decoupled RHS); biased in proportion to the off-diagonal entries of J_f otherwise. Requires prior.xi diagonal.
  • "cumulative" -- legacy multiplicative scheme; scalar sigma from the full H P_pred H.T post-multiplies the whole step. Non-Markovian; robust on multi-scale by accident (inflated P inherited from earlier transitions). Prefer "diagonal_ekf0".
  • "none" -- propagate sigma=1, still report sigma_hat^2 for post-hoc rescaling / diagnostics.

All four modes work with :class:JointPrior / :class:PrecondJointPrior: the diagonal modes scale only the state block per-component and require prior._prior_x.xi to be diagonal.

For the diagonal modes the per-step entries of result.sigma_sqr_seq are length-d arrays; for the scalar modes they are Python floats.

'dynamic'
sigma_in_error SigmaInError

How sigma_hat^2 enters the local-error estimate.

  • "running_mean" (default) -- use the cumulative running mean of past accepted-step sigma_hat^2 values in the error formula (sigma_hat^2 is still per-step in the Q calibration). The running mean has variance 2/(nd), so estimator noise vanishes after a few accepted steps -- h becomes much smoother and reject counts typically drop by several times. For "diagonal" modes the running mean is per-component.
  • "per_step" -- use the just-computed sigma_hat^2 (Bosch et al. 2021 original recipe). Honest but inherits the chi-squared estimator noise (relative std sqrt(2/d)), which propagates into step-size decisions and causes visible h oscillations in flat regions.
'running_mean'
min_sigma_sqr float

Lower bound applied to the per-step sigma_hat^2 (and to each component of the per-component vector in diagonal modes) before it is baked into Q_step_sqr. Default 0.0 preserves the unclamped behavior. Use a small positive value (e.g. 1e-30) on problems where the trivial zero-residual fixed point would otherwise collapse the state-block diffusion to zero and propagate NaN.

0.0
calibrate bool | None

Deprecated. True maps to calibration="dynamic", False maps to calibration="none". Pass calibration directly instead.

None
max_steps int

Hard cap on iterations (rejected + accepted) as a safety valve against infinite loops.

100000

Returns:

Type Description
AdaptiveLoopResult

class:AdaptiveLoopResult with the accepted trajectory.

Raises:

Type Description
RuntimeError

If the controller proposes a step below h_min or the iteration cap is exceeded.

ValueError

If tspan is non-increasing or calibration / sigma_in_error is unknown.

ode_filters.filters.ode_filter_loop

ode_filters.filters.ode_filter_loop.ekf1_sqr_loop(mu_0: Array, Sigma_0_sqr: Array, prior: BasePrior, measure: BaseODEInformation, tspan: tuple[float, float], N: int) -> LoopResult

Run a square-root EKF over N observation steps.

Parameters:

Name Type Description Default
mu_0 Array

Initial state mean estimate.

required
Sigma_0_sqr Array

Initial state covariance (square-root form).

required
prior BasePrior

Prior model (e.g., IWP or PrecondIWP).

required
measure BaseODEInformation

Measurement model (e.g., ODEInformation or subclass).

required
tspan tuple[float, float]

Time interval (t_start, t_end).

required
N int

Number of filter steps.

required

Returns:

Type Description
LoopResult

Tuple of 9 arrays and a scalar log marginal likelihood.

ode_filters.filters.ode_filter_loop.ekf1_sqr_loop_dynamic(mu_0: Array, Sigma_0_sqr: Array, prior: BasePrior, measure: BaseODEInformation, tspan: tuple[float, float], N: int, *, calibration: str = 'dynamic', min_sigma_sqr: float = 0.0, calibrate: bool | None = None) -> DynamicLoopResult

Fixed-step square-root EKF with online per-step diffusion calibration.

Same fixed grid as :func:ekf1_sqr_loop, but at every step the per-step quasi-MLE diffusion is estimated and -- depending on calibration -- baked into the current step's Q_h before propagation. See :func:~ode_filters.filters.ekf1_sqr_adaptive_loop for the full description of the modes; this loop supports the same names:

  • "dynamic" (default, scalar): Bosch et al. 2021 Eq. 32. P_pred = A P_prev A.T + sigma_hat^2_n * Q_h.
  • "diagonal_ekf0": per-component sigma_hat^2_i from (E_1 Q E_1.T)_ii (exact-diagonal denominator). Genuine MLE in a block-diagonal sub-model for component i. Recommended default for multi-component problems.
  • "diagonal": same per-component formula but with the EK1 Jacobian H_t = E_1 - J_f E_0 in the denominator. Heuristic -- equivalent to "diagonal_ekf0" when J_f is block-diagonal in components, biased in proportion to J_f's off-diagonal entries otherwise.
  • "none": propagate sigma=1, report sigma_hat^2 for post-hoc rescaling (combine with :func:posthoc_mle_sigma_sqr / :func:rescale_sqr_seq).

"cumulative" is intentionally not exposed on the fixed-step loop: the cumulative scheme's appeal is robust adaptive-step control, and at fixed step there's nothing to gain over post-hoc rescaling.

The returned log_likelihood reflects the calibrated posterior -- each step's contribution uses the calibrated Pz_sqr. With calibration="none" it matches the log-likelihood from :func:ekf1_sqr_loop; in the calibrated modes it is the log-likelihood under the corresponding generative model and is not directly comparable across modes.

Parameters:

Name Type Description Default
mu_0 Array

Initial state mean.

required
Sigma_0_sqr Array

Initial state covariance (square-root form).

required
prior BasePrior

Prior model (e.g. IWP).

required
measure BaseODEInformation

Measurement model (e.g. ODEInformation).

required
tspan tuple[float, float]

Time interval (t_start, t_end).

required
N int

Number of filter steps.

required
calibration str

"dynamic" (default), "diagonal", or "none". All modes work with :class:JointPrior / :class:PrecondJointPrior; diagonal modes scale the state block only and require prior._prior_x.xi to be diagonal.

'dynamic'
min_sigma_sqr float

Lower bound on the per-step sigma_hat^2 (and on each component in diagonal modes) before it is baked into Q_step_sqr. Default 0.0 preserves unclamped behavior.

0.0
calibrate bool | None

Deprecated; True -> "dynamic", False -> "none".

None

Returns:

Type Description
DynamicLoopResult

Tuple of 10 sequences plus the scalar log-marginal-likelihood:

DynamicLoopResult

``(m_seq, P_seq_sqr, m_pred_seq, P_pred_seq_sqr, G_back_seq,

DynamicLoopResult

d_back_seq, P_back_seq_sqr, mz_seq, Pz_seq_sqr, sigma_sqr_seq,

DynamicLoopResult

log_likelihood). In the diagonal modessigma_sqr_seq``

DynamicLoopResult

entries are length-d_ode arrays; otherwise Python floats.

ode_filters.filters.ode_filter_loop.rts_sqr_smoother_loop(m_N: Array, P_N_sqr: Array, G_back_seq: Array, d_back_seq: Array, P_back_seq_sqr: Array, N: int) -> tuple[Array, Array]

Run a Rauch-Tung-Striebel smoother over N steps.

Parameters:

Name Type Description Default
m_N Array

Final filtered state mean.

required
P_N_sqr Array

Final filtered state covariance (square-root form).

required
G_back_seq Array

Backward pass gain sequence from filter.

required
d_back_seq Array

Backward pass offset sequence from filter.

required
P_back_seq_sqr Array

Backward pass covariance sequence (square-root form) from filter.

required
N int

Number of smoothing steps.

required

Returns:

Type Description
tuple[Array, Array]

Tuple of smoothed state means and covariances (square-root form).

ode_filters.filters.ode_filter_loop.ekf1_sqr_loop_preconditioned(mu_0: Array, P_0_sqr: Array, prior: BasePrior, measure: BaseODEInformation, tspan: tuple[float, float], N: int) -> tuple[Array, Array, Array, Array, Array, Array, Array, Array, Array, Array, Array, Array, float]

Run a preconditioned square-root EKF over N observation steps.

Parameters:

Name Type Description Default
mu_0 Array

Initial state mean estimate.

required
P_0_sqr Array

Initial state covariance (square-root form).

required
prior BasePrior

Prior model (typically PrecondIWP).

required
measure BaseODEInformation

Measurement model (e.g., ODEInformation or subclass).

required
tspan tuple[float, float]

Time interval (t_start, t_end).

required
N int

Number of filter steps.

required

Returns:

Type Description
Array

Tuple of 11 arrays, preconditioning matrix, and scalar log marginal

Array

likelihood.

ode_filters.filters.ode_filter_loop.ekf1_sqr_loop_preconditioned_dynamic(mu_0: Array, P_0_sqr: Array, prior: BasePrior, measure: BaseODEInformation, tspan: tuple[float, float], N: int, *, calibration: str = 'dynamic', min_sigma_sqr: float = 0.0, calibrate: bool | None = None) -> PrecondDynamicLoopResult

Preconditioned fixed-step square-root EKF with online calibration.

Same algorithm as :func:ekf1_sqr_loop_dynamic, executed in the preconditioned (bar) space. The calibration estimators are invariant under preconditioning: H Q H.T = (H T) Q_bar (H T).T, so sigma_hat^2 computed in bar space matches the original-space value exactly.

Supported calibration modes: "dynamic", "diagonal", "diagonal_ekf0", "none" -- see :func:ekf1_sqr_adaptive_loop for the full description. ("cumulative" is not exposed for fixed-step loops.)

For the diagonal modes the bar-space time-block Q_bar_time is obtained from prior._Q_bar if present (the PrecondIWP convention, where Q_bar_time is h-independent); otherwise this function raises AttributeError. Other preconditioned priors that expose a different bar-space convention are not supported in the diagonal modes.

Parameters:

Name Type Description Default
mu_0 Array

Initial state mean (original space).

required
P_0_sqr Array

Initial state covariance (square-root form, original space).

required
prior BasePrior

Preconditioned prior (typically PrecondIWP).

required
measure BaseODEInformation

Measurement model (e.g., ODEInformation).

required
tspan tuple[float, float]

Time interval (t_start, t_end).

required
N int

Number of filter steps.

required
calibration str

"dynamic" (default), "diagonal", "diagonal_ekf0", or "none". All modes work with :class:JointPrior / :class:PrecondJointPrior; diagonal modes scale the state block only and require prior._prior_x.xi to be diagonal.

'dynamic'
min_sigma_sqr float

Lower bound applied to per-step sigma_hat^2 before baking into Q_step_sqr_bar. Default 0.0 preserves unclamped behavior.

0.0
calibrate bool | None

Deprecated; True -> "dynamic", False -> "none".

None

Returns:

Type Description
PrecondDynamicLoopResult

Tuple of 12 sequences, the preconditioning matrix T_h, and the

PrecondDynamicLoopResult

scalar log-marginal-likelihood. sigma_sqr_seq entries are

PrecondDynamicLoopResult

length-d_ode arrays for the diagonal modes; floats otherwise.

ode_filters.filters.ode_filter_loop.rts_sqr_smoother_loop_preconditioned(m_N: Array, P_N_sqr: Array, m_N_bar: Array, P_N_sqr_bar: Array, G_back_seq_bar: Array, d_back_seq_bar: Array, P_back_seq_sqr_bar: Array, N: int, T_h: Array) -> tuple[Array, Array]

Run a preconditioned Rauch-Tung-Striebel smoother over N steps.

Parameters:

Name Type Description Default
m_N Array

Final filtered state mean (original space).

required
P_N_sqr Array

Final filtered state covariance (square-root form, original space).

required
m_N_bar Array

Final filtered state mean (preconditioned space).

required
P_N_sqr_bar Array

Final filtered state covariance (square-root form, preconditioned space).

required
G_back_seq_bar Array

Backward pass gain sequence (preconditioned).

required
d_back_seq_bar Array

Backward pass offset sequence (preconditioned).

required
P_back_seq_sqr_bar Array

Backward pass covariance sequence (square-root form, preconditioned).

required
N int

Number of smoothing steps.

required
T_h Array

Preconditioning transformation matrix.

required

Returns:

Type Description
tuple[Array, Array]

Smoothed state means and covariances (square-root form, original space).

ode_filters.filters.ode_filter_loop.ekf1_sqr_loop_sequential(mu_0: Array, Sigma_0_sqr: Array, prior: BasePrior, measure: BaseODEInformation, tspan: tuple[float, float], N: int, *, observations: list[Measurement] | None = None) -> SeqLoopResult

Run a sequential square-root EKF over N steps.

At each step the ODE + Conservation update is applied first. Then, if observations are provided, an observation update is applied using the ODE-updated state as the prior.

Parameters:

Name Type Description Default
mu_0 Array

Initial state mean estimate.

required
Sigma_0_sqr Array

Initial state covariance (square-root form).

required
prior BasePrior

Prior model (e.g., IWP or PrecondIWP).

required
measure BaseODEInformation

Measurement model (e.g., ODEInformation or subclass).

required
tspan tuple[float, float]

Time interval (t_start, t_end).

required
N int

Number of filter steps.

required
observations list[Measurement] | None

Optional list of :class:Measurement constraints. When provided, build_obs_at_time is called at each step to construct an (H, c, R_sqr) tuple for the observation update.

None

Returns:

Type Description
SeqLoopResult

Tuple of 11 lists and two scalars

SeqLoopResult

(log_likelihood_ode, log_likelihood_obs).

ode_filters.filters.ode_filter_loop.ekf1_sqr_loop_preconditioned_sequential(mu_0: Array, P_0_sqr: Array, prior: BasePrior, measure: BaseODEInformation, tspan: tuple[float, float], N: int, *, observations: list[Measurement] | None = None) -> tuple

Run a preconditioned sequential square-root EKF over N steps.

Parameters:

Name Type Description Default
mu_0 Array

Initial state mean estimate.

required
P_0_sqr Array

Initial state covariance (square-root form).

required
prior BasePrior

Prior model (typically PrecondIWP).

required
measure BaseODEInformation

Measurement model (e.g., ODEInformation or subclass).

required
tspan tuple[float, float]

Time interval (t_start, t_end).

required
N int

Number of filter steps.

required
observations list[Measurement] | None

Optional list of :class:Measurement constraints.

None

Returns:

Type Description
tuple

Tuple of arrays, preconditioning matrix, and two scalars

tuple

(log_likelihood_ode, log_likelihood_obs).

ode_filters.filters.ode_filter_loop.ekf1_sqr_loop_sequential_scan(mu_0: Array, Sigma_0_sqr: Array, prior: BasePrior, measure: BaseODEInformation, tspan: tuple[float, float], N: int, obs_model: ObsModel | None = None) -> SeqScanLoopResult

Run a sequential square-root EKF using jax.lax.scan.

At each step the ODE + Conservation update is applied first. Then, if obs_model is provided, an observation update is applied using the ODE-updated state as the prior. Inactive observation steps are masked via jax.lax.select so that the state and log-likelihood are unaffected.

Parameters:

Name Type Description Default
mu_0 Array

Initial state mean estimate.

required
Sigma_0_sqr Array

Initial state covariance (square-root form).

required
prior BasePrior

Prior model (e.g., IWP).

required
measure BaseODEInformation

Measurement model (ODE + Conservation only, no :class:Measurement constraints bundled in).

required
tspan tuple[float, float]

Time interval (t_start, t_end).

required
N int

Number of filter steps.

required
obs_model ObsModel | None

Pre-computed observation data from :func:prepare_observations, or None for ODE-only.

None

Returns:

Type Description
SeqScanLoopResult

Tuple of 11 arrays and two scalar log-likelihoods

SeqScanLoopResult

(log_likelihood_ode, log_likelihood_obs).

ode_filters.filters.ode_filter_loop.ekf1_sqr_loop_preconditioned_sequential_scan(mu_0: Array, P_0_sqr: Array, prior: BasePrior, measure: BaseODEInformation, tspan: tuple[float, float], N: int, obs_model: ObsModel | None = None) -> SeqPrecondScanLoopResult

Run a preconditioned sequential square-root EKF using jax.lax.scan.

Parameters:

Name Type Description Default
mu_0 Array

Initial state mean estimate.

required
P_0_sqr Array

Initial state covariance (square-root form).

required
prior BasePrior

Prior model (typically PrecondIWP).

required
measure BaseODEInformation

Measurement model (ODE + Conservation only).

required
tspan tuple[float, float]

Time interval (t_start, t_end).

required
N int

Number of filter steps.

required
obs_model ObsModel | None

Pre-computed observation data from :func:prepare_observations, or None for ODE-only.

None

Returns:

Type Description
SeqPrecondScanLoopResult

Tuple of 13 arrays, preconditioning matrix, and two scalar

SeqPrecondScanLoopResult

log-likelihoods (log_likelihood_ode, log_likelihood_obs).

ode_filters.filters.ode_filter_loop.ekf1_sqr_loop_dynamic_scan(mu_0: Array, Sigma_0_sqr: Array, prior: BasePrior, measure: BaseODEInformation, tspan: tuple[float, float], N: int, *, calibration: str = 'dynamic', min_sigma_sqr: float = 0.0, calibrate: bool | None = None) -> DynamicScanLoopResult

jax.lax.scan variant of :func:ekf1_sqr_loop_dynamic.

Same per-step semantics, executed inside a single jax.lax.scan for full JIT/compile-once efficiency. Returns stacked Array sequences rather than lists. Supports calibration in {"dynamic", "diagonal", "diagonal_ekf0", "none"}; in the diagonal modes sigma_sqr_seq has shape (N, d) rather than (N,).

Parameters:

Name Type Description Default
mu_0 Array

Initial state mean.

required
Sigma_0_sqr Array

Initial state covariance (square-root form).

required
prior BasePrior

Prior (e.g. IWP).

required
measure BaseODEInformation

Measurement model.

required
tspan tuple[float, float]

Time interval (t_start, t_end).

required
N int

Number of filter steps.

required
calibration str

"dynamic" (default), "diagonal", "diagonal_ekf0", or "none". All modes work with :class:JointPrior / :class:PrecondJointPrior; diagonal modes require prior._prior_x.xi to be diagonal.

'dynamic'
min_sigma_sqr float

Lower bound applied to per-step sigma_hat^2 (or each component in diagonal modes) before baking into Q_step_sqr. Default 0.0 preserves unclamped behavior.

0.0
calibrate bool | None

Deprecated; True -> "dynamic", False -> "none".

None

Returns:

Type Description
DynamicScanLoopResult

class:DynamicScanLoopResult -- 10 arrays plus the scalar

DynamicScanLoopResult

log-marginal-likelihood. sigma_sqr_seq has shape (N, d_ode)

DynamicScanLoopResult

in the diagonal modes, (N,) otherwise.

ode_filters.filters.ode_filter_loop.ekf1_sqr_loop_preconditioned_dynamic_scan(mu_0: Array, P_0_sqr: Array, prior: BasePrior, measure: BaseODEInformation, tspan: tuple[float, float], N: int, *, calibration: str = 'dynamic', min_sigma_sqr: float = 0.0, calibrate: bool | None = None) -> PrecondDynamicScanLoopResult

jax.lax.scan variant of :func:ekf1_sqr_loop_preconditioned_dynamic.

Operates in preconditioned (bar) space. Calibration is invariant under preconditioning: H Q H.T = (HT) Q_bar (HT).T, so all four modes behave identically to the original-space variants.

Parameters:

Name Type Description Default
mu_0 Array

Initial state mean (original space).

required
P_0_sqr Array

Initial state covariance (square-root form, original space).

required
prior BasePrior

Preconditioned prior (e.g. PrecondIWP).

required
measure BaseODEInformation

Measurement model.

required
tspan tuple[float, float]

Time interval (t_start, t_end).

required
N int

Number of filter steps.

required
calibration str

"dynamic" (default), "diagonal", "diagonal_ekf0", or "none". All modes work with :class:PrecondJointPrior; diagonal modes require prior._prior_x.xi to be diagonal.

'dynamic'
min_sigma_sqr float

Lower bound applied to per-step sigma_hat^2 before baking into Q_step_sqr_bar. Default 0.0.

0.0
calibrate bool | None

Deprecated; True -> "dynamic", False -> "none".

None

Returns:

Type Description
PrecondDynamicScanLoopResult

class:PrecondDynamicScanLoopResult -- 12 arrays, T_h, and the

PrecondDynamicScanLoopResult

scalar log-marginal-likelihood. sigma_sqr_seq shape is

PrecondDynamicScanLoopResult

(N, d_ode) in diagonal modes, (N,) otherwise.

ode_filters.filters.ode_filter_step

ode_filters.filters.ode_filter_step.ekf1_sqr_filter_step(A_t: Array, b_t: Array, Q_t_sqr: Array, m_prev: Array, P_prev_sqr: Array, measure: BaseODEInformation, t: float = 0.0) -> FilterStepResult

Perform a single square-root EKF prediction and update step.

Parameters:

Name Type Description Default
A_t Array

State transition matrix for current step.

required
b_t Array

Drift vector for current step.

required
Q_t_sqr Array

Square-root of process noise covariance.

required
m_prev Array

Previous state mean estimate.

required
P_prev_sqr Array

Previous state covariance (square-root form).

required
measure BaseODEInformation

Measurement model (e.g., ODEInformation or subclass).

required
t float

Current time (default 0.0).

0.0

Returns:

Type Description
FilterStepResult

Tuple of 4 tuples containing prediction, backward pass, and update results.

ode_filters.filters.ode_filter_step.rts_sqr_smoother_step(G_back: Array, d_back: Array, P_back_sqr: Array, m_s: Array, P_s_sqr: Array) -> tuple[Array, Array]

Perform a single Rauch-Tung-Striebel backward smoothing step.

Parameters:

Name Type Description Default
G_back Array

Backward pass gain matrix.

required
d_back Array

Backward pass offset vector.

required
P_back_sqr Array

Backward pass covariance (square-root form).

required
m_s Array

Smoothed state mean from next time step.

required
P_s_sqr Array

Smoothed state covariance (square-root form) from next time step.

required

Returns:

Type Description
tuple[Array, Array]

Tuple of previous smoothed mean and covariance (square-root form).

ode_filters.filters.ode_filter_step.ekf1_sqr_filter_step_preconditioned(A_bar: Array, b_bar: Array, Q_sqr_bar: Array, T_t: Array, m_prev_bar: Array, P_prev_sqr_bar: Array, measure: BaseODEInformation, t: float = 0.0) -> PreconditionedFilterStepResult

Perform a single preconditioned square-root EKF step.

Parameters:

Name Type Description Default
A_bar Array

Stepsize-independent state transition matrix.

required
b_bar Array

Stepsize-independent drift vector.

required
Q_sqr_bar Array

Square-root of stepsize-independent process noise covariance.

required
T_t Array

Preconditioning transformation matrix for current step.

required
m_prev_bar Array

Previous state mean estimate (preconditioned space).

required
P_prev_sqr_bar Array

Previous state covariance (square-root form, preconditioned space).

required
measure BaseODEInformation

Measurement model (e.g., ODEInformation or subclass).

required
t float

Current time (default 0.0).

0.0

Returns:

Type Description
PreconditionedFilterStepResult

Tuple of 5 tuples with preconditioned and original-space results.

ode_filters.filters.ode_filter_step.rts_sqr_smoother_step_preconditioned(G_back_bar: Array, d_back_bar: Array, P_back_sqr_bar: Array, m_s_bar: Array, P_s_sqr_bar: Array, T_t: Array) -> tuple[tuple[Array, Array], tuple[Array, Array]]

Perform a single preconditioned Rauch-Tung-Striebel backward smoothing step.

Parameters:

Name Type Description Default
G_back_bar Array

Backward pass gain matrix (preconditioned space).

required
d_back_bar Array

Backward pass offset vector (preconditioned space).

required
P_back_sqr_bar Array

Backward pass covariance (square-root form, preconditioned space).

required
m_s_bar Array

Smoothed state mean from next step (preconditioned space).

required
P_s_sqr_bar Array

Smoothed state covariance (square-root form, preconditioned space).

required
T_t Array

Preconditioning transformation matrix.

required

Returns:

Type Description
tuple[Array, Array]

Tuple containing:

tuple[Array, Array]
  • (mean, covariance) in preconditioned space
tuple[tuple[Array, Array], tuple[Array, Array]]
  • (mean, covariance) in original space

ode_filters.filters.ode_filter_step.ekf1_sqr_filter_step_sequential(A_t: Array, b_t: Array, Q_t_sqr: Array, m_prev: Array, P_prev_sqr: Array, measure: BaseODEInformation, t: float = 0.0, *, obs: tuple[Array, Array, Array] | None = None) -> SeqFilterStepResult

Perform a sequential square-root EKF step: ODE update then observation.

Unlike the joint update in ekf1_sqr_filter_step, this function first incorporates ODE + Conservation information, then applies an optional observation update using the ODE-updated state as the prior.

Parameters:

Name Type Description Default
A_t Array

State transition matrix for current step.

required
b_t Array

Drift vector for current step.

required
Q_t_sqr Array

Square-root of process noise covariance.

required
m_prev Array

Previous state mean estimate.

required
P_prev_sqr Array

Previous state covariance (square-root form).

required
measure BaseODEInformation

Measurement model (e.g., ODEInformation or subclass).

required
t float

Current time (default 0.0).

0.0
obs tuple[Array, Array, Array] | None

Optional observation tuple (H_obs, c_obs, R_obs_sqr). When provided, the observation update uses the ODE-updated state (m_ode, P_ode_sqr) as its prior.

None

Returns:

Type Description
SeqFilterStepResult

Tuple of 5 tuples: prediction, backward pass, ODE observation

SeqFilterStepResult

marginal, observation marginal, and final updated state.

ode_filters.filters.ode_filter_step.ekf1_sqr_filter_step_preconditioned_sequential(A_bar: Array, b_bar: Array, Q_sqr_bar: Array, T_t: Array, m_prev_bar: Array, P_prev_sqr_bar: Array, measure: BaseODEInformation, t: float = 0.0, *, obs: tuple[Array, Array, Array] | None = None) -> SeqPrecondFilterStepResult

Perform a sequential preconditioned square-root EKF step.

Parameters:

Name Type Description Default
A_bar Array

Stepsize-independent state transition matrix.

required
b_bar Array

Stepsize-independent drift vector.

required
Q_sqr_bar Array

Square-root of stepsize-independent process noise covariance.

required
T_t Array

Preconditioning transformation matrix for current step.

required
m_prev_bar Array

Previous state mean estimate (preconditioned space).

required
P_prev_sqr_bar Array

Previous state covariance (square-root form, preconditioned).

required
measure BaseODEInformation

Measurement model (e.g., ODEInformation or subclass).

required
t float

Current time (default 0.0).

0.0
obs tuple[Array, Array, Array] | None

Optional observation tuple (H_obs, c_obs, R_obs_sqr).

None

Returns:

Type Description
SeqPrecondFilterStepResult

Tuple of 6 tuples with preconditioned and original-space results.

ode_filters.filters.ode_filter_step.ekf1_sqr_filter_step_sequential_scan(A_t: Array, b_t: Array, Q_t_sqr: Array, m_prev: Array, P_prev_sqr: Array, measure: BaseODEInformation, t: float, H_obs: Array, c_obs: Array, R_obs_sqr: Array, obs_active: Array) -> tuple[tuple[Array, Array], tuple[Array, Array, Array], tuple[Array, Array], tuple[Array, Array], tuple[Array, Array]]

Scan-compatible sequential square-root EKF step.

Always executes the observation update with fixed-shape arrays. When obs_active is False, the ODE-only state is selected via jnp.where so that observations have no effect.

Parameters:

Name Type Description Default
A_t Array

State transition matrix.

required
b_t Array

Drift vector.

required
Q_t_sqr Array

Square-root process noise covariance.

required
m_prev Array

Previous state mean.

required
P_prev_sqr Array

Previous state covariance (square-root form).

required
measure BaseODEInformation

Measurement model (ODE + Conservation only).

required
t float

Current time.

required
H_obs Array

Observation Jacobian, shape [obs_dim, state_dim].

required
c_obs Array

Observation offset, shape [obs_dim].

required
R_obs_sqr Array

Square-root observation noise, shape [obs_dim, obs_dim].

required
obs_active Array

Scalar boolean — whether observations are active.

required

Returns:

Type Description
tuple[Array, Array]

Tuple of 5 tuples: prediction, backward pass, ODE observation

tuple[Array, Array, Array]

marginal, observation marginal, and final updated state.

ode_filters.filters.ode_filter_step.ekf1_sqr_filter_step_preconditioned_sequential_scan(A_bar: Array, b_bar: Array, Q_sqr_bar: Array, T_t: Array, m_prev_bar: Array, P_prev_sqr_bar: Array, measure: BaseODEInformation, t: float, H_obs: Array, c_obs: Array, R_obs_sqr: Array, obs_active: Array) -> tuple[tuple[Array, Array], tuple[Array, Array, Array], tuple[Array, Array], tuple[Array, Array], tuple[Array, Array], tuple[Array, Array]]

Scan-compatible preconditioned sequential square-root EKF step.

Parameters:

Name Type Description Default
A_bar Array

Stepsize-independent state transition matrix.

required
b_bar Array

Stepsize-independent drift vector.

required
Q_sqr_bar Array

Square-root stepsize-independent process noise.

required
T_t Array

Preconditioning transformation matrix.

required
m_prev_bar Array

Previous state mean (preconditioned space).

required
P_prev_sqr_bar Array

Previous state covariance (preconditioned, square-root).

required
measure BaseODEInformation

Measurement model (ODE + Conservation only).

required
t float

Current time.

required
H_obs Array

Observation Jacobian, shape [obs_dim, state_dim].

required
c_obs Array

Observation offset, shape [obs_dim].

required
R_obs_sqr Array

Square-root observation noise, shape [obs_dim, obs_dim].

required
obs_active Array

Scalar boolean — whether observations are active.

required

Returns:

Type Description
tuple[tuple[Array, Array], tuple[Array, Array, Array], tuple[Array, Array], tuple[Array, Array], tuple[Array, Array], tuple[Array, Array]]

Tuple of 6 tuples with preconditioned and original-space results.