Unscented Kalman Inversion

EnsembleKalmanProcesses.UnscentedType

An Unscented Kalman Inversion process.

mutable struct Unscented{FT<:AbstractFloat, IT<:Int64} <: EnsembleKalmanProcesses.Process

Fields

  • u_mean::Any: iterable of vectors of length N_parameters containing the parameter mean at each UKI iteration; taken prior to the prediction step and therefore not equal to the sigma-ensemble mean

  • uu_cov::Any: iterable of matrices of size (N_parameters, N_parameters) containing the parameter covariance at each UKI iteration; taken prior to the prediction step and therefore not equal to the sigma-ensemble covariance

  • obs_pred::Any: iterable of vectors of length N_y containing the predicted observation mean at each UKI iteration

  • c_weights::AbstractVecOrMat{FT} where FT<:AbstractFloat: sigma-point weights used to shift the mean; vector of length N_ens for symmetric sigma points or matrix of size (N_parameters, N_ens) for simplex sigma points

  • mean_weights::AbstractVector{FT} where FT<:AbstractFloat: quadrature weights used to reconstruct the mean from the sigma ensemble

  • cov_weights::AbstractVector{FT} where FT<:AbstractFloat: quadrature weights used to reconstruct the covariance from the sigma ensemble

  • N_ens::Int64: number of sigma particles: 2N_parameters + 1 for symmetric or N_parameters + 2 for simplex

  • Σ_ω::AbstractMatrix{FT} where FT<:AbstractFloat: covariance of the artificial evolution noise added during the prediction step

  • Σ_ν_scale::AbstractFloat: scaling factor for the artificial observation noise covariance

  • α_reg::AbstractFloat: regularization parameter controlling shrinkage toward the prior mean (0 < α_reg ≤ 1)

  • r::AbstractVector{FT} where FT<:AbstractFloat: regularization reference vector; defaults to the prior mean

  • update_freq::Int64: frequency at which the evolution covariance Σ_ω is updated; 0 disables updates

  • impose_prior::Bool: flag to use augmented-system Tikhonov regularization (Chada et al. 2020, Huang et al. 2022), which imposes the prior during inversion

  • prior_mean::Any: prior mean used for regularization; defaults to the initial mean

  • prior_cov::Any: prior covariance used for regularization; defaults to the initial covariance

  • iter::Int64: current iteration number

Constructors

Unscented(
    u_mean,
    uu_cov,
    obs_pred,
    c_weights,
    mean_weights,
    cov_weights,
    N_ens,
    Σ_ω,
    Σ_ν_scale,
    α_reg,
    r,
    update_freq,
    impose_prior,
    prior_mean,
    prior_cov,
    iter
)

defined at /home/runner/work/EnsembleKalmanProcesses.jl/EnsembleKalmanProcesses.jl/src/UnscentedKalmanInversion.jl:17.

Unscented(
    u0_mean,
    uu0_cov;
    α_reg,
    update_freq,
    modified_unscented_transform,
    impose_prior,
    prior_mean,
    prior_cov,
    sigma_points
)

defined at /home/runner/work/EnsembleKalmanProcesses.jl/EnsembleKalmanProcesses.jl/src/UnscentedKalmanInversion.jl:75.

Unscented(prior; kwargs...)

defined at /home/runner/work/EnsembleKalmanProcesses.jl/EnsembleKalmanProcesses.jl/src/UnscentedKalmanInversion.jl:200.

source
EnsembleKalmanProcesses.TransformUnscentedType

A square-root-transform variant of the Unscented Kalman Inversion (UTKI) process.

TransformUnscented replaces the classical UKI analysis step with a computationally efficient square-root/Woodbury formulation, reducing the cost of the covariance update from O(Npar³) to O(Nens · N_par²). It shares all constructors with Unscented; see that type for the full constructor reference.

mutable struct TransformUnscented{FT<:AbstractFloat, IT<:Int64} <: EnsembleKalmanProcesses.Process

