EnsembleKalmanProcess

Primary objects and functions

EnsembleKalmanProcesses.EnsembleKalmanProcessType
EnsembleKalmanProcess{FT <: AbstractFloat, IT <: Int, P <: Process}

Structure that is used in Ensemble Kalman processes.

Fields

  • u::Array{EnsembleKalmanProcesses.DataContainers.DataContainer{FT}} where FT<:AbstractFloat: array of stores for parameters (u), each of size [N_par × N_ens]

  • observation_series::ObservationSeries: Container for the observation(s) - and minibatching mechanism

  • N_ens::Int64: ensemble size

  • g::Array{EnsembleKalmanProcesses.DataContainers.DataContainer{FT}} where FT<:AbstractFloat: Array of stores for forward model outputs, each of size [N_obs × N_ens]

  • error_metrics::Dict: Dict of error metric

  • scheduler::EnsembleKalmanProcesses.LearningRateScheduler: Scheduler to calculate the timestep size in each EK iteration

  • accelerator::EnsembleKalmanProcesses.Accelerator: accelerator object that informs EK update steps, stores additional state variables as needed

  • Δt::Vector{FT} where FT<:AbstractFloat: stored vector of timesteps used in each EK iteration

  • update_groups::AbstractVector: vector of update groups, defining which parameters should be updated by which data

  • process::EnsembleKalmanProcesses.Process: the particular EK process (Inversion etc.)

  • rng::Random.AbstractRNG: Random number generator object (algorithm + seed) used for sampling and noise, for reproducibility. Defaults to Random.GLOBAL_RNG.

  • failure_handler::FailureHandler: struct storing failsafe update directives

  • localizer::EnsembleKalmanProcesses.Localizers.Localizer: Localization kernel, implemented for (Inversion, SparseInversion, Unscented)

  • nan_tolerance::AbstractFloat: Fraction of allowable NaNs in model output before an ensemble member is considered a failed member and handled by failure handler

  • nan_row_values::Union{Nothing, AbstractVector}: Default-value vector to impute over entire-NaN rows of data (get_obs(ekp) used if value nothing)

  • verbose::Bool: Whether to print diagnostics for each EK iteration

Generic constructor

EnsembleKalmanProcess(
    params::AbstractMatrix{FT},
    observation_series::OS,
    obs_noise_cov::Union{AbstractMatrix{FT}, UniformScaling{FT}},
    process::P;
    scheduler = DefaultScheduler(1),
    Δt = FT(1),
    rng::AbstractRNG = Random.GLOBAL_RNG,
    failure_handler_method::FM = IgnoreFailures(),
    nan_tolerance = 0.1,
    nan_row_values = nothing,
    localization_method::LM = NoLocalization(),
    verbose::Bool = false,
) where {FT <: AbstractFloat, P <: Process, FM <: FailureHandlingMethod, LM <: LocalizationMethod, OS <: ObservationSeries}

Inputs:

  • params :: Initial parameter ensemble
  • observation_series :: Container for observations (and possible minibatching)
  • process :: Algorithm used to evolve the ensemble
  • scheduler :: Adaptive timestep calculator
  • Δt :: Initial time step or learning rate
  • rng :: Random number generator
  • failure_handler_method :: Method used to handle particle failures
  • nan_tolerance :: Fraction of allowable NaNs in ensemble member before considered failure (0.1 by default)
  • nan_row_values :: Default-value vector to impute over entire-NaN rows of data (get_obs(ekp) used if value nothing)
  • localization_method :: Method used to localize sample covariances
  • verbose :: Whether to print diagnostic information

Other constructors:

EnsembleKalmanProcess(
    u,
    observation_series,
    N_ens,
    g,
    error_metrics,
    scheduler,
    accelerator,
    Δt,
    update_groups,
    process,
    rng,
    failure_handler,
    localizer,
    nan_tolerance,
    nan_row_values,
    verbose
)

defined at /home/runner/work/EnsembleKalmanProcesses.jl/EnsembleKalmanProcesses.jl/src/EnsembleKalmanProcess.jl:205.

EnsembleKalmanProcess(
    params,
    observation_series,
    process,
    configuration;
    update_groups,
    rng,
    nan_tolerance,
    nan_row_values,
    verbose
)

defined at /home/runner/work/EnsembleKalmanProcesses.jl/EnsembleKalmanProcesses.jl/src/EnsembleKalmanProcess.jl:240.

EnsembleKalmanProcess(
    params,
    observation_series,
    process;
    scheduler,
    accelerator,
    failure_handler_method,
    localization_method,
    Δt,
    update_groups,
    rng,
    nan_tolerance,
    nan_row_values,
    verbose
)

