EnsembleKalmanProcess
Primary objects and functions
EnsembleKalmanProcesses.EnsembleKalmanProcess — TypeEnsembleKalmanProcess{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:209.
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:244.
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:342.
EnsembleKalmanProcess(params, observation, args; kwargs...)defined at /home/runner/work/EnsembleKalmanProcesses.jl/EnsembleKalmanProcesses.jl/src/EnsembleKalmanProcess.jl:401.
EnsembleKalmanProcess(
params,
obs,
obs_noise_cov,
args;
kwargs...
)defined at /home/runner/work/EnsembleKalmanProcesses.jl/EnsembleKalmanProcesses.jl/src/EnsembleKalmanProcess.jl:411.
EnsembleKalmanProcess(
observation_series,
process;
kwargs...
)defined at /home/runner/work/EnsembleKalmanProcesses.jl/EnsembleKalmanProcesses.jl/src/UnscentedKalmanInversion.jl:256.
EnsembleKalmanProcess(observation, process; kwargs...)defined at /home/runner/work/EnsembleKalmanProcesses.jl/EnsembleKalmanProcesses.jl/src/UnscentedKalmanInversion.jl:267.
EnsembleKalmanProcess(
obs_mean,
obs_noise_cov,
process;
kwargs...
)defined at /home/runner/work/EnsembleKalmanProcesses.jl/EnsembleKalmanProcesses.jl/src/UnscentedKalmanInversion.jl:278.
EnsembleKalmanProcesses.construct_initial_ensemble — Functionconstruct_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
Constructs the initial ensemble for the Unscented or TransformUnscented process. Returned with parameters as columns in unconstrained space by default (constrain by settingconstrained=true`)
NOTE: This function is created just to see what the initial sigma_ensemble will be without constructing the EKP object. Do not pass the initial ensemble into the EnsembleKalmanProcess object.
EnsembleKalmanProcesses.update_ensemble! — Functionupdate_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!.
update_ensemble!(
ekp::EnsembleKalmanProcess{FT, IT, II<:Inversion},
g::AbstractArray{FT, 2},
process::Inversion,
u_idx::Vector{Int64},
g_idx::Vector{Int64};
deterministic_forward_map,
failed_ens,
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 aN_obs × N_ensarray (i.e data are columms).process:: Type of the EKP.u_idx:: indices of u to update (seeUpdateGroup)g_idx:: indices of g,y,Γ with which to update u (seeUpdateGroup)deterministic_forward_map:: Whether outputgcomes from a deterministic model.failed_ens:: Indices of failed particles. If nothing, failures are computed as columns ofgwith NaN entries.
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_ensarray (i.e data are columms). - process :: Type of the EKP.
- u_idx :: indices of u to update (see
UpdateGroup) - g_idx :: indices of g,y,Γ with which to update u (see
UpdateGroup) - failed_ens :: Indices of failed particles. If nothing, failures are computed as columns of
gwith NaN entries.
update_ensemble!(
ekp::EnsembleKalmanProcess{FT, IT, GaussNewtonInversion},
g::AbstractMatrix{FT},
process::GaussNewtonInversion;
deterministic_forward_map::Bool = true,
failed_ens = nothing,
) where {FT, IT}Updates the ensemble according to an GaussNewtonInversion process. The specific update is given by Algorithm 4.1 of Chada, Chen, Sanz-Alonso (2021), https://doi.org/10.3934/fods.2021011, a.k.a Iterated Ensemble Kalman Filter with Statistical Linearization (IEKF-SL)
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). - process :: Type of the EKP.
- deterministicforwardmap :: Whether output
gcomes from a deterministic model. - failed_ens :: Indices of failed particles. If nothing, failures are computed as columns of
gwith NaN entries.
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 aN_obs × N_ensarray (i.e data are columms).process:: Type of the EKP.u_idx:: indices of u to update (seeUpdateGroup)g_idx:: indices of g,y,Γ with which to update u (seeUpdateGroup)deterministic_forward_map:: Whether outputgcomes from a deterministic model.failed_ens:: Indices of failed particles. If nothing, failures are computed as columns ofgwith NaN entries.
update_ensemble!(
ekp::EnsembleKalmanProcess{FT, IT, Sampler{FT,ST}},
g::AbstractMatrix{FT},
process::Sampler{FT, ST};
failed_ens = nothing,
) where {FT, IT, ST}Updates the ensemble according to a Sampler process.
Inputs:
ekp:: The EnsembleKalmanProcess to update.g:: Model outputs, they need to be stored as aN_obs × N_ensarray (i.e data are columms).process:: Type of the EKP.u_idx:: indices of u to update (seeUpdateGroup)g_idx:: indices of g,y,Γ with which to update u (seeUpdateGroup)failed_ens:: Indices of failed particles. If nothing, failures are computed as columns ofgwith NaN entries.
update_ensemble!(
uki::EnsembleKalmanProcess{FT<:AbstractFloat, IT<:Int64, TU<:TransformUnscented},
g::AbstractArray{FT<:AbstractFloat, 2},
process::TransformUnscented,
u_idx::Vector{Int64},
g_idx::Vector{Int64};
group_idx,
failed_ens,
kwargs...
) -> Any
Updates the ensemble according to an TransformUnscented process.
Inputs:
uki:: The EnsembleKalmanProcess to update.g_in:: Model outputs, they need to be stored as aN_obs × N_ensarray (i.e data are columms).process:: Type of the EKP.u_idx:: indices of u to update (seeUpdateGroup)g_idx:: indices of g,y,Γ with which to update u (seeUpdateGroup)group_idx:: the label of the update group (1 is "first update this iteration")failed_ens:: Indices of failed particles. If nothing, failures are computed as columns ofgwith NaN entries.
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 aN_obs × N_ensarray (i.e data are columms).process:: Type of the EKP.u_idx:: indices of u to update (seeUpdateGroup)g_idx:: indices of g,y,Γ with which to update u (seeUpdateGroup)group_idx:: the label of the update group (1 is "first update this iteration")failed_ens:: Indices of failed particles. If nothing, failures are computed as columns ofgwith NaN entries.
Getter functions
EnsembleKalmanProcesses.get_u — Functionget_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 — Functionget_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_ϕ — Functionget_ϕ(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 — Methodget_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 — Methodget_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 — Methodget_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 — Methodget_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 — Methodget_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 — Methodget_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 — Functionget_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 — Functionget_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 — Functionget_ϕ_final(ekp::EnsembleKalmanProcess; return_array=true)Get the constrained parameters at the last iteration.
EnsembleKalmanProcesses.get_u_mean — Functionget_u_mean(ekp::EnsembleKalmanProcess, iteration::IT) where {IT <: Integer}Returns the mean unconstrained parameter at the given iteration.
get_u_mean(uki::EnsembleKalmanProcess{FT, IT, Unscented}, iteration::IT)Returns the mean unconstrained parameter at the requested iteration.
EnsembleKalmanProcesses.get_u_cov — Functionget_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, IT, Unscented}, iteration::IT)Returns the unconstrained parameter covariance at the requested iteration.
EnsembleKalmanProcesses.get_g_mean — Functionget_g_mean(ekp::EnsembleKalmanProcess, iteration::IT) where {IT <: Integer}Returns the mean forward map evaluation at the given iteration.
EnsembleKalmanProcesses.get_ϕ_mean — Functionget_ϕ_mean(prior::ParameterDistribution, ekp::EnsembleKalmanProcess, iteration::IT)Returns the constrained transform of the mean unconstrained parameter at the given iteration.
EnsembleKalmanProcesses.get_u_mean_final — Functionget_u_mean_final(ekp::EnsembleKalmanProcess)Get the mean unconstrained parameter at the last iteration.
EnsembleKalmanProcesses.get_u_cov_final — Functionget_u_cov_final(ekp::EnsembleKalmanProcess)Get the mean unconstrained parameter covariance at the last iteration.
EnsembleKalmanProcesses.get_g_mean_final — Functionget_g_mean_final(ekp::EnsembleKalmanProcess)Get the mean forward model evaluation at the last iteration.
EnsembleKalmanProcesses.get_ϕ_mean_final — Functionget_ϕ_mean_final(prior::ParameterDistribution, ekp::EnsembleKalmanProcess)Get the constrained transform of the mean unconstrained parameter at the last iteration.
EnsembleKalmanProcesses.get_N_iterations — Functionget_N_iterations(ekp::EnsembleKalmanProcess)Get number of times update has been called (equals size(g), or size(u)-1).
EnsembleKalmanProcesses.get_N_ens — Functionget_N_ens(ekp::EnsembleKalmanProcess)Return N_ens field of EnsembleKalmanProcess.
EnsembleKalmanProcesses.get_accelerator — Functionget_accelerator(ekp::EnsembleKalmanProcess)Return accelerator field of EnsembleKalmanProcess.
EnsembleKalmanProcesses.get_scheduler — Functionget_scheduler(ekp::EnsembleKalmanProcess)Return scheduler field of EnsembleKalmanProcess.
EnsembleKalmanProcesses.get_process — Functionget_process(ekp::EnsembleKalmanProcess)Return process field of EnsembleKalmanProcess.
EnsembleKalmanProcesses.get_rng — Functionget_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 — Functionget_Δt(ekp::EnsembleKalmanProcess)Return Δt field of EnsembleKalmanProcess.
EnsembleKalmanProcesses.get_failure_handler — Functionget_failure_handler(ekp::EnsembleKalmanProcess)Return failure_handler field of EnsembleKalmanProcess.
EnsembleKalmanProcesses.get_localizer — Functionget_localizer(ekp::EnsembleKalmanProcess)Return localizer field of EnsembleKalmanProcess.
EnsembleKalmanProcesses.get_localizer_type — Functionget_localizer_type(ekp::EnsembleKalmanProcess)Return first parametric type of the localizer field of EnsembleKalmanProcess.
EnsembleKalmanProcesses.get_nan_tolerance — Functionget_nan_tolerance(ekp::EnsembleKalmanProcess)Return nan_tolerance field of EnsembleKalmanProcess.
EnsembleKalmanProcesses.get_nan_row_values — Functionget_nan_row_values(ekp::EnsembleKalmanProcess)Return nan_row_values field of EnsembleKalmanProcess.
EnsembleKalmanProcesses.list_update_groups_over_minibatch — Functionlist_update_groups_over_minibatch(ekp::EnsembleKalmanProcess)Return ugroups and ggroups for the current minibatch, i.e. the subset of
Error metrics
EnsembleKalmanProcesses.compute_average_rmse — Functioncompute_average_rmse(ekp::EnsembleKalmanProcess) -> Any
Computes the average covariance-weighted error over the forward model ensemble, normalized by the dimension: 1/dim(y) * tr((g_i - y)' * Γ⁻¹ * (g_i - y)). The error is retrievable as get_error_metrics(ekp)["avg_rmse"]
For Unscented processes it doesn't make sense to average RMSE at sigma points, so it is evaluated at the mean only, where it is exactly equal to sqrt(compute_loss_at_mean(uki))
EnsembleKalmanProcesses.compute_loss_at_mean — Functioncompute_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 — Functioncompute_average_unweighted_rmse(
ekp::EnsembleKalmanProcess
) -> Any
Computes the average unweighted error over the forward model ensemble, normalized by the dimension: 1/dim(y) * tr((g_i - y)' * (g_i - y)). The error is retrievable as get_error_metrics(ekp)["unweighted_avg_rmse"]
For Unscented processes it doesn't make sense to average unweighted RMSE at sigma points, so it is evaluated at the mean only, where it is exactly equal to sqrt(compute_unweighted_loss_at_mean(uki))
EnsembleKalmanProcesses.compute_unweighted_loss_at_mean — Functioncompute_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 — Functioncompute_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 — Functioncompute_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! — Functioncompute_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 — Functionget_error_metrics(ekp::EnsembleKalmanProcess)Returns the stored error_metrics, created with compute_error!
EnsembleKalmanProcesses.get_error — Functionget_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 — Functionlmul_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 — Functionlmul_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! — Functionlmul_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! — Functionlmul_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 — Typestruct 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 — Typestruct 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 — Typestruct 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 — Typestruct 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! — Functioncalculate_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
Failure and NaN handling
EnsembleKalmanProcesses.FailureHandler — TypeFailureHandler{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:82.
FailureHandler(process, method)defined at /home/runner/work/EnsembleKalmanProcesses.jl/EnsembleKalmanProcesses.jl/src/EnsembleKalmanInversion.jl:95.
FailureHandler(process, method)defined at /home/runner/work/EnsembleKalmanProcesses.jl/EnsembleKalmanProcesses.jl/src/EnsembleKalmanSampler.jl:62.
FailureHandler(process, method)FailureHandler(process, method)FailureHandler(process, method)defined at /home/runner/work/EnsembleKalmanProcesses.jl/EnsembleKalmanProcesses.jl/src/GaussNewtonKalmanInversion.jl:28.
FailureHandler(process, method)defined at /home/runner/work/EnsembleKalmanProcesses.jl/EnsembleKalmanProcesses.jl/src/GaussNewtonKalmanInversion.jl:41.
FailureHandler(process, method)FailureHandler(process, method)FailureHandler(process, method)defined at /home/runner/work/EnsembleKalmanProcesses.jl/EnsembleKalmanProcesses.jl/src/UnscentedKalmanInversion.jl:316.
FailureHandler(process, method)defined at /home/runner/work/EnsembleKalmanProcesses.jl/EnsembleKalmanProcesses.jl/src/UnscentedKalmanInversion.jl:335.
FailureHandler(process, method)FailureHandler(process, method)EnsembleKalmanProcesses.SampleSuccGauss — TypeSampleSuccGauss <: FailureHandlingMethodFailure handling method that substitutes failed ensemble members by new samples from the empirical Gaussian distribution defined by the updated successful ensemble.
EnsembleKalmanProcesses.IgnoreFailures — TypeFailure handling method that ignores forward model failures
EnsembleKalmanProcesses.sample_empirical_gaussian — Functionsample_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 — Functionimpute_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 — Functionget_prior_mean(
process::Inversion
) -> Union{Nothing, AbstractVector}
Returns the stored prior_mean from the Inversion process
get_prior_mean(
process::TransformInversion
) -> Union{Nothing, AbstractVector}
Returns the stored prior_mean from the TransformInversion process
EnsembleKalmanProcesses.get_prior_cov — Functionget_prior_cov(
process::Inversion
) -> Union{Nothing, LinearAlgebra.UniformScaling, AbstractMatrix}
Returns the stored prior_cov from the Inversion process
get_prior_cov(
process::TransformInversion
) -> Union{Nothing, LinearAlgebra.UniformScaling, AbstractMatrix}
Returns the stored prior_cov from the TransformInversion process
EnsembleKalmanProcesses.get_impose_prior — Functionget_impose_prior(process::Inversion) -> Bool
Returns the stored impose_prior from the Inversion process
get_impose_prior(process::TransformInversion) -> Bool
Returns the stored impose_prior from the TransformInversion process
EnsembleKalmanProcesses.get_buffer — Functionget_buffer(p::TransformInversion) -> AbstractVector
Returns the stored buffer from the TransformInversion process
get_buffer(p::TransformUnscented) -> AbstractVector
Returns the stored buffer from the TransformUnscented process
EnsembleKalmanProcesses.get_default_multiplicative_inflation — Functionget_default_multiplicative_inflation(
process::Inversion
) -> AbstractFloat
Returns the stored default_multiplicative_inflation from the Inversion process
get_default_multiplicative_inflation(
p::TransformInversion
) -> AbstractFloat
Returns the stored default_multiplicative_inflation from the TransformInversion process