Fields

  • u_mean::Any: iterable of vectors of length N_parameters containing the parameter mean at each UKI iteration; taken prior to the prediction step and therefore not equal to the sigma-ensemble mean

  • uu_cov::Any: iterable of matrices of size (N_parameters, N_parameters) containing the parameter covariance at each UKI iteration; taken prior to the prediction step and therefore not equal to the sigma-ensemble covariance

  • obs_pred::Any: iterable of vectors of length N_y containing the predicted observation mean at each UKI iteration

  • c_weights::AbstractVecOrMat{FT} where FT<:AbstractFloat: sigma-point weights used to shift the mean; vector of length N_ens for symmetric sigma points or matrix of size (N_parameters, N_ens) for simplex sigma points

  • mean_weights::AbstractVector{FT} where FT<:AbstractFloat: quadrature weights used to reconstruct the mean from the sigma ensemble

  • cov_weights::AbstractVector{FT} where FT<:AbstractFloat: quadrature weights used to reconstruct the covariance from the sigma ensemble

  • N_ens::Int64: number of sigma particles: 2N_parameters + 1 for symmetric or N_parameters + 2 for simplex

  • Σ_ω::AbstractMatrix{FT} where FT<:AbstractFloat: covariance of the artificial evolution noise added during the prediction step

  • Σ_ν_scale::AbstractFloat: scaling factor for the artificial observation noise covariance

  • α_reg::AbstractFloat: regularization parameter controlling shrinkage toward the prior mean (0 < α_reg ≤ 1)

  • r::AbstractVector{FT} where FT<:AbstractFloat: regularization reference vector; defaults to the prior mean

  • update_freq::Int64: frequency at which the evolution covariance Σ_ω is updated; 0 disables updates

  • impose_prior::Bool: flag to use augmented-system Tikhonov regularization (Chada et al. 2020, Huang et al. 2022), which imposes the prior during inversion

  • prior_mean::Any: prior mean used for regularization; defaults to the initial mean

  • prior_cov::Any: prior covariance used for regularization; defaults to the initial covariance

  • iter::Int64: current iteration number

  • buffer::AbstractVector: used to store matrices: buffer[1] = Y' * Γinv, buffer[2] = Y' * Γinv * Y

Constructors

TransformUnscented(process)

defined at /home/runner/work/EnsembleKalmanProcesses.jl/EnsembleKalmanProcesses.jl/src/UnscentedKalmanInversion.jl:217.

TransformUnscented(u0_mean, uu0_cov; kwargs...)

defined at /home/runner/work/EnsembleKalmanProcesses.jl/EnsembleKalmanProcesses.jl/src/UnscentedKalmanInversion.jl:246.

TransformUnscented(prior; kwargs...)

defined at /home/runner/work/EnsembleKalmanProcesses.jl/EnsembleKalmanProcesses.jl/src/UnscentedKalmanInversion.jl:257.

TransformUnscented(
    u_mean,
    uu_cov,
    obs_pred,
    c_weights,
    mean_weights,
    cov_weights,
    N_ens,
    Σ_ω,
    Σ_ν_scale,
    α_reg,
    r,
    update_freq,
    impose_prior,
    prior_mean,
    prior_cov,
    iter,
    buffer
)

defined at /home/runner/work/EnsembleKalmanProcesses.jl/EnsembleKalmanProcesses.jl/src/UnscentedTransformKalmanInversion.jl:23.

source
EnsembleKalmanProcesses.Gaussian_2dFunction
Gaussian_2d(
    u_mean::AbstractArray{FT<:AbstractFloat, 1},
    uu_cov::AbstractArray{FT<:AbstractFloat, 2},
    Nx::Int64,
    Ny::Int64;
    xx,
    yy
) -> Tuple{Vector, Vector, Matrix{T} where T<:AbstractFloat}

Evaluate a 2-dimensional Gaussian density on a regular grid and return the grid arrays and density values.

Arguments

  • u_mean: mean of the 2-d Gaussian (length-2 vector).
  • uu_cov: 2×2 covariance matrix.
  • Nx: number of grid points along the first dimension.
  • Ny: number of grid points along the second dimension.
  • xx: optional custom grid along the first dimension; computed from u_mean and uu_cov if nothing.
  • yy: optional custom grid along the second dimension; computed from u_mean and uu_cov if nothing.

Returns (xx, yy, Z) where xx and yy are the 1-d grid vectors and Z is the Nx × Ny density matrix.

source
EnsembleKalmanProcesses.construct_sigma_ensembleFunction
construct_sigma_ensemble(
    process::Union{TransformUnscented, Unscented},
    x_mean::AbstractArray{FT<:AbstractFloat, 1},
    x_cov::AbstractArray{FT<:AbstractFloat, 2}
) -> Any

Return the sigma-point ensemble matrix of size (N_x, N_ens) built from mean x_mean and covariance x_cov.

The layout and number of sigma points depend on the sigma_points scheme stored in process ("symmetric" or "simplex"). When Cholesky factorization of x_cov fails, an SVD-based fallback is used.

source
EnsembleKalmanProcesses.construct_meanFunction
construct_mean(
    uki::EnsembleKalmanProcess{FT<:AbstractFloat, IT<:Int64, UorTU<:Union{TransformUnscented, Unscented}},
    x::Union{AbstractArray{FT<:AbstractFloat, 1}, AbstractArray{FT<:AbstractFloat, 2}};
    mean_weights
) -> Any

Return the weighted mean of an ensemble x using the UKI quadrature weights.