defined at /home/runner/work/EnsembleKalmanProcesses.jl/EnsembleKalmanProcesses.jl/src/EnsembleKalmanProcess.jl:338.

EnsembleKalmanProcess(params, observation, args; kwargs...)

defined at /home/runner/work/EnsembleKalmanProcesses.jl/EnsembleKalmanProcesses.jl/src/EnsembleKalmanProcess.jl:397.

EnsembleKalmanProcess(
    params,
    obs,
    obs_noise_cov,
    args;
    kwargs...
)

defined at /home/runner/work/EnsembleKalmanProcesses.jl/EnsembleKalmanProcesses.jl/src/EnsembleKalmanProcess.jl:407.

EnsembleKalmanProcess(
    observation_series,
    process;
    kwargs...
)

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

EnsembleKalmanProcess(observation, process; kwargs...)

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

EnsembleKalmanProcess(
    obs_mean,
    obs_noise_cov,
    process;
    kwargs...
)

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

source
EnsembleKalmanProcesses.construct_initial_ensembleFunction
construct_initial_ensemble(
    rng::Random.AbstractRNG,
    prior::EnsembleKalmanProcesses.ParameterDistributions.ParameterDistribution,
    N_ens::Int64;
    constrained
) -> Any

Construct the initial parameters, by sampling N_ens samples from specified prior distribution. Returned with parameters as columns in unconstrained space by default (constrain by setting constrained=true)

Note: Unscented and TransformUnscented processes require different arguments than the other processes

source
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
EnsembleKalmanProcesses.update_ensemble!Function
update_ensemble!(
    ekp::EnsembleKalmanProcess,
    g_in::AbstractMatrix;
    multiplicative_inflation,
    additive_inflation,
    additive_inflation_cov,
    s,
    Δt_new,
    ekp_kwargs...
) -> Any

Updates the ensemble according to an Inversion process. Inputs:

  • ekp :: The EnsembleKalmanProcess to update.
  • g :: Model outputs, they need to be stored as a N_obs × N_ens array (i.e data are columms).
  • multiplicative_inflation :: Flag indicating whether to use multiplicative inflation.
  • additive_inflation :: Flag indicating whether to use additive inflation.
  • additiveinflationcov :: specifying an additive inflation matrix (default is the prior covariance) assumed positive semi-definite If false (default), parameter covariance from the current iteration is used.
  • s :: Scaling factor for time step in inflation step.
  • ekpkwargs :: Keyword arguments to pass to standard ekp updateensemble!.
source

Getter functions

EnsembleKalmanProcesses.get_uFunction
get_u(ekp::EnsembleKalmanProcess, iteration::IT; return_array=true) where {IT <: Integer}

Returns the unconstrained parameters at the given iteration. Returns a DataContainer object unless return_array is true.

source
get_u(ekp::EnsembleKalmanProcess; return_array=true)

Returns the unconstrained parameters from all iterations. The outer dimension is given by the number of iterations, and the inner objects are DataContainer objects unless return_array is true.

source
EnsembleKalmanProcesses.get_gFunction
get_g(ekp::EnsembleKalmanProcess, iteration::IT; return_array=true) where {IT <: Integer}

Returns the forward model evaluations at the given iteration. Returns a DataContainer object unless return_array is true.

source
get_g(ekp::EnsembleKalmanProcess; return_array=true)

Returns the forward model evaluations from all iterations. The outer dimension is given by the number of iterations, and the inner objects are DataContainer objects unless return_array is true.

source
EnsembleKalmanProcesses.get_ϕFunction
get_ϕ(prior::ParameterDistribution, ekp::EnsembleKalmanProcess, iteration::IT; return_array=true)

Returns the constrained parameters at the given iteration.

source
get_ϕ(prior::ParameterDistribution, ekp::EnsembleKalmanProcess; return_array=true)

Returns the constrained parameters from all iterations. The outer dimension is given by the number of iterations, and the inner objects are DataContainer objects unless return_array is true.

source
EnsembleKalmanProcesses.get_obsMethod
get_obs(ekp::EnsembleKalmanProcess; build) -> Any

Get the observation from the current batch in ObservationSeries build=false: returns a vector of vectors, build=true: returns a concatenated vector,

source
EnsembleKalmanProcesses.get_obsMethod
get_obs(ekp::EnsembleKalmanProcess, iteration; build) -> Any

Get the observation for iteration from the ObservationSeries build=false: returns a vector of vectors, build=true: returns a concatenated vector,

source
EnsembleKalmanProcesses.get_obs_noise_covMethod
get_obs_noise_cov(ekp::EnsembleKalmanProcess; build) -> Any

Get the obs_noise_cov from the current batch in ObservationSeries build=false:, returns a vector of blocks, build=true: returns a block matrix,

