API Reference

Types

Main Types

QuadraticKalman.QKDataType
QKData{T<:Real,N}

A data structure for holding an M x T_bar array of real values (Y) plus two integer fields (M, T_bar) derived from its dimensions. M is the dimension of the measurement vector, 'Yt', and `Tbar` is the number of time periods minus 1 (due to the autoregressive structure of the model).

Fields

  • Y::AbstractArray{T, N} : The underlying data array (N-dimensional, element type <: Real).
  • M::Int : The first dimension of Y if Y is at least 2D. If N == 1, we define M = 1.
  • T_bar::Int : One less than the “second dimension” in a 2D case, or length(Y) - 1 for a 1D vector.
source
QuadraticKalman.QKModelType
QKModel{T<:Real, T2<:Real}

Main structure containing all parameters and moments needed for the quadratic Kalman filter.

This composite model encapsulates every component necessary to specify a quadratic state-space model. It bundles together the standard state evolution parameters, measurement parameters with quadratic terms, augmented state parameters for richer dynamics, and the unconditional moments that summarize the state behavior over time. This design enables integrated filtering and smoothing procedures within a unified framework.

Fields

  • state::StateParams{T}: Parameters governing the state process, defined by the equation Xt = μ + Φ X{t-1} + Ω εt, where μ is the initial state mean, Φ is the state transition matrix, and Ω scales the process noise εt.
  • meas::MeasParams{T}: Parameters for the measurement model, given by Yt = A + B Xt + α Y{t-1} + ∑{i=1}^M (Xt' Ci Xt) + D εt, including the intercept A, linear loading B, autoregressive term α, quadratic terms involving matrices C_i, and the measurement noise scaling D.
  • aug_state::AugStateParams{T,T2}: Parameters for the augmented state space which extend the state representation. This component includes transformed state means, transitions, and additional matrices that capture nonlinear or auxiliary features of the state process. The use of a secondary numeric type T2 facilitates compatibility with automatic differentiation.
  • moments::Moments{T}: The unconditional (or stationary) moments of both the state and the augmented state. This includes the long-run mean and covariance for the state dynamics as well as the augmented state, which are critical for initialization and diagnostic evaluation of the model.

Type Parameters

  • T: The primary numeric type used for most parameters (e.g., Float64). It must be a subtype of Real, ensuring both numerical precision and compatibility with standard arithmetic operations.
  • T2: A potentially different numeric type used specifically for parameters like Lambda in AugStateParams, often employed to leverage automatic differentiation (AD) techniques.
source
QuadraticKalman.StateParamsType
StateParams{T<:Real}

A structure representing the parameters that govern the evolution of the state in a state-space model. This type encapsulates the fundamental components required to describe the behavior of the state process, including its initial mean, transition dynamics, noise characteristics, and the resulting covariance.

Fields

  • N::Int The dimensionality of the state vector. This parameter specifies the number of state variables.
  • mu::AbstractVector{T} The initial state mean vector. It must have a length equal to N.
  • Phi::AbstractMatrix{T} The state transition matrix, which models how the state evolves over time. This matrix should be of size N×N.
  • Omega::AbstractMatrix{T} The state noise matrix, used to scale the impact of the stochastic noise on the state evolution. It must be of size N×N.
  • Sigma::AbstractMatrix{T} The state covariance matrix, typically computed as Omega * Omega'. This matrix quantifies the uncertainty in the state.

Details

The state evolution of a typical state-space model is represented by the equation: Xₜ = mu + Phi * Xₜ₋₁ + Omega * εₜ, where εₜ represents a white noise process. This structure is a critical component in facilitating both the filtering and smoothing processes within the QuadraticKalman framework, ensuring that the model's dynamics are accurately captured and that its stability conditions can be properly validated.

This structure is also designed to integrate smoothly with automatic differentiation tools, taking advantage of Julia's type system to provide both precision and performance in numerical computations.

source
QuadraticKalman.MeasParamsType
MeasParams{T<:Real}

A container for the measurement equation parameters in a quadratic state-space model.

This structure holds all measurement-related parameters essential for specifying the measurement equation:

Yₜ = A + B * Xₜ + α * Yₜ₋₁ + ∑₍ᵢ₌₁₎ᴹ (Xₜ' * Cᵢ * Xₜ) + noise

where: • M: The number of measurement variables. • A: A vector of intercept terms (length M). • B: An M×N matrix mapping the state vector (of dimension N) to the measurement space. • C: A vector of M matrices, each of size N×N, representing quadratic measurement parameters. • D: An M×M matrix scaling the measurement noise. • α: An M×M matrix capturing the autoregressive component in the measurement equation. • V: An auxiliary M×M matrix computed as V = D * D', representing the covariance structure of the measurement noise.

All fields should be provided as concrete matrices or vectors (or will be derived from UniformScaling objects as needed), ensuring consistency and compatibility with downstream filtering and smoothing routines. The use of the @with_kw macro facilitates clear initialization and automatic field assignment.

source
QuadraticKalman.MomentsType
Moments{T<:Real}

Unconditional Moments Structure for the Quadratic Kalman Filter.

This structure encapsulates the long-run (stationary) moments of both the state and the augmented state. These moments include the mean and covariance estimates, which are critical for initializing the filter and for conducting diagnostic evaluations of the model dynamics.

Fields:

  • state_mean::AbstractVector{T}: Unconditional (stationary) mean vector of the state.
  • state_cov::AbstractMatrix{T}: Unconditional covariance matrix of the state.
  • aug_mean::AbstractVector{T}: Unconditional mean vector of the augmented state.
  • aug_cov::AbstractMatrix{T}: Unconditional covariance matrix of the augmented state.
source
QuadraticKalman.AugStateParamsType
AugStateParams{T<:Real, T2<:Real}

An augmented state parameter container designed for quadratic measurement models. This type extends the conventional state-space representation to incorporate quadratic measurement features, enabling advanced filtering and smoothing algorithms to effectively handle non-linear measurement equations.

Fields:

  • mu_aug::AbstractVector{T}: The augmented state mean vector, which integrates the original state mean with additional terms arising from quadratic components.
  • Phi_aug::AbstractMatrix{T}: The augmented state transition matrix. It extends the traditional state transition dynamics to include quadratic interactions.
  • B_aug::AbstractMatrix{T}: The augmented measurement matrix that relates the extended state vector to the observed measurements, accounting for both linear and quadratic effects.
  • H_aug::AbstractMatrix{T}: The augmented selection matrix used in mapping the original state space to the augmented space, facilitating the extraction of relevant subcomponents.
  • G_aug::AbstractMatrix{T}: The augmented duplication matrix, which assists in preserving the symmetry properties of quadratic forms when processing covariance or moment adjustments.
  • Lambda::AbstractMatrix{T2}: A core structural matrix that captures the key quadratic interactions within the model. Its specific formulation supports the reconstruction of quadratic measures.
  • L1::AbstractMatrix{T}: An auxiliary matrix used for computing first-order moment corrections in the augmented state representation.
  • L2::AbstractMatrix{T}: An auxiliary matrix involved in the computation of second-order moment adjustments, contributing to the accurate determination of state covariances.
  • L3::AbstractMatrix{T}: An auxiliary matrix designed to support higher-order moment computations, often necessary for fine-tuning the filtering process.
  • P::Int: The total augmented state dimension, typically defined as the sum of the original state dimension and the square of the state dimension.

This structure is pivotal for models that incorporate quadratic measurement equations, allowing for the direct integration of quadratic terms into the state estimation process. It facilitates the derivation of augmented transition and measurement matrices, which are essential for achieving improved filtering and smoothing performance in non-linear state-space models.

source
QuadraticKalman.QKFOutputType
QKFOutput{T<:Real}

Combined container for both filtering and smoothing outputs from the Quadratic Kalman algorithm.

This type encapsulates the results of the forward filtering pass—where state estimates and associated covariances are computed recursively based solely on past and present observations—and the subsequent backward smoothing pass that refines these estimates by incorporating future observations. The unified structure provides a clear and convenient interface for diagnostic analysis, visualization, and further model-based inference tasks.

Fields

  • filter::FilterOutput{T}: Contains the outputs of the filtering process, such as the filtered log-likelihood, state estimates, and covariance matrices at each time step. These results represent the model’s estimates obtained in real time as the data was observed.
  • smoother::SmootherOutput{T}: Contains the outputs of the smoothing process, which refines and improves upon the filter results by leveraging information from the entire observation sequence. This typically includes the smoothed state estimates and corresponding covariance matrices, providing a more accurate reconstruction of the underlying state dynamics.

Details

Using the QKFOutput structure, users can conveniently access both the instantaneous (filtering) and retrospectively improved (smoothing) estimates, making it easier to perform post-hoc analysis, diagnostics, or forecasting.

source
QuadraticKalman.FilterOutputType
FilterOutput{T<:Real}

A container for the outputs produced by the Quadratic Kalman Filter applied to state-space models with quadratic measurement equations. This structure organizes and stores all key filtering results, facilitating subsequent analysis, diagnostics, or visualization.

Fields: • llt::Vector{T} A vector containing the log-likelihood values computed at each time step. • Ztt::Matrix{T} A matrix representing the filtered state estimates at the current time step. • Ptt::Array{T,3} A 3-dimensional array containing the error covariance matrices corresponding to the filtered state estimates. • Yttm1::Union{Vector{T}, Matrix{T}} The one-step-ahead (t-1) predicted measurements; it can be a vector for univariate cases or a matrix for multivariate cases. • Mttm1::Array{T,3} A 3-dimensional array holding auxiliary statistics from the filtering process. • Kt::Array{T,3} A 3-dimensional array of Kalman gain matrices computed at each filter update. • Zttm1::Matrix{T} A matrix of the one-step-ahead state predictions (prior estimates). • Pttm1::Array{T,3} A 3-dimensional array representing the error covariance matrices for the one-step-ahead state predictions.

Usage: This structure is used to encapsulate all relevant outputs from the Quadratic Kalman Filter, ensuring that users can easily access and work with the filtered estimates, prediction errors, and associated metrics that describe the performance of the filtering process.

source
QuadraticKalman.SmootherOutputType
SmootherOutput{T<:Real}

Container for outputs from the Quadratic Kalman Smoother. This structure encapsulates the results produced by the smoothing algorithm, which refines state estimates by incorporating information from both past and future observations. The smoother typically yields more accurate state estimates and associated uncertainty quantification than the filter alone.

Fields

  • Z_smooth::Matrix{T}: A matrix containing the smoothed state estimates. The matrix dimensions are P × (T̄+1), where P represents the dimension of the state vector and T̄+1 denotes the number of time steps.
  • P_smooth::Array{T,3}: A three-dimensional array holding the smoothed covariance matrices. Each slice P_smooth[:, :, t] corresponds to the covariance estimate for the state at time step t, with overall dimensions P × P × (T̄+1).

Details

The smoothed states and covariances are calculated using the Quadratic Kalman Smoother, which enhances the filtering results by backward recursion. This allows for improved state estimation by considering future as well as past measurements. The structure is essential for diagnostic checks and subsequent analyses in applications where state estimation uncertainty plays a critical role.

source

Plotting Types

QuadraticKalman.KalmanFilterTruthPlotType
KalmanFilterTruthPlot{T,M<:FilterOutput{<:Real}}

A type for plotting Kalman filter results against true states.

Fields

  • true_states::T: The true state values to compare against
  • results::M: The output from running the Kalman filter (FilterOutput)

This type is used internally by the plotting recipes to create comparison plots between the true states and the Kalman filter estimates, including confidence intervals.

source
QuadraticKalman.KalmanSmootherTruthPlotType
KalmanSmootherTruthPlot{T,M<:SmootherOutput{<:Real}}

A type for plotting Kalman smoother results against true states.

Fields

  • true_states::T: An N×T_bar matrix representing the actual underlying state values, X_t, that the model is attempting to filter. Each column corresponds to the state at a specific time step.
  • results::M: The output from running the Kalman smoother, provided as a SmootherOutput, which includes the smoothed state estimates and associated covariance matrices.

This type is used internally by the plotting recipes to create comparison plots between the true states and the Kalman smoother estimates. The plots typically display:

  • The true state trajectories (as represented by the N×T_bar matrix of X_t values).
  • The smoother’s state estimates along with confidence intervals,

thus facilitating a visual diagnostic of the smoother's performance in recovering the true state dynamics.

source
QuadraticKalman.KalmanFilterPlotType
KalmanFilterPlot{M<:FilterOutput{<:Real}}

An abstraction for constructing detailed visualizations of Kalman filter outcomes. This type is designed to encapsulate all necessary components for plotting the evolution of state estimates, along with the corresponding uncertainty measures, over time. It serves as a central piece within the plotting recipes framework, enabling users to create diagnostic graphics that compare filtered state trajectories against observed data, thereby providing deep insights into the performance and reliability of the Kalman filter.

Fields

  • results::M: The filtering results computed from a Kalman filter run, expected to be an instance of FilterOutput. This container typically includes:
    • Filtered state estimates (Z_tt): The real-time estimates of the states.
    • Covariance matrices (P_tt): The corresponding error covariances associated with the state estimates.
    • Log-likelihood values (ll_t): Metrics that indicate the goodness-of-fit at each time step.
    • Supplementary outputs such as one-step-ahead predictions and other auxiliary statistics necessary for an in-depth analysis of the filtering process.

Details

This type is intended for internal use by the plotting recipes module to standardize and streamline the visualization process. Its role is to facilitate the generation of clear and informative plots that not only display the state estimates but also visually represent the associated uncertainty through confidence intervals. This makes it easier to assess the dynamic behavior of the filter, detect potential issues, and diagnose problems related to model fit and convergence.

The comprehensive documentation provided here ensures that users extending or interfacing with the plotting system will have a complete understanding of the type's purpose and structure, thereby enhancing the overall diagnostic workflow.

source
QuadraticKalman.KalmanSmootherPlotType
KalmanSmootherPlot{M<:SmootherOutput{<:Real}}

A container for generating visualizations of Kalman smoother outputs. This type encapsulates the results of the smoother, which typically include smoothed state estimates and the corresponding covariance matrices, providing a unified interface for plotting these outputs along with their confidence intervals.

Fields

  • results::M: An instance of SmootherOutput that contains:
    • Z_smooth: A matrix of the smoothed state estimates.
    • P_smooth: A three-dimensional array of the smoothed state covariances.
    These outputs are generated by executing a backward recursion that refines state estimates using future information.

Usage

This type is used internally by the plotting recipes to produce diagnostic plots that display the smoother’s performance. The standard visualizations enable users to visually compare the smoothed state trajectories against true state values (or other benchmarks), with uncertainty represented as confidence intervals. This facilitates a detailed analysis of the model’s estimation efficiency and convergence properties.

Details

By encapsulating the smoother results, KalmanSmootherPlot simplifies the integration of smoothing diagnostics into the plotting framework. The associated plots typically illustrate:

  • The estimated state trajectories over time.
  • The uncertainty in these estimates, often depicted using shaded regions corresponding to a 95% confidence interval.

This enhanced visualization aids in assessing the accuracy and reliability of the smoothing procedure.

source

Functions

Filtering and Smoothing

QuadraticKalman.qkf_filterFunction
qkf_filter(data::QKData{T1,1}, model::QKModel{T,T2})

Run the Quadratic Kalman Filter (QKF) on a time series of length , returning a new set of result arrays (out-of-place).

Description

This function implements the same quadratic Kalman filter recursion as qkf_filter!, but instead of updating arrays in-place, it allocates new arrays for predictions, updates, and outputs. This can be simpler to use in contexts where you don't want to mutate or reuse data and model, but it may be less memory-efficient for large-scale problems.

At each time step, it performs:

  1. State Prediction (predict_Z_tt / predict_P_tt)
  2. Measurement Prediction (predict_Y_tt / predict_M_ttm1)
  3. Kalman Gain Computation (compute_K_t)
  4. State & Covariance Update (update_Z_tt, update_P_tt)
  5. PSD Correction (correct_Zₜₜ)
  6. Log-Likelihood computation.

Arguments

  • data::QKData{T1,1} Same structure as in qkf_filter!, with fields:
    • Y::Vector{T1}, T̄::Int, M::Int.
  • model::QKModel{T,T2} Same parameter structure as in qkf_filter!, with fields:
    • N::Int, P::Int, μ̃ᵘ, Σ̃ᵘ, etc.

Return

A named tuple with fields:

  • ll_t::Vector{Float64} Per-time-step log-likelihoods (size = ).
  • Z_tt::Array{T,3} The filtered state at each time step (dimensions (T, P, T̄+1) in your usage).
  • P_tt::Array{T,4} The filtered state covariance array.
  • Y_ttm1::Vector{T} Predicted measurement for each step.
  • M_ttm1::Array{T,4} Predicted measurement covariances.
  • K_t::Array{T,4} The Kalman gain for each time step.
  • Z_ttm1::Array{T,3} One-step-ahead predicted states.
  • P_ttm1::Array{T,4} One-step-ahead predicted covariances.

Details

  1. Initialization:
    • Creates new arrays for Z_tt and P_tt and sets the initial state to aug_mean and aug_cov.
  2. Time Loop:
    • Prediction: predict_Z_tt, predict_P_tt.
    • Measurement: predict_Y_tt, predict_M_ttm1.
    • Gain & Update: compute_K_t, update_Z_tt, update_P_tt.
    • Correction: correct_Z_tt for PSD.
    • Likelihood: compute_loglik.
  3. No In-Place Mutation:
    • Each step returns fresh arrays; original inputs are not modified.

Example

data = QKData(Y, M=measurement_dim, T̄=length(Y)-1)
model = QKModel(...)
result = qkf_filter(data, model)

@show result.ll_t
@show result.Z_tt[:, end]   # final state
source
QuadraticKalman.qkf_filter!Function
qkf_filter!(data::QKData{T1,1}, model::QKModel{T,T2})

Run the Quadratic Kalman Filter (QKF) on a time series of length , modifying the result in-place.

Description

This function implements a Kalman-like recursive filter where the state vector Zₜ includes not only the usual mean component xₜ but also terms for the second-moment (x xᵀ)ₜ, making it a quadratic extension. At each time step, it performs:

  1. State Prediction (predict_Z_ttm1! / predict_P_ttm1!)
  2. Measurement Prediction (predict_Y_ttm1! / predict_M_ttm1!)
  3. Kalman Gain Computation (compute_K_t!)
  4. State & Covariance Update (update_Z_tt!, update_P_tt!)
  5. PSD Correction (correct_Z_tt!)
  6. Log-Likelihood computation for the current innovation.

Unlike the non-mutating version (qkf_filter), this function reuses and mutates internal arrays and data structures in-place, which can improve performance and reduce memory allocations.

Arguments

  • data::QKData{T1,1} A structure holding:

    • Y::Vector{T1} of length T̄+1, containing observations. Typically, Y[1] is an initial placeholder and Y[2..end] are the actual measurements.
    • T_bar::Int the total number of time steps (excluding index 0).
    • M::Int the dimension of the measurement at each time step.
  • model::QKModel{T,T2} A parameter structure holding:

    • N::Int: State dimension (for the mean part).
    • P::Int: Dimension of the augmented "quadratic" state (P = N + N(N+1)/2).
    • μ_aug, Σ_aug: The unconditional mean and covariance used for initialization.
    • Additional model matrices or functions (e.g., Φ_aug, B_aug, A, V) accessed via subroutines.

Return

A named tuple with fields:

  • ll_t::Vector{Float64} The per-time-step log-likelihoods (size = ).
  • Z_tt::Array{T,3} The filtered state at each time step. Dimensions: (T, P, T̄+1) in your specific code (or (P, T̄+1) in a more generic version).
  • Pₜₜ::Array{T,4} The filtered state covariance array. Dimensions often (T, P, P, T̄+1) in your code.
  • Yₜₜ₋₁::Vector{T} The predicted measurement at each time step (size = ).
  • M_ttm1::Array{T,4} The predicted measurement covariance, dimensions (T, M, M, T̄).
  • K_t::Array{T,4} The Kalman gain for each time step, (T, P, M, T̄).
  • Zₜₜ₋₁::Array{T,3} One-step-ahead predicted states.
  • Pₜₜ₋₁::Array{T,4} One-step-ahead predicted covariances.
  • Σ_ttm1::Array{T,4} Any intermediate covariance terms used for prediction.

Details

  1. Initialization:
    • Z_tt[:, 1] .= μ̃ᵘ and P_tt[:,:,1] .= Σ̃ᵘ.
  2. Recursive Steps (for t in 1:T̄):
    • Prediction: predict_Z_ttm1! / predict_P_ttm1!.
    • Measurement: predict_Y_ttm1! / predict_M_ttm1!.
    • Gain & Update: compute_K_t!, then update_Z_tt! / update_P_tt!.
    • Correction: correct_Z_tt! clamps negative eigenvalues for PSD.
    • Likelihood: compute_loglik! appends the log-likelihood.
  3. Positive Semidefinite Correction:
    • Negative eigenvalues introduced by approximation are set to zero.

Example

data = QKData(Y, M=measurement_dim, T̄=length(Y)-1)
model = QKModel(...)
result = qkf_filter!(data, model)

@show result.ll_t
@show result.Z_tt[:, end]   # final state
source
qkf_filter!(data::QKData{T1,1}, model::QKModel{T,T2})

Run the Quadratic Kalman Filter (QKF) on a time series of length .

Description

This function implements a Kalman-like recursive filter where the state vector Z_t includes not only the usual mean component xₜ but also terms for the second-moment (x xᵀ)ₜ, making it a quadratic extension. At each time step, it performs:

  1. State Prediction (predict_Z_ttm1 / predict_P_ttm1)
  2. Measurement Prediction (predict_Y_ttm1 / predict_M_ttm1)
  3. Kalman Gain Computation (compute_K_t)
  4. State & Covariance Update (update_Zₜₜ!, update_Pₜₜ!)
  5. PSD Correction: Ensures the implied covariance is positive semidefinite by clamping negative eigenvalues via correct_Zₜₜ!.
  6. Log-Likelihood computation for the current innovation.

The filter stores and returns the entire history of filtered states, covariances, predicted measurements, and related arrays.

Arguments

  • data::QKData{T1,1} A structure holding:

    • Y::Vector{T1} of length T_bar+1, which contains observations (Y[1] is unused or some initial placeholder, and Y[2..end] are the actual measurements).
    • T_bar::Int the total number of time steps (excluding index 0).
    • M::Int the dimension of the measurement at each time step.
  • model::QKModel{T,T2} A parameter structure holding:

    • N::Int: State dimension (for the mean part).
    • P::Int: Dimension of the augmented "quadratic" state vector (P = N + N(N+1)/2).
    • aug_mean, aug_cov: The unconditional mean and covariance used for initialization.
    • Additional model matrices or functions (e.g., Phi_aug, B_aug, A, V), typically accessed via separate predict/update subroutines.

Return

A named tuple with fields:

  • ll_t::Vector{Float64} The per-time-step log-likelihoods of the innovations (size = T_bar).

  • Z_tt::Matrix{T} The updated ("filtered") state at each time step. Dimensions: (P, T_bar+1), where column k corresponds to time index k-1.

  • P_tt::Array{T,3} The updated ("filtered") state covariance array (or the augmented second-moment representation) at each time step. Dimensions: (P, P, T_bar+1.

  • Y_ttm1::Vector{T} The predicted measurement at each time step (size = T_bar).

  • M_ttm1::Array{T,3} The predicted measurement covariance at each time step (dimensions: (M, M, T_bar)).

  • K_t::Array{T,3} The Kalman gain at each time step (P, M, T_bar).

  • Z_ttm1::Matrix{T} The one-step-ahead predicted state (dimensions: (P, T_bar)).

  • P_ttm1::Array{T,4} The one-step-ahead predicted covariance (dimensions: (P, P, T_bar)).

  • Σ_ttm1::Array{T,4} Any intermediate state-dependent covariance terms added during prediction (dimensions: (P, P, T_bar)).

Details

  1. Initialization:

    • Z_tt[:,1] .= aug_mean and P_tt[:,:,1] .= aug_cov, representing the state at time 0.
  2. Time Loop (for t in 1:T_bar):

    • Predict step:
      • predict_Z_ttm1 sets Z_ttm1[:,t] from Z_tt[:,t].
      • predict_P_ttm1 sets P_ttm1[:,:,t] from P_tt[:,:,t].
    • Measurement step:
      • predict_Y_ttm1 computes Y_ttm1[t] from Z_ttm1[:,t].
      • predict_M_ttm1 updates M_ttm1[:,:,t].
    • Gain & Update:
      • compute_K_t obtains the Kalman gain K_t[:,:,t].
      • update_Z_tt and update_P_tt yield the next Z_tt[:,t+1] and P_tt[:,:,t+1].
    • Correction:
      • correct_Z_tt clamps negative eigenvalues in the implied covariance portion of Z_tt to ensure PSD.
    • Likelihood:
      • compute_loglik appends/logs the innovation-based log-likelihood in ll_t[t].
  3. Positive Semidefinite Correction:

    • Because the filter tracks [X Xᵀ] - X Xᵀᵀ as a covariance block, it can become indefinite due to linearization or numerical issues. We enforce PSD by zeroing out negative eigenvalues in correct_Z_tt.
    • This step is not strictly differentiable at eigenvalues crossing zero. If you need AD through the filter, consider an alternative correction or a custom adjoint.

Example

```julia data = QKData(Y, M=measurementdim, Tbar=length(Y)-1) model = QKModel( ... ) # set up your model

result = qkf_filter!(data, model)

@show result.llt @show result.Ztt[:, end] # final state vector

source
QuadraticKalman.qkf_smootherFunction
qkf_smoother( Z, P, Z_pred, P_pred, T_bar, H_aug, G_aug, Φ_aug, dof ) -> (Z_smooth, P_smooth)

Out-of-place version of the QKF smoother. Returns new arrays rather than overwriting the input ones.

Description

Identical to qkfsmoother! in logic, but it allocates new arrays Zsmooth and P_smooth for the smoothed results. This is often simpler for AD frameworks that do not allow in-place mutation of arrays.

Returns

  • Z_smooth::Matrix{T}: (P × (T_bar+1)) smoothed states
  • P_smooth::Array{T,3}: (P × P × (T_bar+1)) smoothed covariances

Example

Z_smooth, P_smooth = qkf_smoother(Z, P, Z_pred, P_pred, T_bar, Hn, Gn, H_aug, Φ_aug, n)
source
qkf_smoother(filter_output::FilterOutput{T}, model::QKModel{T,T2}) where {T<:Real, T2<:Real}

Performs out-of-place backward smoothing for the Quadratic Kalman Filter (QKF) using filtering outputs. This function refines the state estimates produced during filtering by incorporating future observations through a backward smoothing pass.

Arguments

  • filteroutput::FilterOutput{T}: A container holding the outputs from the filtering phase, including: • Ztt: A matrix of filtered augmented state estimates. • Ptt: An array of filtered state covariances. • Zttm1: A matrix containing the one-step-ahead predicted states. • P_ttm1: An array containing the one-step-ahead predicted state covariances.
  • model::QKModel{T,T2}: A model specification that provides the necessary parameters for the smoothing process. The model includes an aug_state field which is unpacked to retrieve: • Haug: The augmented measurement selection matrix. • Gaug: The augmented duplication matrix used for handling quadratic forms. • Phi_aug: The augmented state transition matrix. Additionally, the state field is used to extract the dimensionality (N) of the state.

Details

This function first creates copies of the filtered state and covariance estimates to prevent modification of the original filtering outputs. It then invokes the in-place smoother routine (qkf_smoother!) on these copies using the one-step-ahead predicted values and the model's parameters. The final smoothed results are wrapped in a SmootherOutput struct and returned as fresh copies, which is particularly important for compatibility with automatic differentiation (AD) workflows.

Returns

  • SmootherOutput: A composite structure containing: • Zsmooth: A matrix of smoothed augmented state estimates. • Psmooth: An array of smoothed state covariance matrices.
source
QuadraticKalman.qkf_smoother!Function
qkf_smoother!(
    Z::AbstractMatrix{T},      # Filtered states   (P × (T_bar+1))
    P::AbstractArray{T, 3},    # Filtered covariances (P × P × (T_bar+1))
    Z_pred::AbstractMatrix{T}, # One-step-ahead predicted states (P × T_bar)
    P_pred::AbstractArray{T,3},
    T_bar::Int,
    Hn::Matrix{T},  Gn::Matrix{T},  H_aug::Matrix{T},  Φ_aug::Matrix{T},
    dof::Int
) where {T<:Real}

Perform in-place backward smoothing for the Quadratic Kalman Filter (QKF).

Description

Given the forward-filtered estimates (Z, P) from t=1..T_bar, plus the one-step-ahead predictions (Z_pred, P_pred) and the special matrices (Haug, Gaug) that handle the dimension reduction via Vech(·)/Vec(·), this function computes Z[:,t] and P[:,:,t] for t = T_bar-1 .. 1 in backward fashion to produce the smoothed estimates (Zₜ|Tbar, PZₜ|Tbar).

Mathematical Form (Backward Pass)

  1. Compute Fₜ = (H̃ₙPᵗ|ᵗᶻH̃ₙ')(H̃ₙΦ̃G̃ₙ)'(H̃ₙPᵗ⁺¹|ᵗᶻH̃ₙ')⁻¹ but implemented via solves (rather than explicit inverses) for numerical stability.

  2. Then update (in H̃ₙ-transformed space): H̃ₙZₜ|ₜ = H̃ₙZₜ|ₜ + Fₜ(H̃ₙZₜ₊₁|ₜ - H̃ₙZₜ₊₁|ₜ)

  3. And similarly for the covariance: (H̃ₙPᵗ|ᵀᶻH̃ₙ') = (H̃ₙPᵗ|ᵗᶻH̃ₙ') + Fₜ[(H̃ₙPᵗ⁺¹|ᵀᶻH̃ₙ') - (H̃ₙPᵗ⁺¹|ᵗᶻH̃ₙ')]Fₜ'

  4. Finally, transform back to get Zₜ|T and Pᵗ|ᵀᶻ in the full (Vec/augmented) space if necessary.

Arguments

  • Z::AbstractMatrix{T}: On entry, Z[:,t] = Zₜ|T for each t. On exit, it will contain the smoothed states Zₜ|T.
  • P::AbstractArray{T,3}: On entry, P[:,:,t] = Pᵗ|ᵀᶻ. On exit, P[:,:,t] = Pᵗ|ᵀᶻ.
  • Z_pred::AbstractMatrix{T}, P_pred::AbstractArray{T,3}: The one-step-ahead predicted states/covariances from the forward pass, i.e. Z_pred[:,t] = Zₜ|ₜ, P_pred[:,:,t] = P^Zₜ|ₜ.
  • T_bar::Int: Total time steps (excluding time 0).
  • Hn::Matrix{T}, Gn::Matrix{T}: The selection/duplication operators for Vec/Vech transforms of block (x xᵀ). Usually size (n(n+1) × n(n+3)/2) or similarly.
  • H_aug::Matrix{T}, Φ_aug::Matrix{T}: The augmented versions (Haug, Gaug) used in the QKF recursion.
  • dof::Int: Dimension parameter (often n or P). Adjust to your model.

Notes

  • This function runs backward from t = T_bar-1 down to t = 1, using the final values (Z[:, T_bar], P[:,:, T_bar]) as the terminal condition (Z_{T_bar|T_bar}, P^Z_{T_bar|T_bar}).
  • If your AD library supports destructive updates, this approach should be AD-friendly; if not, consider the out-of-place version qkf_smoother.

Example

Suppose you already ran the forward filter, so you have: Z, P, Zpred, Ppred, plus your special matrices.

qkf_smoother!(Z, P, Z_pred, P_pred, T_bar, Hn, Gn, H_aug, Φ_aug, n)
source
qkf_smoother!(filter_output::FilterOutput{T}, model::QKModel{T,T2}) where {T<:Real, T2<:Real}

Performs in-place backward smoothing for a Quadratic Kalman Filter (QKF) using the outputs obtained during filtering.

This function refines the filtered state estimates and covariance matrices by incorporating future observations, thereby producing a set of smoothed estimates. It operates by first extracting the filtered states (Ztt) and their associated covariances (Ptt), as well as the one-step-ahead predictions (Zttm1 and Pttm1) from the given FilterOutput structure. It then unpacks the augmented state parameters (Haug, Gaug, Phi_aug) along with the state dimension (N) from the QKModel. The function makes local copies of the filtered estimates to avoid any modification of the original filtering results, and then calls the lower-level in-place smoothing routine to update these copies using the one-step-ahead predictions and the model parameters. The smoothed states and covariances are finally encapsulated into a SmootherOutput structure, which is returned.

Parameters:

  • filteroutput: A FilterOutput instance containing: • Ztt :: Matrix of filtered augmented state estimates for time steps 0 to T̄. • Ptt :: Array of filtered state covariance matrices. • Zttm1 :: Matrix of one-step-ahead predicted augmented states. • P_ttm1 :: Array of one-step-ahead predicted covariance matrices.
  • model: A QKModel instance providing the necessary model parameters for smoothing, including: • augstate: A structure containing: - Haug :: The augmented measurement selection matrix. - Gaug :: The augmented duplication matrix for handling quadratic forms. - Phiaug:: The augmented state transition matrix. • state: The dimension (N) of the state vector.

Returns:

  • A SmootherOutput structure containing: • Zsmooth :: Matrix of smoothed augmented state estimates. • Psmooth :: Array of smoothed state covariance matrices.
source
QuadraticKalman.qkfFunction
qkf(model::QKModel{T,T2}, data::QKData{T1,N}) where {T1<:Real, T<:Real, T2<:Real, N}

Execute the full Quadratic Kalman Filter (QKF) process by combining both the filtering and backward smoothing stages. This function first applies the filtering routine to generate real-time state estimates and then refines these estimates using the smoother to incorporate information from future observations. The resulting output is a composite object that encapsulates both the filtered and smoothed results.

Parameters:

  • model::QKModel{T,T2}: A model specification that includes the system dynamics, state transition parameters, and the augmented state representation. This object provides all necessary configurations to perform the QKF.
  • data::QKData{T1,N}: A data container comprising the observed time series measurements, formatted appropriately for the QKF. The parameter N indicates the dimension of the state vector, and T1 denotes the numerical type of the data.

Process:

  1. Filtering Stage: The function invokes qkf_filter with the provided data and model to compute the filtered state estimates and error covariances.
  2. Smoothing Stage: It then calls qkf_smoother to perform backward smoothing on the filtered results, enhancing the state estimates by leveraging future information.

Returns:

  • QKFOutput{T}: A composite output object that bundles both filtering and smoothing results. This output includes: • filteroutput: The results from the filtering phase, containing state estimates and covariances. • smootheroutput: The refined state estimates obtained after applying the smoother.

This combined output is useful for subsequent analysis, diagnostics, and visualization of the QKF performance.

source
qkf(data::QKData{T1,N}, model::QKModel{T,T2}) where {T1<:Real, T<:Real, T2<:Real, N}

This function offers a backwards-compatible interface to the Quadratic Kalman Filter (QKF) by accepting the data and model arguments in a reversed order compared to the standard interface. Typically, the preferred order is to pass the model first followed by the data. This method ensures that legacy code that follows the old argument order still functions correctly.

Parameters:

  • data::QKData{T1,N}: An instance containing the observed time series data formatted for the QKF. The data structure organizes the measurements and any associated time indices to be processed by the filter.
  • model::QKModel{T,T2}: An instance containing the model specifications which include the system dynamics, the augmented state representation, and all relevant parameters required to perform the filtering and smoothing operations.

Returns:

  • QKFOutput{T}: A composite output object that includes both the filtering results and the smoothed state estimates. This output encapsulates all intermediate steps and final results computed by invoking the standard qkf function with the proper argument order.

Note: This reversed argument order function is maintained solely for backwards compatibility. Internally, it simply calls qkf(model, data) to ensure that the original processing logic remains unchanged.

source

Convenience Functions

QuadraticKalman.get_measurementFunction
get_measurement(data::QKData{T,N}, t::Int) where {T<:Real, N}

Extract measurement at time t from QKData. For vector data returns a scalar, for matrix data returns a vector.

Arguments

  • data::QKData{T,N}: Data structure
  • t::Int: Time index (1-based)

Returns

  • For N=1: Single measurement Y[t]
  • For N=2: Vector of measurements Y[:,t]

Throws

  • ArgumentError if t is out of bounds
source
QuadraticKalman.qkf_negloglikFunction
qkf_negloglik(params::Vector{T}, data::QKData, N::Int, M::Int) where T<:Real -> Real

Compute the negative log-likelihood for a Quadratic Kalman Filter model given parameters and data.

Arguments

  • params::Vector{T}: Vector of model parameters to be converted into a QKModel
  • data::QKData: Data container with observations and optional initial conditions
  • N::Int: Dimension of the state vector
  • M::Int: Dimension of the measurement vector

Returns

The negative log-likelihood value computed by:

  1. Converting parameters to a QKModel
  2. Running the Kalman filter
  3. Taking the negative sum of the per-period log-likelihoods

Note

This function is typically used as an objective function for maximum likelihood estimation, where the goal is to minimize the negative log-likelihood.

For optimization, you may want to wrap this function with N, M and data specified:

negloglik(params) = qkf_negloglik(params, data, N, M)
source
QuadraticKalman.model_to_paramsFunction
model_to_params(model::QKModel{T, T2}) where {T<:Real, T2}

Convert a QKModel object into a vector of unconstrained parameters.

The ordering of the parameters is as follows:

  1. State parameters:

    • mu (length N)
    • Phi (N×N, stored columnwise)
    • Unconstrained parameters for the state noise scaling factor Ω: For each row i = 1,...,N and column j = 1,...,i (i.e. the lower–triangular part):
      • If i == j: the parameter is log(Ω[i,i])
      • Else: the parameter is Ω[i,j]
  2. Measurement parameters:

    • A (length M)
    • B (M×N, stored columnwise)
    • C (a vector of M matrices; each matrix is N×N and is flattened columnwise)
    • Unconstrained parameters for the measurement noise scaling factor D: For each row i = 1,...,M and column j = 1,...,i:
      • If i == j: the parameter is log(D[i,i])
      • Else: the parameter is D[i,j]
    • alpha (M×M, stored columnwise)

Returns

A vector of unconstrained parameters that, when passed to params_to_model, reconstructs the original QKModel.

source
QuadraticKalman.params_to_modelFunction
params_to_model(params::Vector{T}, N::Int, M::Int) where T<:Real -> QKModel

Convert a parameter vector into a QKModel object with state and measurement parameters.

Arguments

  • params::Vector{T}: A vector of unconstrained parameters.
  • N::Int: Dimension of the state vector.
  • M::Int: Dimension of the measurement vector.

Returns

A QKModel object containing:

  • State parameters: (mu, Phi, Omega) where the state equation is Xₜ = μ + Φ Xₜ₋₁ + Omega εₜ. Here, Omega is constructed as Omega = Dstate * Dstate′, with D_state a lower–triangular matrix (of size N×N) whose diagonal entries are positive.
  • Measurement parameters: (A, B, C, D, α) where the measurement equation is Yₜ = A + B Xₜ + α Yₜ₋₁ + ∑₍ᵢ₌₁₎ᴹ Xₜ′ Cᵢ Xₜ + D εₜ. Here, D is constructed as D = Dmeas * Dmeas′, with D_meas a lower–triangular matrix (of size M×M) whose diagonal entries are positive.
  • Augmented state parameters and model moments (computed via helper functions).

Parameter vector layout

The parameter vector is assumed to contain:

  1. State parameters:

    • First N entries: state mean mu.
    • Next N^2 entries: entries of Phi (stored columnwise).
    • Next N(N+1)/2 entries: unconstrained parameters for D_state (used to form Omega).
  2. Measurement parameters:

    • Next M entries: A.
    • Next M×N entries: entries of B (reshaped as an M×N matrix).
    • Next M×N^2 entries: entries for C. (Interpreted as M matrices of size N×N.)
    • Next M(M+1)/2 entries: unconstrained parameters for D_meas (used to form D).
    • Final M×M entries: entries of α (reshaped as an M×M matrix).

Total expected length:

N + N^2 + N(N+1)/2  +  M + M×N + M×N^2 + M(M+1)/2 + M^2
source

Plotting API

QuadraticKalman.kalman_filter_truth_plotFunction
kalman_filter_truth_plot(X, results)

Create a plot comparing true states with Kalman filter estimates.

Arguments

  • X: Matrix of true state values (N×T)
  • results: FilterOutput object containing filtered state estimates

Returns

A plot showing true states, filtered estimates, and confidence intervals

source
QuadraticKalman.kalman_smoother_truth_plotFunction
kalman_smoother_truth_plot(X, results)

Create a plot comparing true states with Kalman smoother estimates.

Arguments

  • X: Matrix of true state values (N×T)
  • results: SmootherOutput object containing smoothed state estimates

Returns

A plot showing true states, smoothed estimates, and confidence intervals

source
QuadraticKalman.kalman_filter_plotFunction
kalman_filter_plot(results)

Create a plot of Kalman filter estimates.

Arguments

  • results: FilterOutput object containing filtered state estimates

Returns

A plot showing filtered estimates and confidence intervals

source
QuadraticKalman.kalman_smoother_plotFunction
kalman_smoother_plot(results)

Create a plot of Kalman smoother estimates.

Arguments

  • results: SmootherOutput object containing smoothed state estimates

Returns

A plot showing smoothed estimates and confidence intervals

source