When x is a matrix, each column is treated as one ensemble member and the result is a vector. When x is a vector, the result is a scalar.

source
EnsembleKalmanProcesses.construct_covFunction
construct_cov(
    uki::EnsembleKalmanProcess{FT<:AbstractFloat, IT<:Int64, UorTU<:Union{TransformUnscented, Unscented}},
    x::Union{AbstractArray{FT<:AbstractFloat, 1}, AbstractArray{FT<:AbstractFloat, 2}};
    ...
) -> Any
construct_cov(
    uki::EnsembleKalmanProcess{FT<:AbstractFloat, IT<:Int64, UorTU<:Union{TransformUnscented, Unscented}},
    x::Union{AbstractArray{FT<:AbstractFloat, 1}, AbstractArray{FT<:AbstractFloat, 2}},
    x_mean::Union{Nothing, AbstractArray{FT<:AbstractFloat, 1}, FT<:AbstractFloat};
    cov_weights
) -> Any

Return the weighted covariance of ensemble x with respect to mean x_mean.

When x_mean is nothing it is computed via construct_mean. When x is a matrix the result is a matrix; when x is a vector the result is a scalar.

source
construct_cov(
    uki::EnsembleKalmanProcess{FT<:AbstractFloat, IT<:Int64, UorTU<:Union{TransformUnscented, Unscented}},
    x::AbstractArray{FT<:AbstractFloat, 2},
    x_mean::AbstractArray{FT<:AbstractFloat, 1},
    obs_mean::AbstractArray{FT<:AbstractFloat, 2},
    y_mean::AbstractArray{FT<:AbstractFloat, 1};
    cov_weights
) -> Any

Return the weighted cross-covariance between ensemble x (with mean x_mean) and ensemble obs_mean (with mean y_mean).

source
EnsembleKalmanProcesses.update_ensemble_prediction!Function
update_ensemble_prediction!(
    process::Union{TransformUnscented, Unscented},
    Δt::AbstractFloat,
    u_idx::Vector{Int64}
) -> Any

Perform the UKI prediction step for parameter indices u_idx and return the new sigma ensemble as a matrix of size (N_parameters, N_ens).

The predicted mean and covariance are propagated by the regularized evolution model (α_reg, r, Σ_ω). If update_freq > 0 and the current iteration is a multiple of update_freq, the evolution covariance Σ_ω is updated from the current parameter covariance.

source
EnsembleKalmanProcesses.update_ensemble_analysis!Function
update_ensemble_analysis!(
    uki::EnsembleKalmanProcess{FT<:Real, IT<:Int64, TU<:TransformUnscented},
    u_p_full::AbstractMatrix,
    g_full::AbstractMatrix,
    u_idx::Vector{Int64},
    g_idx::Vector{Int64}
)

Perform the UTKI analysis (mean and covariance update) step in-place.

Uses a Woodbury/square-root formulation: constructs the parameter perturbation matrix X and observation perturbation matrix Y from the sigma ensemble, then solves for the updated mean and covariance without forming the full N_par × N_par gain matrix. When process.impose_prior is true the system is augmented with the prior constraints before the update.

Arguments

  • uki: The EnsembleKalmanProcess being updated.
  • u_p_full: Full parameter ensemble matrix (N_par × N_ens).
  • g_full: Full predicted-observation matrix (N_obs × N_ens).
  • u_idx: Parameter indices belonging to the current UpdateGroup.
  • g_idx: Observation indices belonging to the current UpdateGroup.
source
update_ensemble_analysis!(
    uki::EnsembleKalmanProcess{FT<:AbstractFloat, IT<:Int64, U<:Unscented},
    u_p_full::AbstractArray{FT<:AbstractFloat, 2},
    g_full::AbstractArray{FT<:AbstractFloat, 2},
    u_idx::Vector{Int64},
    g_idx::Vector{Int64}
) -> Any

Perform the UKI analysis step for an Unscented process, updating the stored parameter mean, covariance, and predicted observations in-place.

g_full is the full (N_obs, N_ens) matrix of predicted observations; only the rows in g_idx and the parameter dimensions in u_idx are updated. When impose_prior=true, uses the augmented-system formulation.

source
EnsembleKalmanProcesses.construct_initial_ensembleMethod
construct_initial_ensemble(
    prior::EnsembleKalmanProcesses.ParameterDistributions.ParameterDistribution,
    process::Union{TransformUnscented, Unscented};
    constrained
) -> Any

Return the initial sigma ensemble for an Unscented or TransformUnscented process.

Parameters are returned as columns in unconstrained space by default; pass constrained=true to transform to constrained space via prior.

Note: this function exists to inspect the initial sigma ensemble without constructing an EnsembleKalmanProcess object. Do not pass the returned ensemble into EnsembleKalmanProcess — the constructor generates it internally.

source