source
EnsembleKalmanProcesses.get_obs_noise_covMethod
get_obs_noise_cov(
    ekp::EnsembleKalmanProcess,
    iteration;
    build
) -> Any

Get the obs_noise_cov for iteration from the ObservationSeries build=false:, returns a vector of blocks, build=true: returns a block matrix,

source
EnsembleKalmanProcesses.get_obs_noise_cov_invMethod
get_obs_noise_cov_inv(
    ekp::EnsembleKalmanProcess;
    build
) -> Any

Get the obs_noise_cov inverse from the current batch in ObservationSeries build=false:, returns a vector of blocks, build=true: returns a block matrix,

source
EnsembleKalmanProcesses.get_obs_noise_cov_invMethod
get_obs_noise_cov_inv(
    ekp::EnsembleKalmanProcess,
    iteration;
    build
) -> Any

Get the obsnoisecov inverse for iteration from the ObservationSeries build=false:, returns a vector of blocks, build=true: returns a block matrix,

source
EnsembleKalmanProcesses.get_u_finalFunction
get_u_final(ekp::EnsembleKalmanProcess; return_array=true)

Get the unconstrained parameters at the last iteration, returning a DataContainer Object if return_array is false.

source
EnsembleKalmanProcesses.get_g_finalFunction
get_g_final(ekp::EnsembleKalmanProcess; return_array=true)

Get forward model outputs at the last iteration, returns a DataContainer Object if return_array is false.

source
EnsembleKalmanProcesses.get_u_meanFunction
get_u_mean(ekp::EnsembleKalmanProcess, iteration::IT) where {IT <: Integer}

Returns the mean unconstrained parameter at the given iteration.

source
get_u_mean(
    uki::EnsembleKalmanProcess{FT<:AbstractFloat, IT<:Int64, UorTU<:Union{TransformUnscented, Unscented}},
    iteration::Int64
) -> Any

Return the mean unconstrained parameter vector at the given iteration.

source
EnsembleKalmanProcesses.get_u_covFunction
get_u_cov(ekp::EnsembleKalmanProcess, iteration::IT) where {IT <: Integer}

Returns the unconstrained parameter sample covariance at the given iteration.

source
get_u_cov(
    uki::EnsembleKalmanProcess{FT<:AbstractFloat, IT<:Int64, UorTU<:Union{TransformUnscented, Unscented}},
    iteration::Int64
) -> Any

Return the unconstrained parameter covariance matrix at the given iteration.

source
EnsembleKalmanProcesses.get_ϕ_meanFunction
get_ϕ_mean(prior::ParameterDistribution, ekp::EnsembleKalmanProcess, iteration::IT)

Returns the constrained transform of the mean unconstrained parameter at the given iteration.

source
EnsembleKalmanProcesses.get_u_priorFunction
get_u_prior(ekp::EnsembleKalmanProcess; return_array=true)

Get the unconstrained parameters as drawn from the prior, returning a DataContainer Object if return_array is false.

source
EnsembleKalmanProcesses.get_rngFunction
get_rng(m::FixedMinibatcher) -> Random.AbstractRNG

gets the rng field from the FixedMinibatcher object

source
get_rng(m::RandomFixedSizeMinibatcher) -> Random.AbstractRNG

gets the rng field from the RandomFixesSizeMinibatcher object

source
get_rng(ekp::EnsembleKalmanProcess)

Return rng field of EnsembleKalmanProcess.

source
EnsembleKalmanProcesses.get_cov_blocksFunction
get_cov_blocks(cov::AbstractMatrix{FT}, p::IT) where {FT <: Real, IT <: Integer}

Given a covariance matrix cov and number of parameters p, returns the matrix blocks corresponding to the u–u covariance, the u–G(u) covariance, and the G(u)–G(u) covariance.

source

Error metrics

EnsembleKalmanProcesses.compute_average_rmseFunction
compute_average_rmse(ekp::EnsembleKalmanProcess) -> Any

