EnsembleKalmanProcess

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]

  • obs_mean::Vector{FT} where FT<:AbstractFloat

    vector of the observed vector size [N_obs]

  • obs_noise_cov::Union{LinearAlgebra.UniformScaling{FT}, AbstractMatrix{FT}} where FT<:AbstractFloat

    covariance matrix of the observational noise, of size [N_obs × N_obs]

  • 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]

  • err::Vector{FT} where FT<:AbstractFloat

    vector of errors

  • 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

  • process::EnsembleKalmanProcesses.Process

    the particular EK process (Inversion or Sampler or Unscented or TransformInversion or SparseInversion)

  • 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, implemented for (Inversion, SparseInversion, Unscented, TransformInversion)

  • localizer::EnsembleKalmanProcesses.Localizers.Localizer

    Localization kernel, implemented for (Inversion, SparseInversion, Unscented)

  • verbose::Bool

    Whether to print diagnostics for each EK iteration

Generic constructor

EnsembleKalmanProcess(
    params::AbstractMatrix{FT},
    obs_mean,
    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(),
    localization_method::LM = NoLocalization(),
    verbose::Bool = false,
) where {FT <: AbstractFloat, P <: Process, FM <: FailureHandlingMethod, LM <: LocalizationMethod}

Inputs:

  • params :: Initial parameter ensemble
  • obs_mean :: Vector of observations
  • obs_noise_cov :: Noise covariance associated with the observations obs_mean
  • 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
  • localization_method :: Method used to localize sample covariances
  • verbose :: Whether to print diagnostic information

Other constructors:

EnsembleKalmanProcess(u, obs_mean, obs_noise_cov, N_ens, g, err, scheduler, accelerator, Δt, process, rng, failure_handler, localizer, verbose)

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

EnsembleKalmanProcess(params, obs_mean, obs_noise_cov, process)

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

EnsembleKalmanProcess(obs_mean, obs_noise_cov, process)

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

source
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:10.

FailureHandler(process, method)

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

FailureHandler(process, method)

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

FailureHandler(process, method)

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

FailureHandler(process, method)

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

FailureHandler(process, method)

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

FailureHandler(process, method)

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

FailureHandler(process, method)

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

FailureHandler(process, method)

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

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.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)

Returns the constrained parameters at the given iteration.

source
get_ϕ(prior::ParameterDistribution, ekp::EnsembleKalmanProcess)

Returns the constrained parameters from all iterations. The outer dimension is given by the number of iterations.

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, IT, Unscented}, iteration::IT)

Returns the mean unconstrained parameter at the requested 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, IT, Unscented}, iteration::IT)

Returns the unconstrained parameter covariance at the requested 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.compute_error!Function
compute_error!(ekp::EnsembleKalmanProcess)

Computes the covariance-weighted error of the mean forward model output, (ḡ - y)'Γ_inv(ḡ - y). The error is stored within the EnsembleKalmanProcess.

source
EnsembleKalmanProcesses.construct_initial_ensembleFunction
construct_initial_ensemble(
    rng::AbstractRNG,
    prior::ParameterDistribution,
    N_ens::IT
) where {IT <: Int}
construct_initial_ensemble(prior::ParameterDistribution, N_ens::IT) where {IT <: Int}

Construct the initial parameters, by sampling N_ens samples from specified prior distribution. Returned with parameters as columns.

source
EnsembleKalmanProcesses.update_ensemble!Function
update_ensemble!(
    ekp::EnsembleKalmanProcess,
    g::AbstractMatrix{FT};
    multiplicative_inflation::Bool = false,
    additive_inflation::Bool = false,
    additive_inflation_cov::MorUS = get_u_cov_prior(ekp),
    s::FT = 0.0,
    ekp_kwargs...,
) where {FT, IT}

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
update_ensemble!(
    ekp::EnsembleKalmanProcess{FT, IT, Inversion},
    g::AbstractMatrix{FT},
    process::Inversion;
    deterministic_forward_map::Bool = true,
    failed_ens = nothing,
) where {FT, IT}

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).
  • process :: Type of the EKP.
  • deterministicforwardmap :: Whether output g comes from a deterministic model.
  • failed_ens :: Indices of failed particles. If nothing, failures are computed as columns of g with NaN entries.
source
update_ensemble!(
    ekp::EnsembleKalmanProcess{FT, IT, TransformInversion},
    g::AbstractMatrix{FT},
    process::TransformInversion;
    failed_ens = nothing,
) where {FT, IT}

Updates the ensemble according to a TransformInversion 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).
  • process :: Type of the EKP.
  • failed_ens :: Indices of failed particles. If nothing, failures are computed as columns of g with NaN entries.
source
update_ensemble!(
    ekp::EnsembleKalmanProcess{FT, IT, SparseInversion{FT}},
    g::AbstractMatrix{FT},
    process::SparseInversion{FT};
    deterministic_forward_map = true,
    failed_ens = nothing,
) where {FT, IT}

Updates the ensemble according to a SparseInversion 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).
  • process :: Type of the EKP.
  • deterministic_forward_map :: Whether output g comes from a deterministic model.
  • failed_ens :: Indices of failed particles. If nothing, failures are computed as columns of g with NaN entries.
source
update_ensemble!(
    ekp::EnsembleKalmanProcess{FT, IT, Sampler{FT}},
    g::AbstractMatrix{FT},
    process::Sampler{FT};
    failed_ens = nothing,
) where {FT, IT}

Updates the ensemble according to a Sampler 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).
  • process :: Type of the EKP.
  • failed_ens :: Indices of failed particles. If nothing, failures are computed as columns of g with NaN entries.
source
update_ensemble!(
    uki::EnsembleKalmanProcess{FT, IT, Unscented},
    g_in::AbstractMatrix{FT},
    process::Unscented;
    failed_ens = nothing,
) where {FT <: AbstractFloat, IT <: Int}

Updates the ensemble according to an Unscented process.

Inputs:

  • uki :: The EnsembleKalmanProcess to update.
  • g_in :: Model outputs, they need to be stored as a N_obs × N_ens array (i.e data are columms).
  • process :: Type of the EKP.
  • failed_ens :: Indices of failed particles. If nothing, failures are computed as columns of g with NaN entries.
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

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, M, 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

  • inv_sqrt_noise::Vector

    the inverse square-root of the noise covariance is stored

  • 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