EnsembleKalmanProcess
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 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::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
orSampler
orUnscented
orTransformInversion
orSparseInversion
)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},
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(),
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 failureslocalization_method
:: Method used to localize sample covariancesverbose
:: Whether to print diagnostic information
Other constructors:
EnsembleKalmanProcess(u, observation_series, N_ens, g, error, scheduler, accelerator, Δt, process, rng, failure_handler, localizer, verbose)
defined at /home/runner/work/EnsembleKalmanProcesses.jl/EnsembleKalmanProcesses.jl/src/EnsembleKalmanProcess.jl:159
.
EnsembleKalmanProcess(params, observation_series, process, configuration)
defined at /home/runner/work/EnsembleKalmanProcesses.jl/EnsembleKalmanProcesses.jl/src/EnsembleKalmanProcess.jl:188
.
EnsembleKalmanProcess(params, observation_series, process)
defined at /home/runner/work/EnsembleKalmanProcesses.jl/EnsembleKalmanProcesses.jl/src/EnsembleKalmanProcess.jl:268
.
EnsembleKalmanProcess(params, observation, args)
defined at /home/runner/work/EnsembleKalmanProcesses.jl/EnsembleKalmanProcesses.jl/src/EnsembleKalmanProcess.jl:312
.
EnsembleKalmanProcess(params, obs, obs_noise_cov, args)
defined at /home/runner/work/EnsembleKalmanProcesses.jl/EnsembleKalmanProcesses.jl/src/EnsembleKalmanProcess.jl:322
.
EnsembleKalmanProcess(observation_series, process)
defined at /home/runner/work/EnsembleKalmanProcesses.jl/EnsembleKalmanProcesses.jl/src/UnscentedKalmanInversion.jl:222
.
EnsembleKalmanProcess(observation, process)
defined at /home/runner/work/EnsembleKalmanProcesses.jl/EnsembleKalmanProcesses.jl/src/UnscentedKalmanInversion.jl:233
.
EnsembleKalmanProcess(obs_mean, obs_noise_cov, process)
defined at /home/runner/work/EnsembleKalmanProcesses.jl/EnsembleKalmanProcesses.jl/src/UnscentedKalmanInversion.jl:244
.
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: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)
FailureHandler(process, method)
FailureHandler(process, method)
FailureHandler(process, method)
FailureHandler(process, method)
defined at /home/runner/work/EnsembleKalmanProcesses.jl/EnsembleKalmanProcesses.jl/src/UnscentedKalmanInversion.jl:255
.
FailureHandler(process, method)
defined at /home/runner/work/EnsembleKalmanProcesses.jl/EnsembleKalmanProcesses.jl/src/UnscentedKalmanInversion.jl:274
.
EnsembleKalmanProcesses.SampleSuccGauss
— Type" SampleSuccGauss <: FailureHandlingMethod
Failure 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.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_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.compute_error!
— Functioncompute_error!(ekp::EnsembleKalmanProcess)
Computes the covariance-weighted error of the mean forward model output, (ḡ - y)'Γ_inv(ḡ - y)
. The error is stored within the EnsembleKalmanProcess
.
EnsembleKalmanProcesses.get_error
— Functionget_error(ekp::EnsembleKalmanProcess)
Returns the mean forward model output error as a function of algorithmic time.
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
— Functiongets the rng
field from the FixedMinibatcher
object
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_failuer_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.construct_initial_ensemble
— Functionconstruct_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.
EnsembleKalmanProcesses.update_ensemble!
— Functionupdate_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!.
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.
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.
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_ens
array (i.e data are columms).process
:: Type of the EKP.deterministic_forward_map
:: Whether outputg
comes from a deterministic model.failed_ens
:: Indices of failed particles. If nothing, failures are computed as columns ofg
with NaN entries.
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.
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_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 ofg
with NaN entries.
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.
Learning Rate Schedulers
EnsembleKalmanProcesses.DefaultScheduler
— Typestruct DefaultScheduler{FT} <: EnsembleKalmanProcesses.LearningRateScheduler
Scheduler 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.LearningRateScheduler
Scheduler 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.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$
EnsembleKalmanProcesses.DataMisfitController
— Typestruct 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",
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