EnsembleKalmanProcess
Primary objects and functions
EnsembleKalmanProcesses.EnsembleKalmanProcess — Type
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 mechanismN_ens::Int64: ensemble sizeg::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 metricscheduler::EnsembleKalmanProcesses.LearningRateScheduler: Scheduler to calculate the timestep size in each EK iterationaccelerator::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 iterationupdate_groups::AbstractVector: vector of update groups, defining which parameters should be updated by which dataprocess::EnsembleKalmanProcesses.Process: the particular EK process (Inversionetc.)rng::Random.AbstractRNG: Random number generator object (algorithm + seed) used for sampling and noise, for reproducibility. Defaults toRandom.GLOBAL_RNG.failure_handler::FailureHandler: struct storing failsafe update directiveslocalizer::EnsembleKalmanProcesses.Localizers.Localizer: Localization kernel, implemented for (Inversion,SparseInversion,Unscented)nan_tolerance::AbstractFloat: Fraction of allowableNaNs in model output before an ensemble member is considered a failed member and handled by failure handlernan_row_values::Union{Nothing, AbstractVector}: Default-value vector to impute over entire-NaN rows of data (get_obs(ekp)used if valuenothing)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 ensembleobservation_series:: Container for observations (and possible minibatching)process:: Algorithm used to evolve the ensemblescheduler:: Adaptive timestep calculatorΔt:: Initial time step or learning raterng:: Random number generatorfailure_handler_method:: Method used to handle particle failuresnan_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 valuenothing)localization_method:: Method used to localize sample covariancesverbose:: 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.
EnsembleKalmanProcesses.construct_initial_ensemble — Function
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
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.
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_ensarray (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!.
Getter functions
EnsembleKalmanProcesses.get_u — Function
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.
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.
EnsembleKalmanProcesses.get_g — Function
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.
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.
EnsembleKalmanProcesses.get_ϕ — Function
get_ϕ(prior::ParameterDistribution, ekp::EnsembleKalmanProcess, iteration::IT; return_array=true)Returns the constrained parameters at the given iteration.
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.
EnsembleKalmanProcesses.get_obs — Method
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,
EnsembleKalmanProcesses.get_obs — Method
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,
EnsembleKalmanProcesses.get_obs_noise_cov — Method
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,
EnsembleKalmanProcesses.get_obs_noise_cov — Method
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,
EnsembleKalmanProcesses.get_obs_noise_cov_inv — Method
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,
EnsembleKalmanProcesses.get_obs_noise_cov_inv — Method
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,
EnsembleKalmanProcesses.get_u_final — Function
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.
EnsembleKalmanProcesses.get_g_final — Function
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.
EnsembleKalmanProcesses.get_ϕ_final — Function
get_ϕ_final(ekp::EnsembleKalmanProcess; return_array=true)Get the constrained parameters at the last iteration.
EnsembleKalmanProcesses.get_u_mean — Function
get_u_mean(ekp::EnsembleKalmanProcess, iteration::IT) where {IT <: Integer}Returns the mean unconstrained parameter at the given iteration.
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.
EnsembleKalmanProcesses.get_u_cov — Function
get_u_cov(ekp::EnsembleKalmanProcess, iteration::IT) where {IT <: Integer}Returns the unconstrained parameter sample covariance at the given iteration.
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.
EnsembleKalmanProcesses.get_g_mean — Function
get_g_mean(ekp::EnsembleKalmanProcess, iteration::IT) where {IT <: Integer}Returns the mean forward map evaluation at the given iteration.
EnsembleKalmanProcesses.get_ϕ_mean — Function
get_ϕ_mean(prior::ParameterDistribution, ekp::EnsembleKalmanProcess, iteration::IT)Returns the constrained transform of the mean unconstrained parameter at the given iteration.
EnsembleKalmanProcesses.get_u_prior — Function
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.
EnsembleKalmanProcesses.get_u_mean_final — Function
get_u_mean_final(ekp::EnsembleKalmanProcess)Get the mean unconstrained parameter at the last iteration.
EnsembleKalmanProcesses.get_u_cov_prior — Function
get_u_cov_prior(ekp::EnsembleKalmanProcess)Returns the unconstrained parameter sample covariance for the initial ensemble.
EnsembleKalmanProcesses.get_u_cov_final — Function
get_u_cov_final(ekp::EnsembleKalmanProcess)Get the mean unconstrained parameter covariance at the last iteration.
EnsembleKalmanProcesses.get_g_mean_final — Function
get_g_mean_final(ekp::EnsembleKalmanProcess)Get the mean forward model evaluation at the last iteration.
EnsembleKalmanProcesses.get_ϕ_mean_final — Function
get_ϕ_mean_final(prior::ParameterDistribution, ekp::EnsembleKalmanProcess)Get the constrained transform of the mean unconstrained parameter at the last iteration.
EnsembleKalmanProcesses.get_N_iterations — Function
get_N_iterations(ekp::EnsembleKalmanProcess)Get number of times update has been called (equals size(g), or size(u)-1).
EnsembleKalmanProcesses.get_N_ens — Function
get_N_ens(ekp::EnsembleKalmanProcess)Return N_ens field of EnsembleKalmanProcess.
EnsembleKalmanProcesses.get_accelerator — Function
get_accelerator(ekp::EnsembleKalmanProcess)Return accelerator field of EnsembleKalmanProcess.
EnsembleKalmanProcesses.get_scheduler — Function
get_scheduler(ekp::EnsembleKalmanProcess)Return scheduler field of EnsembleKalmanProcess.
EnsembleKalmanProcesses.get_process — Function
get_process(ekp::EnsembleKalmanProcess)Return process field of EnsembleKalmanProcess.
EnsembleKalmanProcesses.get_rng — Function
get_rng(m::FixedMinibatcher) -> Random.AbstractRNG
gets the rng field from the FixedMinibatcher object
get_rng(m::RandomFixedSizeMinibatcher) -> Random.AbstractRNG
gets the rng field from the RandomFixesSizeMinibatcher object
get_rng(ekp::EnsembleKalmanProcess)Return rng field of EnsembleKalmanProcess.
EnsembleKalmanProcesses.get_Δt — Function
get_Δt(ekp::EnsembleKalmanProcess)Return Δt field of EnsembleKalmanProcess.
EnsembleKalmanProcesses.get_algorithm_time — Function
get_algorithm_time(ekp::EnsembleKalmanProcess) -> Any
Get the accumulated Δt over iterations, also known as algorithm- or pseudo-time.
EnsembleKalmanProcesses.get_observation_series — Function
get_observation_series(ekp::EnsembleKalmanProcess)Return observation_series field of EnsembleKalmanProcess.
EnsembleKalmanProcesses.get_failure_handler — Function
get_failure_handler(ekp::EnsembleKalmanProcess)Return failure_handler field of EnsembleKalmanProcess.
EnsembleKalmanProcesses.get_localizer — Function
get_localizer(ekp::EnsembleKalmanProcess)Return localizer field of EnsembleKalmanProcess.
EnsembleKalmanProcesses.get_localizer_type — Function
get_localizer_type(ekp::EnsembleKalmanProcess)Return first parametric type of the localizer field of EnsembleKalmanProcess.
EnsembleKalmanProcesses.get_nan_tolerance — Function
get_nan_tolerance(ekp::EnsembleKalmanProcess)Return nan_tolerance field of EnsembleKalmanProcess.
EnsembleKalmanProcesses.get_nan_row_values — Function
get_nan_row_values(ekp::EnsembleKalmanProcess)Return nan_row_values field of EnsembleKalmanProcess.
EnsembleKalmanProcesses.list_update_groups_over_minibatch — Function
list_update_groups_over_minibatch(ekp::EnsembleKalmanProcess)Return ugroups and ggroups for the current minibatch, i.e. the subset of
EnsembleKalmanProcesses.get_cov_blocks — Function
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.
Error metrics
EnsembleKalmanProcesses.compute_average_rmse — Function
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"]
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.
EnsembleKalmanProcesses.compute_loss_at_mean — Function
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
EnsembleKalmanProcesses.compute_average_unweighted_rmse — Function
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"]
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.
EnsembleKalmanProcesses.compute_unweighted_loss_at_mean — Function
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"]
EnsembleKalmanProcesses.compute_bayes_loss_at_mean — Function
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
EnsembleKalmanProcesses.compute_crps — Function
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).
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_rmsecomputed withcompute_average_rmse(ekp)losscomputed withcompute_loss_at_mean(ekp)unweighed_avg_rmsecomputed withcompute_average_unweighted_rmse(ekp)unweighted_losscomputed withcompute_unweighted_loss_at_mean(ekp)bayes_losscomputed withcompute_bayes_loss_at_mean(ekp)crpscomputed withcompute_crps(ekp)
EnsembleKalmanProcesses.get_error_metrics — Function
get_error_metrics(ekp::EnsembleKalmanProcess)Returns the stored error_metrics, created with compute_error!
EnsembleKalmanProcesses.get_error — Function
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 withcompute_bayes_loss_at_mean(ekp) - If prior not provided to the process:
get_error_metrics(ekp)["loss"], the loss computed withcompute_loss_at_mean(ekp)
EnsembleKalmanProcesses.lmul_obs_noise_cov — Function
lmul_obs_noise_cov(
ekp::EnsembleKalmanProcess,
X::AbstractVecOrMat
) -> Any
Convenience function to multiply X by obs_noise_cov on the left without building the full matrix.
EnsembleKalmanProcesses.lmul_obs_noise_cov_inv — Function
lmul_obs_noise_cov_inv(
ekp::EnsembleKalmanProcess,
X::AbstractVecOrMat
) -> Any
Convenience function to multiply X by obs_noise_cov_inv on the left without building the full matrix.
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_idxfrom the unbuiltobs_noise_cov_inv - select local indices
[:,local_idx]of this block - multiply with the corresponding global indices
X[global_idx,:]
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_idxfrom the unbuiltobs_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()
Learning Rate Schedulers
EnsembleKalmanProcesses.DefaultScheduler — Type
struct DefaultScheduler{FT} <: EnsembleKalmanProcesses.LearningRateSchedulerScheduler containing a default constant step size, users can override this temporarily within update_ensemble!.
Δt_default::Any: step size
EnsembleKalmanProcesses.MutableScheduler — Type
struct MutableScheduler{FT} <: EnsembleKalmanProcesses.LearningRateSchedulerScheduler containing a mutable constant step size, users can override this permanently within update_ensemble!.
Δt_mutable::Vector: mutable step size
EnsembleKalmanProcesses.EKSStableScheduler — Type
struct EKSStableScheduler{FT} <: EnsembleKalmanProcesses.LearningRateSchedulerScheduler 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$
EnsembleKalmanProcesses.DataMisfitController — Type
struct DataMisfitController{FT, S} <: EnsembleKalmanProcesses.LearningRateSchedulerScheduler 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 iterationterminate_at::Any: the algorithm time for termination, default: 1.0on_terminate::Any: the action on termination, default: "stop",
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
EnsembleKalmanProcesses.posdef_correct — Function
posdef_correct(mat::AbstractMatrix; tol) -> Any
Makes square matrix mat positive definite, by symmetrizing and bounding the minimum eigenvalue below by tol
Failure and NaN handling
EnsembleKalmanProcesses.FailureHandler — Type
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)FailureHandler(process, method)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)FailureHandler(process, method)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)FailureHandler(process, method)EnsembleKalmanProcesses.SampleSuccGauss — Type
SampleSuccGauss <: FailureHandlingMethodFailure handling method that substitutes failed ensemble members by new samples from the empirical Gaussian distribution defined by the updated successful ensemble.
EnsembleKalmanProcesses.IgnoreFailures — Type
Failure handling method that ignores forward model failures
EnsembleKalmanProcesses.sample_empirical_gaussian — Function
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.
EnsembleKalmanProcesses.split_indices_by_success — Function
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.
EnsembleKalmanProcesses.impute_over_nans — Function
impute_over_nans(
g::AbstractMatrix,
nan_tolerance,
nan_row_values::AbstractVector;
verbose
) -> Any
Imputation of "reasonable values" over NaNs in the following manner
- Detect failures: check if any column contains NaNs exceeding the fraction
nan_tolerance, such members are flagged as failures - 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.
- Impute a value for row with all NaNs: Of the admissible columns, the value of the observation itself
yis imputed
Process-specific
EnsembleKalmanProcesses.get_prior_mean — Function
get_prior_mean(
process::Inversion
) -> Union{Nothing, AbstractVector}
Return the stored prior_mean from the Inversion process.
get_prior_mean(
process::TransformInversion
) -> Union{Nothing, AbstractVector}
Return the stored prior_mean from the TransformInversion process.
get_prior_mean(
process::GaussNewtonInversion
) -> AbstractVector
Return the stored prior_mean from the GaussNewtonInversion process.
get_prior_mean(process::SparseInversion)
Return nothing; SparseInversion does not store a prior mean.
get_prior_mean(
process::Sampler
) -> Vector{FT} where FT<:AbstractFloat
Return the prior mean vector stored in process.
get_prior_mean(
process::Union{TransformUnscented, Unscented}
) -> Any
Return the prior mean stored in an Unscented or TransformUnscented process.
EnsembleKalmanProcesses.get_prior_cov — Function
get_prior_cov(
process::Inversion
) -> Union{Nothing, LinearAlgebra.UniformScaling, AbstractMatrix}
Return the stored prior_cov from the Inversion process.
get_prior_cov(
process::TransformInversion
) -> Union{Nothing, LinearAlgebra.UniformScaling, AbstractMatrix}
Return the stored prior_cov from the TransformInversion process.
get_prior_cov(
process::GaussNewtonInversion
) -> Union{LinearAlgebra.UniformScaling, AbstractMatrix}
Return the stored prior_cov from the GaussNewtonInversion process.
get_prior_cov(process::SparseInversion)
Return nothing; SparseInversion does not store a prior covariance.
get_prior_cov(
process::Sampler
) -> Union{LinearAlgebra.UniformScaling{FT}, AbstractMatrix{FT}} where FT<:AbstractFloat
Return the prior covariance matrix stored in process.
get_prior_cov(
process::Union{TransformUnscented, Unscented}
) -> Any
Return the prior covariance stored in an Unscented or TransformUnscented process.
EnsembleKalmanProcesses.get_impose_prior — Function
get_impose_prior(process::Inversion) -> Bool
Return the stored impose_prior from the Inversion process.
get_impose_prior(process::TransformInversion) -> Bool
Return the stored impose_prior from the TransformInversion process.
EnsembleKalmanProcesses.get_buffer — Function
get_buffer(p::TransformInversion) -> AbstractVector
Return the stored buffer from the TransformInversion process.
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.
EnsembleKalmanProcesses.get_default_multiplicative_inflation — Function
get_default_multiplicative_inflation(
process::Inversion
) -> AbstractFloat
Return the stored default_multiplicative_inflation from the Inversion process.
get_default_multiplicative_inflation(
p::TransformInversion
) -> AbstractFloat
Return the stored default_multiplicative_inflation from the TransformInversion process.
Update Groups
EnsembleKalmanProcesses.UpdateGroup — Type
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 UpdateGroupFields
u_group::Vector{Int64}: vector of parameter indices forming a partition of1:input_dimwith otherUpdateGroupsg_group::Vector{Int64}: vector of data indices within1:output_dimgroup_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.
EnsembleKalmanProcesses.get_u_group — Function
get_u_group(group::UpdateGroup) -> Vector{Int64}
Return the parameter index vector stored in group.
EnsembleKalmanProcesses.get_g_group — Function
get_g_group(group::UpdateGroup) -> Vector{Int64}
Return the data index vector stored in group.
EnsembleKalmanProcesses.get_group_id — Function
get_group_id(group::UpdateGroup) -> Dict
Return the group identification dictionary stored in group.
EnsembleKalmanProcesses.update_group_consistency — Function
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
EnsembleKalmanProcesses.create_update_groups — Function
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, ...)
Accelerators
EnsembleKalmanProcesses.DefaultAccelerator — Type
struct DefaultAccelerator <: EnsembleKalmanProcesses.AcceleratorDefault accelerator provides no acceleration, runs traditional EKI
EnsembleKalmanProcesses.ConstantNesterovAccelerator — Type
mutable struct ConstantNesterovAccelerator{FT<:AbstractFloat} <: EnsembleKalmanProcesses.AcceleratorAccelerator 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
EnsembleKalmanProcesses.NesterovAccelerator — Type
mutable struct NesterovAccelerator{FT<:AbstractFloat} <: EnsembleKalmanProcesses.AcceleratorAccelerator 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 θ
EnsembleKalmanProcesses.FirstOrderNesterovAccelerator — Type
mutable struct FirstOrderNesterovAccelerator{FT<:AbstractFloat} <: EnsembleKalmanProcesses.AcceleratorAccelerator 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
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.
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.
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.
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.
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.
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.
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.