Computes the average covariance-weighted error over the forward model ensemble, normalized by the dimension: mean(sqrt( 1/dim(y) * tr((g_i - y)' * Γ⁻¹ * (g_i - y)))). The error is retrievable as get_error_metrics(ekp)["avg_rmse"]

source
compute_average_rmse(
    uki::EnsembleKalmanProcess{FT<:AbstractFloat, IT<:Int64, UorTU<:Union{TransformUnscented, Unscented}}
) -> Any

Return the RMSE for an Unscented process evaluated at the mean, equal to sqrt(compute_loss_at_mean(uki)).

For Unscented processes, averaging RMSE at individual sigma points is not meaningful.

source
EnsembleKalmanProcesses.compute_loss_at_meanFunction
compute_loss_at_mean(ekp::EnsembleKalmanProcess) -> Any

Computes the covariance-weighted error of the mean of the forward model output, normalized by the dimension 1/dim(y) * (ḡ - y)' * Γ⁻¹ * (ḡ - y). The error is retrievable as get_error_metrics(ekp)["loss"] or returned from get_error(ekp) if a prior is not provided to the process

source
EnsembleKalmanProcesses.compute_average_unweighted_rmseFunction
compute_average_unweighted_rmse(
    ekp::EnsembleKalmanProcess
) -> Any

Computes the average unweighted error over the forward model ensemble, normalized by the dimension: mean(sqrt( 1/dim(y) * tr((g_i - y)' * (g_i - y)))). The error is retrievable as get_error_metrics(ekp)["unweighted_avg_rmse"]

source
compute_average_unweighted_rmse(
    uki::EnsembleKalmanProcess{FT<:AbstractFloat, IT<:Int64, UorTU<:Union{TransformUnscented, Unscented}}
) -> Any

Return the unweighted RMSE for an Unscented process evaluated at the mean, equal to sqrt(compute_unweighted_loss_at_mean(uki)).

For Unscented processes, averaging RMSE at individual sigma points is not meaningful.

source
EnsembleKalmanProcesses.compute_unweighted_loss_at_meanFunction
compute_unweighted_loss_at_mean(
    ekp::EnsembleKalmanProcess
) -> Any

Computes the unweighted error of the mean of the forward model output, normalized by the dimension 1/dim(y) * (ḡ - y)' * (ḡ - y). The error is retrievable as get_error_metrics(ekp)["unweighted_loss"]

source
EnsembleKalmanProcesses.compute_bayes_loss_at_meanFunction
compute_bayes_loss_at_mean(
    ekp::EnsembleKalmanProcess
) -> Any

Computes the bayes loss of the mean of the forward model output, normalized by dimensions (1/dim(y)*dim(u)) * [(ḡ - y)' * Γ⁻¹ * (ḡ - y) + (̄u - m)' * C⁻¹ * (̄u - m)]. If the prior is not provided to the process on creation of EKP, then m and C are estimated from the initial ensemble.

The error is retrievable as get_error_metrics(ekp)["bayes_loss"] or returned from get_error(ekp) if a prior is provided to the process

source
EnsembleKalmanProcesses.compute_crpsFunction
compute_crps(ekp::EnsembleKalmanProcess) -> Any

Computes a Gaussian approximation of CRPS (continuous rank probability score) of the ensemble with the observation (performing through a whitening by C^GG, see e.g., Zheng, Sun, 2025, https://arxiv.org/abs/2410.09133).

source
EnsembleKalmanProcesses.compute_error!Function
compute_error!(ekp::EnsembleKalmanProcess) -> Any

Computes a variety of error metrics and stores this in EnsembleKalmanProcess. (retrievable with get_error_metrics(ekp)) currently available:

  • avg_rmse computed with compute_average_rmse(ekp)
  • loss computed with compute_loss_at_mean(ekp)
  • unweighed_avg_rmse computed with compute_average_unweighted_rmse(ekp)
  • unweighted_loss computed with compute_unweighted_loss_at_mean(ekp)
  • bayes_loss computed with compute_bayes_loss_at_mean(ekp)
  • crps computed with compute_crps(ekp)
source
EnsembleKalmanProcesses.get_errorFunction
get_error(ekp::EnsembleKalmanProcess)

[For back compatability] Returns the relevant loss function that is minimized over EKP iterations.

  • If prior provided to the process: get_error_metrics(ekp)["bayes_loss"], the loss computed with compute_bayes_loss_at_mean(ekp)
  • If prior not provided to the process: get_error_metrics(ekp)["loss"], the loss computed with compute_loss_at_mean(ekp)
source
EnsembleKalmanProcesses.lmul_obs_noise_cov!Function
lmul_obs_noise_cov!(
    out::AbstractMatrix,
    ekp::EnsembleKalmanProcess,
    X::AbstractVecOrMat,
    idx_set::AbstractVector
)

Convenience function to multiply X by obs_noise_cov on the left without building the full matrix, and storing the result in the first argument out, and at a set of indices idx_triple = [(block_idx, local_idx, global_idx), ...]. Here, for each triple the multiplication will:

  • select block block_idx from the unbuilt obs_noise_cov_inv
  • select local indices [:,local_idx] of this block
  • multiply with the corresponding global indices X[global_idx,:]
source
EnsembleKalmanProcesses.lmul_obs_noise_cov_inv!Function
lmul_obs_noise_cov_inv!(
    out::AbstractMatrix,
    ekp::EnsembleKalmanProcess,
    X::AbstractVecOrMat,
    idx_set::AbstractVector
)

Convenience function to multiply X by obs_noise_cov_inv on the left without building the full matrix, and storing the result in the first argument out, and at a set of indices idx_triple = [(block_idx, local_idx, global_idx), ...]. Here, for each triple the multiplication will:

  • select block block_idx from the unbuilt obs_noise_cov_inv
  • select local indices [:,local_idx] of this block
  • multiply with the corresponding global indices X[global_idx,:]

The primary use case is in update_ensemble! for TransformInversion()

source

Learning Rate Schedulers

EnsembleKalmanProcesses.DefaultSchedulerType
struct DefaultScheduler{FT} <: EnsembleKalmanProcesses.LearningRateScheduler

Scheduler containing a default constant step size, users can override this temporarily within update_ensemble!.

  • Δt_default::Any: step size
source
EnsembleKalmanProcesses.MutableSchedulerType
struct MutableScheduler{FT} <: EnsembleKalmanProcesses.LearningRateScheduler

Scheduler containing a mutable constant step size, users can override this permanently within update_ensemble!.

  • Δt_mutable::Vector: mutable step size
source
EnsembleKalmanProcesses.EKSStableSchedulerType
struct EKSStableScheduler{FT} <: EnsembleKalmanProcesses.LearningRateScheduler

Scheduler known to be stable for EKS, In particular, $\Delta t = \frac{\alpha}{\|U\| + \varepsilon}$ where $U = (G(u) - \bar{G(u)})^T\Gamma^{-1}(G(u) - y)$. Cannot be overriden.

  • numerator::Any: the numerator $\alpha$

  • nugget::Any: the nugget term $\varepsilon$

source
EnsembleKalmanProcesses.DataMisfitControllerType
struct DataMisfitController{FT, S} <: EnsembleKalmanProcesses.LearningRateScheduler

Scheduler from Iglesias, Yang, 2021, Based on Bayesian Tempering. Terminates at T=1 by default, and at this time, ensemble spread provides a (more) meaningful approximation of posterior uncertainty In particular, for parameters $\theta_j$ at step $n$, to calculate the next timestep $\Delta t_n = \min\left(\max\left(\frac{J}{2\Phi}, \sqrt{\frac{J}{2\langle \Phi, \Phi \rangle}}\right), 1-\sum^{n-1}_i t_i\right)$ where $\Phi_j = \|\Gamma^{-\frac{1}{2}}(G(\theta_j) - y)\|^2$. Cannot be overriden by user provided timesteps. By default termination returns true from update_ensemble! and

  • if on_terminate == "stop", stops further iteration.
  • if on_terminate == "continue_fixed", continues iteration with the final timestep fixed
  • if on_terminate == "continue", continues the algorithm (though no longer compares to $1-\sum^{n-1}_i t_i$)

The user may also change the T with terminate_at keyword.

  • iteration::Vector{Int64}: the current iteration

  • terminate_at::Any: the algorithm time for termination, default: 1.0

  • on_terminate::Any: the action on termination, default: "stop",

source
EnsembleKalmanProcesses.calculate_timestep!Function
calculate_timestep!(
    ekp::EnsembleKalmanProcess,
    g::AbstractMatrix,
    Δt_new::Union{Nothing, AbstractFloat}
) -> Any

Calculates next timestep by pushing to ekp.Δt, !isnothing(return_value) implies termination condition has been met

source

Failure and NaN handling

EnsembleKalmanProcesses.FailureHandlerType
FailureHandler{P <: Process, FM <: FailureHandlingMethod}

Structure defining the failure handler method used in the EnsembleKalmanProcess.

Fields

  • failsafe_update::Function: Failsafe algorithmic update equation

Constructors

FailureHandler(process, method)

defined at /home/runner/work/EnsembleKalmanProcesses.jl/EnsembleKalmanProcesses.jl/src/EnsembleKalmanInversion.jl:97.

FailureHandler(process, method)

defined at /home/runner/work/EnsembleKalmanProcesses.jl/EnsembleKalmanProcesses.jl/src/EnsembleKalmanInversion.jl:110.

FailureHandler(process, method)

defined at /home/runner/work/EnsembleKalmanProcesses.jl/EnsembleKalmanProcesses.jl/src/EnsembleKalmanSampler.jl:100.

FailureHandler(process, method)

defined at /home/runner/work/EnsembleKalmanProcesses.jl/EnsembleKalmanProcesses.jl/src/EnsembleTransformKalmanInversion.jl:106.

FailureHandler(process, method)

defined at /home/runner/work/EnsembleKalmanProcesses.jl/EnsembleKalmanProcesses.jl/src/EnsembleTransformKalmanInversion.jl:118.

FailureHandler(process, method)

defined at /home/runner/work/EnsembleKalmanProcesses.jl/EnsembleKalmanProcesses.jl/src/GaussNewtonKalmanInversion.jl:52.

FailureHandler(process, method)

defined at /home/runner/work/EnsembleKalmanProcesses.jl/EnsembleKalmanProcesses.jl/src/GaussNewtonKalmanInversion.jl:65.

FailureHandler(process, method)

defined at /home/runner/work/EnsembleKalmanProcesses.jl/EnsembleKalmanProcesses.jl/src/SparseEnsembleKalmanInversion.jl:66.

FailureHandler(process, method)

defined at /home/runner/work/EnsembleKalmanProcesses.jl/EnsembleKalmanProcesses.jl/src/SparseEnsembleKalmanInversion.jl:78.

FailureHandler(process, method)

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

FailureHandler(process, method)

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

FailureHandler(process, method)

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

FailureHandler(process, method)

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

source
EnsembleKalmanProcesses.SampleSuccGaussType
SampleSuccGauss <: FailureHandlingMethod

Failure handling method that substitutes failed ensemble members by new samples from the empirical Gaussian distribution defined by the updated successful ensemble.

source
EnsembleKalmanProcesses.sample_empirical_gaussianFunction
sample_empirical_gaussian(
    rng::AbstractRNG,
    u::AbstractMatrix{FT},
    n::IT;
    inflation::Union{FT, Nothing} = nothing,
) where {FT <: Real, IT <: Int}

Returns n samples from an empirical Gaussian based on point estimates u, adding inflation if the covariance is singular.

source
EnsembleKalmanProcesses.split_indices_by_successFunction
 split_indices_by_success(g::AbstractMatrix{FT}) where {FT <: Real}

Returns the successful/failed particle indices given a matrix with output vectors stored as columns. Failures are defined for particles containing at least one NaN output element.

source
EnsembleKalmanProcesses.impute_over_nansFunction
impute_over_nans(
    g::AbstractMatrix,
    nan_tolerance,
    nan_row_values::AbstractVector;
    verbose
) -> Any

Imputation of "reasonable values" over NaNs in the following manner

  1. Detect failures: check if any column contains NaNs exceeding the fraction nan_tolerance, such members are flagged as failures
  2. Impute values in rows with few NaNs: Of the admissible columns, any NaNs are replaced by finite values of the ensemble-mean (without NaNs) over the row.
  3. Impute a value for row with all NaNs: Of the admissible columns, the value of the observation itself y is imputed
source

Process-specific

EnsembleKalmanProcesses.get_prior_meanFunction
get_prior_mean(
    process::Inversion
) -> Union{Nothing, AbstractVector}

Return the stored prior_mean from the Inversion process.

source
get_prior_mean(
    process::TransformInversion
) -> Union{Nothing, AbstractVector}

Return the stored prior_mean from the TransformInversion process.

source
get_prior_mean(
    process::GaussNewtonInversion
) -> AbstractVector

Return the stored prior_mean from the GaussNewtonInversion process.

source
get_prior_mean(process::SparseInversion)

Return nothing; SparseInversion does not store a prior mean.

source
get_prior_mean(
    process::Sampler
) -> Vector{FT} where FT<:AbstractFloat

Return the prior mean vector stored in process.

source
get_prior_mean(
    process::Union{TransformUnscented, Unscented}
) -> Any

Return the prior mean stored in an Unscented or TransformUnscented process.

source
EnsembleKalmanProcesses.get_prior_covFunction
get_prior_cov(
    process::Inversion
) -> Union{Nothing, LinearAlgebra.UniformScaling, AbstractMatrix}

Return the stored prior_cov from the Inversion process.

source
get_prior_cov(
    process::TransformInversion
) -> Union{Nothing, LinearAlgebra.UniformScaling, AbstractMatrix}

Return the stored prior_cov from the TransformInversion process.

source
get_prior_cov(
    process::GaussNewtonInversion
) -> Union{LinearAlgebra.UniformScaling, AbstractMatrix}

Return the stored prior_cov from the GaussNewtonInversion process.

source
get_prior_cov(process::SparseInversion)

Return nothing; SparseInversion does not store a prior covariance.

source
get_prior_cov(
    process::Sampler
) -> Union{LinearAlgebra.UniformScaling{FT}, AbstractMatrix{FT}} where FT<:AbstractFloat

Return the prior covariance matrix stored in process.

source
get_prior_cov(
    process::Union{TransformUnscented, Unscented}
) -> Any

Return the prior covariance stored in an Unscented or TransformUnscented process.

source
EnsembleKalmanProcesses.get_impose_priorFunction
get_impose_prior(process::Inversion) -> Bool

Return the stored impose_prior from the Inversion process.

source
get_impose_prior(process::TransformInversion) -> Bool

Return the stored impose_prior from the TransformInversion process.

source
EnsembleKalmanProcesses.get_bufferFunction
get_buffer(p::TransformInversion) -> AbstractVector

Return the stored buffer from the TransformInversion process.

source
get_buffer(p::TransformUnscented) -> AbstractVector

Return the computation buffer stored in the TransformUnscented process.

The buffer is a two-element vector: buffer[1] holds Y' * Γ⁻¹ and buffer[2] holds Y' * Γ⁻¹ * Y, both pre-allocated to avoid repeated allocation during the analysis step.

source
EnsembleKalmanProcesses.get_default_multiplicative_inflationFunction
get_default_multiplicative_inflation(
    process::Inversion
) -> AbstractFloat

Return the stored default_multiplicative_inflation from the Inversion process.

source
get_default_multiplicative_inflation(
    p::TransformInversion
) -> AbstractFloat

Return the stored default_multiplicative_inflation from the TransformInversion process.

source

Update Groups

EnsembleKalmanProcesses.UpdateGroupType

Container of indices indicating which parameters (u_group) are updated using which data (g_group).

Provide an array of UpdateGroups to partition the parameter space. This partitioning assumes conditional independence between the different u_groups.

struct UpdateGroup

Fields

  • u_group::Vector{Int64}: vector of parameter indices forming a partition of 1:input_dim with other UpdateGroups

  • g_group::Vector{Int64}: vector of data indices within 1:output_dim

  • group_id::Dict: mapping of parameter index range to data index range, used for group identification

Constructors

UpdateGroup(u_group, g_group, group_id)
UpdateGroup(u_group, g_group, group_id)

defined at /home/runner/work/EnsembleKalmanProcesses.jl/EnsembleKalmanProcesses.jl/src/UpdateGroup.jl:24.

UpdateGroup(u_group, g_group)

defined at /home/runner/work/EnsembleKalmanProcesses.jl/EnsembleKalmanProcesses.jl/src/UpdateGroup.jl:35.

source
EnsembleKalmanProcesses.update_group_consistencyFunction
update_group_consistency(
    groups::AbstractVector,
    input_dim::Int64,
    output_dim::Int64
) -> Bool

Check the consistency of sizing and partitioning of the UpdateGroup array if it contains indices No consistency check if u,g has strings internally

source
EnsembleKalmanProcesses.create_update_groupsFunction
create_update_groups(
    prior::EnsembleKalmanProcesses.ParameterDistributions.ParameterDistribution,
    observation::Observation,
    group_identifiers::Dict
) -> Vector{Any}

To construct a list of update-groups populated by indices of parameter distributions and indices of observations, from a dictionary of group_identifiers = Dict(..., group_k_input_names => group_k_output_names, ...)

source

Accelerators

EnsembleKalmanProcesses.ConstantNesterovAcceleratorType
mutable struct ConstantNesterovAccelerator{FT<:AbstractFloat} <: EnsembleKalmanProcesses.Accelerator

Accelerator that adapts a variant of Nesterov's momentum method for EKI. In this variant the momentum parameter is a constant value $λ$, the default of 0.9.

  • λ::AbstractFloat: constant momentum parameter λ ∈ (0, 1)

  • u_prev::Array{FT} where FT<:AbstractFloat: ensemble parameter matrix from the previous iteration, used to compute the momentum update

source
EnsembleKalmanProcesses.NesterovAcceleratorType
mutable struct NesterovAccelerator{FT<:AbstractFloat} <: EnsembleKalmanProcesses.Accelerator

Accelerator that adapts a variant of Nesterov's momentum method for EKI using the recurrence $\theta_{k+1} = \frac{-\theta_k^2 + \sqrt{\theta_k^4 + 4\theta_k^2}}{2}$ with $\theta_0 = 1$, giving momentum coefficient $\theta_{k+1}(\theta_k^{-1} - 1)$ at each iteration. This recurrence is also described in Su, Boyd, Candes (2014) and is the recommended variant. Stores a previous ensemble state for momentum computation; this is distinct from the state returned as the ensemble value.

  • u_prev::Array{FT} where FT<:AbstractFloat: ensemble parameter matrix from the previous iteration, used to compute the momentum update

  • θ_prev::AbstractFloat: previous value of the Nesterov momentum parameter θ

source
EnsembleKalmanProcesses.FirstOrderNesterovAcceleratorType
mutable struct FirstOrderNesterovAccelerator{FT<:AbstractFloat} <: EnsembleKalmanProcesses.Accelerator

Accelerator that adapts a variant of Nesterov's momentum method for EKI. In this variant the momentum parameter at iteration $k$ is given by $1-\frac{r}{k+2}$, a first-order representation of the asymptotic schedule $1-\frac{r}{k}-\mathcal{O}(k^{-2})$ advocated in Su, Boyd, Candes (2014). Stores a previous ensemble state for momentum computation; this is distinct from the state returned as the ensemble value.

  • r::AbstractFloat: decay exponent r controlling the momentum schedule; the momentum coefficient at iteration k is 1 - r/(k+2)

  • u_prev::Array{FT} where FT<:AbstractFloat: ensemble parameter matrix from the previous iteration, used to compute the momentum update

source
EnsembleKalmanProcesses.accelerate!Function
accelerate!(
    ekp::EnsembleKalmanProcess{FT<:AbstractFloat, IT<:Int64, P<:EnsembleKalmanProcesses.Process, LRS<:EnsembleKalmanProcesses.LearningRateScheduler, DefaultAccelerator},
    u::AbstractMatrix
) -> Array{EnsembleKalmanProcesses.DataContainers.DataContainer{FT}, 1} where FT<:AbstractFloat

Push ensemble state u unchanged into ekp, performing no momentum acceleration.

source
accelerate!(
    ekp::EnsembleKalmanProcess{FT<:AbstractFloat, IT<:Int64, P<:EnsembleKalmanProcesses.Process, LRS<:EnsembleKalmanProcesses.LearningRateScheduler, NesterovAccelerator{FT<:AbstractFloat}},
    u::AbstractMatrix
) -> Any

Apply the NesterovAccelerator momentum update to ensemble state u and push the accelerated state into ekp.

The momentum coefficient follows the recurrence described in NesterovAccelerator; see also the lecture notes for background.

source
accelerate!(
    uki::EnsembleKalmanProcess{FT<:AbstractFloat, IT<:Int64, P<:Unscented, LRS<:EnsembleKalmanProcesses.LearningRateScheduler, NesterovAccelerator{FT<:AbstractFloat}},
    u::AbstractMatrix
) -> Any

Apply the NesterovAccelerator momentum update to UKI ensemble state u and push the accelerated state into uki.

Performs the same momentum step as the generic method, then reconstructs the ensemble mean and covariance of the accelerated sigma points by inverting the UKI prediction transformation, and overwrites the stored u_mean and uu_cov.

source
accelerate!(
    ekp::EnsembleKalmanProcess{FT<:AbstractFloat, IT<:Int64, P<:EnsembleKalmanProcesses.Process, LRS<:EnsembleKalmanProcesses.LearningRateScheduler, ConstantNesterovAccelerator{FT<:AbstractFloat}},
    u::AbstractMatrix
) -> Any

Apply the ConstantNesterovAccelerator momentum update to ensemble state u and push the accelerated state into ekp.

The momentum coefficient is the constant λ stored in the accelerator.

source
accelerate!(
    uki::EnsembleKalmanProcess{FT<:AbstractFloat, IT<:Int64, P<:Unscented, LRS<:EnsembleKalmanProcesses.LearningRateScheduler, ConstantNesterovAccelerator{FT<:AbstractFloat}},
    u::AbstractMatrix
) -> Any

Apply the ConstantNesterovAccelerator momentum update to UKI ensemble state u and push the accelerated state into uki.

Performs the same constant-momentum step as the generic method, then reconstructs the ensemble mean and covariance of the accelerated sigma points by inverting the UKI prediction transformation, and overwrites the stored u_mean and uu_cov.

source
accelerate!(
    ekp::EnsembleKalmanProcess{FT<:AbstractFloat, IT<:Int64, P<:EnsembleKalmanProcesses.Process, LRS<:EnsembleKalmanProcesses.LearningRateScheduler, FirstOrderNesterovAccelerator{FT<:AbstractFloat}},
    u::AbstractMatrix
) -> Any

Apply the FirstOrderNesterovAccelerator momentum update to ensemble state u and push the accelerated state into ekp.

The momentum coefficient at iteration k is $1 - r/(k+2)$, where r is the decay exponent stored in the accelerator.

source
accelerate!(
    uki::EnsembleKalmanProcess{FT<:AbstractFloat, IT<:Int64, P<:Unscented, LRS<:EnsembleKalmanProcesses.LearningRateScheduler, FirstOrderNesterovAccelerator{FT<:AbstractFloat}},
    u::AbstractMatrix
) -> Any

Apply the FirstOrderNesterovAccelerator momentum update to UKI ensemble state u and push the accelerated state into uki.

Performs the same iteration-dependent momentum step as the generic method, then reconstructs the ensemble mean and covariance of the accelerated sigma points by inverting the UKI prediction transformation, and overwrites the stored u_mean and uu_cov.

source