EnsembleKalmanProcesses.ParameterDistributions.ConstraintType
Constraint{T} <: ConstraintType

Class describing a 1D bijection between constrained and unconstrained spaces. Included parametric types for T:

  • NoConstraint
  • BoundedBelow
  • BoundedAbove
  • Bounded

Fields

  • constrained_to_unconstrained::Function: A map from constrained domain -> (-Inf,Inf)

  • c_to_u_jacobian::Function: The jacobian of the map from constrained domain -> (-Inf,Inf)

  • unconstrained_to_constrained::Function: Map from (-Inf,Inf) -> constrained domain

  • bounds::Union{Nothing, Dict}: Dictionary of values used to build the Constraint (e.g. "lowerbound" or "upperbound")

EnsembleKalmanProcesses.ParameterDistributions.GaussianRandomFieldInterfaceType
struct GaussianRandomFieldInterface <: FunctionParameterDistributionType

GaussianRandomFieldInterface object based on a GRF package. Only a ND->1D output-dimension field interface is implemented.

Fields

  • gaussian_random_field::Any: GRF object, containing the mapping from the field of coefficients to the discrete function

  • package::EnsembleKalmanProcesses.ParameterDistributions.GaussianRandomFieldsPackage: the choice of GRF package

  • distribution::EnsembleKalmanProcesses.ParameterDistributions.ParameterDistribution: the distribution of the coefficients that we shall compute with

EnsembleKalmanProcesses.ParameterDistributions.ParameterDistributionType
ParameterDistribution

Structure to hold a parameter distribution, always stored as an array of distributions internally.

Fields

  • distribution::AbstractVector{PDType} where PDType<:EnsembleKalmanProcesses.ParameterDistributions.ParameterDistributionType: Vector of parameter distributions, defined in unconstrained space

  • constraint::AbstractVector{CType} where CType<:EnsembleKalmanProcesses.ParameterDistributions.ConstraintType: Vector of constraints defining transformations between constrained and unconstrained space

  • name::AbstractVector{ST} where ST<:AbstractString: Vector of parameter names

Constructors

ParameterDistribution(distribution, constraint, name)
ParameterDistribution(param_dist_dict)
ParameterDistribution(distribution, constraint, name)
ParameterDistribution(
    distribution_samples,
    constraint,
    name;
    params_are_columns
)
EnsembleKalmanProcesses.ParameterDistributions.ParameterDistributionMethod
ParameterDistribution(distribution_samples::AbstractMatrix,
                      constraint::Union{ConstraintType,AbstractVector{ConstraintType}},
                      name::AbstractString;
    params_are_columns::Bool = true)

constructor of a Samples ParameterDistribution from a matrix distribution_samples of parameters stored as columns by defaut, (array of) constraint, name.

EnsembleKalmanProcesses.ParameterDistributions.ParameterDistributionMethod
ParameterDistribution(distribution::ParameterDistributionType,
                               constraint::Union{ConstraintType,AbstractVector{ConstraintType}},
                               name::AbstractString)

constructor of a ParameterDistribution from a single distribution, (array of) constraint, name. these can used to build another ParameterDistribution

EnsembleKalmanProcesses.ParameterDistributions.ParameterDistributionMethod
ParameterDistribution(param_dist_dict::Union{Dict,AbstractVector})

Constructor taking in a Dict or array of Dicts. Each dict must contain the key-val pairs:

  • "distribution" - a distribution of ParameterDistributionType
  • "constraint" - constraint(s) given as a ConstraintType or array of ConstraintTypes with length equal to the dims of the distribution
  • "name" - a name of the distribution as a String.
EnsembleKalmanProcesses.ParameterDistributions.SamplesType
Samples{FT <: Real} <: ParameterDistributionType

A distribution comprised of only samples, stored as columns of parameters.

Fields

  • distribution_samples::AbstractMatrix{FT} where FT<:Real: Samples defining an empirical distribution, stored as columns
EnsembleKalmanProcesses.ParameterDistributions.VectorOfParameterizedType
VectorOfParameterized <: ParameterDistributionType

A distribution built from an array of Parametrized distributions. A utility to help stacking of distributions where a multivariate equivalent doesn't exist.

Fields

  • distribution::AbstractVector{DT} where DT<:Distributions.Distribution: A vector of parameterized distributions
Base.lengthMethod
length(c<:ConstraintType)

A constraint has length 1.

Base.ndimsMethod
ndims(grfi::GaussianRandomFieldInterface, function_parameter_opt = "dof")

Provides a relevant number of dimensions in different circumstances, If function_parameter_opt =

  • "dof" : returns n_dofs(grfi), the degrees of freedom in the function
  • "eval" : returns n_eval_pts(grfi), the number of discrete evaluation points of the function
  • "constraint": returns 1, the number of constraints in the evaluation space
Base.ndimsMethod
ndims(d<:ParametrizedDistributionType)

The number of dimensions of the parameter space

Base.sizeMethod
size(c<:ConstraintType)

A constraint has size 1.

Distributions.logpdfMethod
logpdf(pd::ParameterDistribution, xarray::Array{<:Real,1})

Obtains the independent logpdfs of the parameter distributions at xarray (non-Samples Distributions only), and returns their sum.

EnsembleKalmanProcesses.ParameterDistributions.batchMethod
batch(pd::ParameterDistribution; function_parameter_opt = "dof")

Returns a list of contiguous [collect(1:i), collect(i+1:j),... ] used to split parameter arrays by distribution dimensions. function_parameter_opt is passed to ndims in the special case of FunctionParameterDistributionTypes.

EnsembleKalmanProcesses.ParameterDistributions.boundedMethod
bounded(lower_bound::Real, upper_bound::Real)

Constructs a Constraint with provided upper and lower bounds, enforced by maps x -> log((x - lower_bound) / (upper_bound - x)) and x -> (upper_bound * exp(x) + lower_bound) / (exp(x) + 1).

EnsembleKalmanProcesses.ParameterDistributions.build_function_sampleMethod
build_function_sample(rng::AbstractRNG, grfi::GaussianRandomFieldInterface, n_draws::Int)

sample function distribution n_draw times on the discrete grid, from the stored prior distributions.

Defaults: n_draw = 1, rng = Random.GLOBAL_RNG, and coeff_vec sampled from the stored prior distribution

EnsembleKalmanProcesses.ParameterDistributions.constrained_gaussianMethod
constrained_gaussian(
    name::AbstractString,
    μ_c::Real,
    σ_c::Real,
    lower_bound::Real,
    upper_bound::Real;
    repeats = 1,
    optim_algorithm::Optim.AbstractOptimizer = NelderMead(),
    optim_kwargs...,
)

Constructor for a 1D ParameterDistribution consisting of a transformed Gaussian, constrained to have support on [lower_bound, upper_bound], with first two moments μ_c and σ_c^2. The moment integrals can't be done in closed form, so we set the parameters of the distribution with numerical optimization.

Note

The intended use case is defining priors set from user expertise for use in inference with adequate data, so for the sake of performance we only require that the optimization reproduce μ_c, σ_c to a loose tolerance (1e-5). Warnings are logged when the optimization fails.

Note

The distribution may be bimodal for σ_c large relative to the width of the bound interval. In extreme cases the distribution becomes concentrated at the bound endpoints. We regard this as a feature, not a bug, and do not warn the user when bimodality occurs.

EnsembleKalmanProcesses.ParameterDistributions.get_distributionMethod
get_distribution(pd::ParameterDistribution)

Returns a Dict of ParameterDistribution distributions, with the parameter names as dictionary keys. For parameters represented by Samples, the samples are returned as a 2D (parameter_dimension x n_samples) array.

EnsembleKalmanProcesses.ParameterDistributions.transform_constrained_to_unconstrainedMethod
transform_constrained_to_unconstrained(d::GaussianRandomFieldInterface, constraint::AbstractVector, x::AbstractMatrix)

Assume x is a matrix with columns as flattened samples of evaluation points. Remove the constraint from constraint to the output space of the function. Note this is the inverse of transform_unconstrained_to_constrained(...,build_flag=false)

EnsembleKalmanProcesses.ParameterDistributions.transform_constrained_to_unconstrainedMethod
transform_constrained_to_unconstrained(d::GaussianRandomFieldInterface, constraint::AbstractVector, x::AbstractVector)

Assume x is a flattened vector of evaluation points. Remove the constraint from constraint to the output space of the function. Note this is the inverse of transform_unconstrained_to_constrained(...,build_flag=false)

EnsembleKalmanProcesses.ParameterDistributions.transform_unconstrained_to_constrainedMethod
transform_unconstrained_to_constrained(d::GaussianRandomFieldInterface, constraint::AbstractVector, x::AbstractMatri)

Optional args: build_flag::Bool = true

Two functions, depending on build_flag If true, assume x is a matrix with columns of coefficient samples. Perform the following 3 maps.

  1. Apply the transformation to map (possibly constrained) parameter samples x into the unconstrained space. Using internally stored constraints (given by the coefficient prior)
  2. Build the unconstrained (flattened) function sample at the evaluation points from these constrained coefficients.
  3. Apply the constraint from constraint to the output space of the function.

If false, Assume x is a matrix with columns as flattened samples of evaluation points. Apply only step 3. above to x.

EnsembleKalmanProcesses.ParameterDistributions.transform_unconstrained_to_constrainedMethod
transform_unconstrained_to_constrained(d::GaussianRandomFieldInterface, constraint::AbstractVector, x::AbstractVector)

Optional Args build_flag::Bool = true

Two functions, depending on build_flag If true, assume x is a vector of coefficients. Perform the following 3 maps.

  1. Apply the transformation to map (possibly constrained) parameter samples x into the unconstrained space. Using internally stored constraints (given by the coefficient prior)
  2. Build the unconstrained (flattened) function sample at the evaluation points from these constrained coefficients.
  3. Apply the constraint from constraint to the output space of the function.

If false, Assume x is a flattened vector of evaluation points. Apply only step 3. above to x.

Statistics.covMethod
cov(pd::ParameterDistribution)

Returns a dense blocked (co)variance of the distributions.

Statistics.meanMethod
mean(pd::ParameterDistribution)

Returns a concatenated mean of the parameter distributions.

Statistics.varMethod
var(pd::ParameterDistribution)

Returns a flattened variance of the distributions

StatsBase.sampleMethod
sample([rng], pd::ParameterDistribution, [n_draws])

Draws n_draws samples from the parameter distributions pd. Returns an array, with parameters as columns. rng is optional and defaults to Random.GLOBAL_RNG. n_draws is optional and defaults to 1. Performed in computational space.

StatsBase.sampleMethod
sample([rng], d::Parameterized, [n_draws])

Draws n_draws samples from the parameter distributions d. Returns an array, with parameters as columns. rng is optional and defaults to Random.GLOBAL_RNG. n_draws is optional and defaults to 1. Performed in computational space.

StatsBase.sampleMethod
sample([rng], d::Samples, [n_draws])

Draws n_draws samples from the parameter distributions d. Returns an array, with parameters as columns. rng is optional and defaults to Random.GLOBAL_RNG. n_draws is optional and defaults to 1. Performed in computational space.

StatsBase.sampleMethod
sample([rng], d::VectorOfParameterized, [n_draws])

Draws n_draws samples from the parameter distributions d. Returns an array, with parameters as columns. rng is optional and defaults to Random.GLOBAL_RNG. n_draws is optional and defaults to 1. Performed in computational space.

Base.sizeMethod
size(dc::DataContainer, idx::IT) where {IT <: Integer}

Returns the size of the stored data. If idx provided, it returns the size along dimension idx.

Base.sizeMethod
size(pdc::PairedDataContainer, idx::IT) where {IT <: Integer}

Returns the sizes of the inputs and ouputs along dimension idx (if provided).

EnsembleKalmanProcesses.Observations.ObservationType
Observation{FT <: AbstractFloat}

Structure that contains the observations

Fields

  • samples::Array{Vector{FT}, 1} where FT<:AbstractFloat: vector of observational samples, each of length sample_dim

  • obs_noise_cov::Union{Nothing, LinearAlgebra.UniformScaling{FT}, AbstractMatrix{FT}, FT} where FT<:AbstractFloat: covariance of the observational noise (assumed to be normally distributed); sample_dim x sample_dim (where sample_dim is the number of elements in each sample), or a scalar if the sample dim is 1. If not supplied, obs_noise_cov is set to a diagonal matrix whose non-zero elements are the variances of the samples, or to a scalar variance in the case of 1d samples. obs_noise_cov is set to nothing if only a single sample is provided.

  • mean::Union{AbstractVector{FT}, FT} where FT<:AbstractFloat: sample mean

  • data_names::Union{String, AbstractVector{String}}: names of the data

EnsembleKalmanProcesses.TOMLInterface.assign_values!Method
assign_values!(
    member, 
    param_array, 
    param_distribution, 
    param_slices, 
    param_dict, 
    names)

Updates param_dict with the values of the given member of the param_array

Args: member - ensemble member (corresponds to column of param_array) param_array - Npar x Nens array of parameter values param_distribution - the parameter distribution underlying param_array param_slices - list of contiguous [collect(1:i), collect(i+1:j),... ] used to split parameter arrays by distribution dimensions param_dict - the dict of parameters to be updated with new parameter values names - array of parameter names

Returns the updated param_dict

EnsembleKalmanProcesses.TOMLInterface.collect_from_exprMethod
collect_from_expr(e, eltype; repeat=false)

Collects distributions or constraints

Args: e - expression containing the distribution or constraint information eltype - string denoting the type of elements that are collected, "d" for distributions, "c" for constraints repeat - true if this distribution or constraint is given as a repeat(...) expression, false otherwise

Returns an array of distributions / constraints, or a single distribution / constraint if only one is present

EnsembleKalmanProcesses.TOMLInterface.construct_constraintMethod
construct_constraint(param_info)

Extracts information on type and arguments of each constraint and uses that information to construct a Constraint.

Args: param_info - dictionary with (at least) a key "constraint", whose value is the parameter's constraint(s) (as parsed from TOML file)

Returns a Constraint

EnsembleKalmanProcesses.TOMLInterface.construct_priorMethod
construct_prior(param_info)

Extracts information on type and arguments of the prior distribution and use that information to construct an actual Distribution

Args: param_info - dictionary with (at least) a key "prior", whose value is the parameter's distribution(s) (as parsed from TOML file)

Returns a distribution of type Parameterized, Samples, or VectorOfParameterized

EnsembleKalmanProcesses.TOMLInterface.generate_subdir_namesMethod
generate_subdir_names(N; prefix="member", mode="all", pad_zeros=3)

Generates N directory names "<prefix>_<i>"; i=1, ..., N

Args: N - number of ensemble members (= number of subdirectories) or for mode=only, the chosen member prefix - prefix used for generation of subdirectory names mode - default =all generates all names, =only generates just the Nth name pad_zeros - amount of digits to pad to Returns a list of directory names

EnsembleKalmanProcesses.TOMLInterface.get_admissible_parametersMethod
get_admissible_parameters(param_dict)

Finds all parameters in param_dict that are admissible for calibration.

Args: param_dict - nested dictionary that has parameter names as keys and the corresponding dictionaries of parameter information as values

Returns an array of the names of all admissible parameters in param_dict. Admissible parameters must have a key "prior" and the value value of this is not set to "fixed". This allows for other parameters to be stored within the TOML file.

EnsembleKalmanProcesses.TOMLInterface.get_parameter_distributionMethod
get_parameter_distribution(param_dict, name)

Constructs a ParameterDistribution for a single parameter

Args: param_dict - nested dictionary that has parameter names as keys and the corresponding dictionary of parameter information (in particular, the parameters' prior distributions and constraints) as values name - parameter name

Returns a ParameterDistribution

EnsembleKalmanProcesses.TOMLInterface.get_parameter_distributionMethod
get_parameter_distribution(param_dict, names)

Constructs a ParameterDistribution for an array of parameters

Args: param_dict - nested dictionary that has parameter names as keys and the corresponding dictionary of parameter information (in particular, the parameters' prior distributions and constraints) as values names - array of parameter names

Returns a ParameterDistribution

EnsembleKalmanProcesses.TOMLInterface.get_parameter_valuesMethod
get_parameter_values(param_dict, names)

Gets parameter values from a parameter dictionary, indexing by name.

Args: param_dict - nested dictionary that has parameter names as keys and the corresponding dictionary of parameter information (in particular, the parameters' values) name - iterable parameter names return_type - return type, default "dict", otherwise "array"

EnsembleKalmanProcesses.TOMLInterface.get_regularizationMethod
get_regularization(param_dict, name)

Returns the regularization information for a single parameter

Args: param_dict - nested dictionary that has parameter names as keys and the corresponding dictionary of parameter information as values name - parameter name

Returns a tuple (<regularizationtype>, <regularizationvalue>), where the regularization type is either "L1" or "L2", and the regularization value is a float. Returns (nothing, nothing) if parameter has no regularization information.

EnsembleKalmanProcesses.TOMLInterface.get_regularizationMethod
get_regularization(param_dict, names)

Returns the regularization information for an array of parameters

Args: param_dict - nested dictionary that has parameter names as keys and the corresponding dictionary of parameter information as values names - array of parameter names

Returns an arary of tuples (<regularizationtype>, <regularizationvalue>), with the ith tuple corresponding to the parameter names[i]. The regularization type is either "L1" or "L2", and the regularization value is a float. Returns (nothing, nothing) for parameters that have no regularization information.

EnsembleKalmanProcesses.TOMLInterface.path_to_ensemble_memberMethod
path_to_ensemble_member(
    base_path,
    iteration,
    member,
    pad_zeros = 3,
)

Obtains the file path to a specified ensemble member. The likely form is base_path/iteration_X/member_Y/ with X,Y padded with zeros. The file path can be reconstructed with: base_path - base path to where EKP parameters are stored member - number of the ensemble member (without zero padding) iteration - iteration of ensemble method (if =nothing then only the load path is used) pad_zeros - amount of digits to pad to

EnsembleKalmanProcesses.TOMLInterface.save_parameter_ensembleMethod
save_parameter_ensemble(
    param_array,
    param_distribution,
    default_param_data,
    save_path,
    save_file,
    iteration
    pad_zeros=3,
apply_constraints=true
)

Saves the parameters in the given param_array to TOML files. The intended use is for saving the ensemble of parameters after each update of an ensemble Kalman process. Each ensemble member (column of param_array) is saved in a separate directory "member<j>" (j=1, ..., Nens). The name of the saved toml file is given by save_file; it is the same for all members. A directory "iteration<iter>" is created in `savepath`, which contains all the "member_<j>" subdirectories.

Args: param_array - array of size Nparam x Nens param_distribution - the parameter distribution underlying param_array apply_constraints - apply the constraints in param_distribution default_param_data - dict of default parameters to be combined and saved with the parameters in param_array into a toml file save_path - path to where the parameters will be saved save_file - name of the toml files to be generated iteration - the iteration of the ensemble Kalman process represented by the given param_array pad_zeros - the amount of zero-padding for the ensemble member number

EnsembleKalmanProcesses.TOMLInterface.write_log_fileMethod
write_log_file(param_dict, file_path)

Writes the parameters in param_dict into a .toml file

Args: param_dict - nested dictionary that has parameter names as keys and the corresponding dictionaries of parameter information as values file_path - path of the file where parameters are saved

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

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

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

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
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$

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
)
EnsembleKalmanProcess(
    params,
    obs_mean,
    obs_noise_cov,
    process;
    scheduler,
    accelerator,
    Δt,
    rng,
    failure_handler_method,
    localization_method,
    verbose
)
EnsembleKalmanProcess(
    obs_mean,
    obs_noise_cov,
    process;
    kwargs...
)
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)
FailureHandler(process, method)
FailureHandler(process, method)
FailureHandler(process, method)
FailureHandler(process, method)
FailureHandler(process, method)
FailureHandler(process, method)
FailureHandler(process, method)
FailureHandler(process, method)
EnsembleKalmanProcesses.FailureHandlerMethod
FailureHandler(process::Inversion, method::SampleSuccGauss)

Provides a failsafe update that

  • updates the successful ensemble according to the EKI update,
  • updates the failed ensemble by sampling from the updated successful ensemble.
EnsembleKalmanProcesses.FailureHandlerMethod
FailureHandler(process::SparseInversion, method::SampleSuccGauss)

Provides a failsafe update that

  • updates the successful ensemble according to the SparseEKI update,
  • updates the failed ensemble by sampling from the updated successful ensemble.
EnsembleKalmanProcesses.FailureHandlerMethod
FailureHandler(process::TransformInversion, method::SampleSuccGauss)

Provides a failsafe update that

  • updates the successful ensemble according to the ETKI update,
  • updates the failed ensemble by sampling from the updated successful ensemble.
EnsembleKalmanProcesses.FailureHandlerMethod
FailureHandler(process::Unscented, method::SampleSuccGauss)

Provides a failsafe update that

  • computes all means and covariances over the successful sigma points,
  • rescales the mean weights and the off-center covariance weights of the successful particles to sum to the same value as the original weight sums.
EnsembleKalmanProcesses.FirstOrderNesterovAcceleratorType
mutable struct FirstOrderNesterovAccelerator{FT<:AbstractFloat} <: EnsembleKalmanProcesses.Accelerator

Accelerator that adapts a variant of Nesterov's momentum method for EKI. In this variant the momentum parameter at iteration $k$, is given by $1-\frac{3}{k+2}$. This variant is a first order repesentation of the desired asymptotic behavior $1-\frac{3}{k}- \mathcal{O}(\frac{1}{k^2})$ advocated in (Su, Boyd, Candes, 2014). Stores a previous state value u_prev for computational purposes (note this is distinct from state returned as "ensemble value")

  • r::AbstractFloat

  • u_prev::Array{FT} where FT<:AbstractFloat

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
EnsembleKalmanProcesses.NesterovAcceleratorType
mutable struct NesterovAccelerator{FT<:AbstractFloat} <: EnsembleKalmanProcesses.Accelerator

Accelerator that adapts a variant of Nesterov's momentum method for EKI. In this variant the momentum parameter is given by $\theta_{k+1}(\frac{1}{\theta_{k}} - 1)$ where $\theta_{k+1} = \frac{-\theta_k^2 + \sqrt{\theta_k^4 + 4 \theta_k^2}}{2}$ and $\theta_0 = 1$. This recurrence, also mentioned in (Su, Boyd, Candes, 2014), is our recommended variant. Stores a previous state value u_prev for computational purposes (note this is distinct from state returned as "ensemble value")

  • u_prev::Array{FT} where FT<:AbstractFloat

  • θ_prev::AbstractFloat

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.

EnsembleKalmanProcesses.SamplerType
Sampler{FT<:AbstractFloat,IT<:Int} <: Process

An ensemble Kalman Sampler process.

Fields

  • prior_mean::Vector{FT} where FT<:AbstractFloat: Mean of Gaussian parameter prior in unconstrained space

  • prior_cov::Union{LinearAlgebra.UniformScaling{FT}, AbstractMatrix{FT}} where FT<:AbstractFloat: Covariance of Gaussian parameter prior in unconstrained space

Constructors

Sampler(prior_mean, prior_cov)
Sampler(prior)
EnsembleKalmanProcesses.SparseInversionType
SparseInversion <: Process

A sparse ensemble Kalman Inversion process

Fields

  • γ::AbstractFloat: upper limit of l1-norm

  • threshold_value::AbstractFloat: threshold below which the norm of parameters is pruned to zero

  • uc_idx::Union{Colon, AbstractVector}: indices of parameters included in the evaluation of l1-norm constraint

  • reg::AbstractFloat: a small regularization value to enhance robustness of convex optimization

Constructors

SparseInversion(γ, threshold_value, uc_idx, reg)
SparseInversion(γ; threshold_value, uc_idx, reg)
SparseInversion(; γ, threshold_value, uc_idx, reg)
EnsembleKalmanProcesses.TransformInversionType
TransformInversion <: Process

An ensemble transform Kalman inversion process.

Fields

  • Γ_inv::Union{LinearAlgebra.UniformScaling{FT}, AbstractMatrix{FT}} where FT<:AbstractFloat: Inverse of the observation error covariance matrix
EnsembleKalmanProcesses.UnscentedType
Unscented{FT<:AbstractFloat, IT<:Int} <: Process

An unscented Kalman Inversion process.

Fields

  • u_mean::Any: an interable of arrays of size N_parameters containing the mean of the parameters (in each uki iteration a new array of mean is added), note - this is not the same as the ensemble mean of the sigma ensemble as it is taken prior to prediction

  • uu_cov::Any: an iterable of arrays of size (N_parameters x N_parameters) containing the covariance of the parameters (in each uki iteration a new array of cov is added), note - this is not the same as the ensemble cov of the sigma ensemble as it is taken prior to prediction

  • obs_pred::Any: an iterable of arrays of size N_y containing the predicted observation (in each uki iteration a new array of predicted observation is added)

  • c_weights::AbstractVecOrMat{FT} where FT<:AbstractFloat: weights in UKI

  • mean_weights::AbstractVector{FT} where FT<:AbstractFloat

  • cov_weights::AbstractVector{FT} where FT<:AbstractFloat

  • N_ens::Int64: number of particles 2N+1 or N+2

  • Σ_ω::AbstractMatrix{FT} where FT<:AbstractFloat: covariance of the artificial evolution error

  • Σ_ν_scale::AbstractFloat: covariance of the artificial observation error

  • α_reg::AbstractFloat: regularization parameter

  • r::AbstractVector{FT} where FT<:AbstractFloat: regularization vector

  • update_freq::Int64: update frequency

  • impose_prior::Bool: using augmented system (Tikhonov regularization with Kalman inversion in Chada et al 2020 and Huang et al (2022)) to regularize the inverse problem, which also imposes prior for posterior estimation.

  • prior_mean::Any: prior mean - defaults to initial mean

  • prior_cov::Any: prior covariance - defaults to initial covariance

  • iter::Int64: current iteration number

Constructors

Unscented(
    u0_mean::AbstractVector{FT},
    uu0_cov::AbstractMatrix{FT};
    α_reg::FT = 1.0,
    update_freq::IT = 0,
    modified_unscented_transform::Bool = true,
    impose_prior::Bool = false,
    prior_mean::Any,
    prior_cov::Any,
    sigma_points::String = symmetric
) where {FT <: AbstractFloat, IT <: Int}

Construct an Unscented Inversion Process.

Inputs:

  • u0_mean: Mean at initialization.
  • uu0_cov: Covariance at initialization.
  • α_reg: Hyperparameter controlling regularization toward the prior mean (0 < α_reg ≤ 1),

default should be 1, without regulariazion.

  • update_freq: Set to 0 when the inverse problem is not identifiable,

namely the inverse problem has multiple solutions, the covariance matrix will represent only the sensitivity of the parameters, instead of posterior covariance information; set to 1 (or anything > 0) when the inverse problem is identifiable, and the covariance matrix will converge to a good approximation of the posterior covariance with an uninformative prior.

  • modified_unscented_transform: Modification of the UKI quadrature given in Huang et al (2021).
  • impose_prior: using augmented system (Tikhonov regularization with Kalman inversion in Chada et al 2020 and Huang et al (2022)) to regularize the inverse problem, which also imposes prior for posterior estimation. If impose_prior == true, prior mean and prior cov must be provided. This is recommended to use, especially when the number of observations is smaller than the number of parameters (ill-posed inverse problems). When this is used, other regularizations are turned off automatically.
  • prior_mean: Prior mean used for regularization.
  • prior_cov: Prior cov used for regularization.
  • sigma_points: String of sigma point type, it can be symmetric with 2N_par+1 ensemble members or simplex with N_par+2 ensemble members.
Unscented(
    u_mean,
    uu_cov,
    obs_pred,
    c_weights,
    mean_weights,
    cov_weights,
    N_ens,
    Σ_ω,
    Σ_ν_scale,
    α_reg,
    r,
    update_freq,
    impose_prior,
    prior_mean,
    prior_cov,
    iter
)
Unscented(
    u0_mean,
    uu0_cov;
    α_reg,
    update_freq,
    modified_unscented_transform,
    impose_prior,
    prior_mean,
    prior_cov,
    sigma_points
)
Unscented(prior; kwargs...)
EnsembleKalmanProcesses.accelerate!Method

State update method for UKI with constant Nesterov Accelerator. Performs identical update as with other methods, but requires reconstruction of mean and covariance of the accelerated positions prior to saving.

EnsembleKalmanProcesses.accelerate!Method

State update method for UKI with first-order Nesterov Accelerator. Performs identical update as with other methods, but requires reconstruction of mean and covariance of the accelerated positions prior to saving.

EnsembleKalmanProcesses.accelerate!Method

State update method for UKI with Nesterov Accelerator. The dependence of the momentum parameter can be found e.g. here "https://www.damtp.cam.ac.uk/user/hf323/M19-OPT/lecture5.pdf" Performs identical update as with other methods, but requires reconstruction of mean and covariance of the accelerated positions prior to saving.

EnsembleKalmanProcesses.accelerate!Method

Performs state update with Nesterov momentum approach. The dependence of the momentum parameter can be found e.g. here "https://www.damtp.cam.ac.uk/user/hf323/M19-OPT/lecture5.pdf"

EnsembleKalmanProcesses.additive_inflation!Method
additive_inflation!(
    ekp::EnsembleKalmanProcess;
    use_prior_cov::Bool = false,
    s::FT = 1.0,
) where {FT <: Real}

Applies additive Gaussian noise to particles. Noise is drawn from normal distribution with 0 mean and scaled parameter covariance. If usepriorcov=false (default), scales parameter covariance matrix from current ekp iteration. Otherwise, scales parameter covariance of initial ensemble. Inputs: - ekp :: The EnsembleKalmanProcess to update. - s :: Scaling factor for time step in additive perturbation. - usepriorcov :: Bool specifying whether to use prior covariance estimate for additive inflation. If false (default), parameter covariance from the current iteration is used.

EnsembleKalmanProcesses.calculate_timestep!Method
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.compute_error!Method
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.

EnsembleKalmanProcesses.construct_covMethod
construct_cov(
    uki::EnsembleKalmanProcess{FT, IT, Unscented},
    x::AbstractMatrix{FT},
    x_mean::AbstractVector{FT},
    obs_mean::AbstractMatrix{FT},
    y_mean::AbstractVector{FT};
    cov_weights = uki.process.cov_weights,
) where {FT <: AbstractFloat, IT <: Int, P <: Process}

Constructs covariance xy_cov from ensemble x and mean x_mean, ensemble obs_mean and mean y_mean.

EnsembleKalmanProcesses.construct_covMethod
construct_cov(
    uki::EnsembleKalmanProcess{FT, IT, Unscented},
    x::AbstractVecOrMat{FT},
    x_mean::Union{FT, AbstractVector{FT}, Nothing} = nothing;
    cov_weights = uki.process.cov_weights,
) where {FT <: AbstractFloat, IT <: Int}

Constructs covariance xx_cov from ensemble x and mean x_mean.

EnsembleKalmanProcesses.construct_initial_ensembleMethod
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.

EnsembleKalmanProcesses.construct_meanMethod
construct_mean(
    uki::EnsembleKalmanProcess{FT, IT, Unscented},
    x::AbstractVecOrMat{FT};
    mean_weights = uki.process.mean_weights,
) where {FT <: AbstractFloat, IT <: Int}

constructs mean x_mean from an ensemble x.

EnsembleKalmanProcesses.construct_sigma_ensembleMethod
construct_sigma_ensemble(
    process::Unscented,
    x_mean::Array{FT},
    x_cov::AbstractMatrix{FT},
) where {FT <: AbstractFloat, IT <: Int}

Construct the sigma ensemble based on the mean x_mean and covariance x_cov.

EnsembleKalmanProcesses.construct_successful_covMethod
construct_successful_cov(
    uki::EnsembleKalmanProcess{FT, IT, Unscented},
    x::AbstractMatrix{FT},
    x_mean::AbstractArray{FT},
    obs_mean::AbstractMatrix{FT},
    y_mean::AbstractArray{FT},
    successful_indices::Union{AbstractVector{IT}, AbstractVector{Any}},
) where {FT <: AbstractFloat, IT <: Int}

Constructs covariance of x and obs_mean - y_mean over successful particles by rescaling the off-center weights over the successful off-center particles.

EnsembleKalmanProcesses.construct_successful_covMethod
construct_successful_cov(
    uki::EnsembleKalmanProcess{FT, IT, Unscented},
    x::AbstractVecOrMat{FT},
    x_mean::Union{AbstractVector{FT}, FT},
    successful_indices::Union{AbstractVector{IT}, AbstractVector{Any}},
) where {FT <: AbstractFloat, IT <: Int}

Constructs variance of x over successful particles by rescaling the off-center weights over the successful off-center particles.

EnsembleKalmanProcesses.construct_successful_meanMethod
construct_successful_mean(
    uki::EnsembleKalmanProcess{FT, IT, Unscented},
    x::AbstractVecOrMat{FT},
    successful_indices::Union{AbstractVector{IT}, AbstractVector{Any}},
) where {FT <: AbstractFloat, IT <: Int}

Constructs mean over successful particles by rescaling the quadrature weights over the successful particles. If the central particle fails in a modified unscented transform, the mean is computed as the ensemble mean over all successful particles.

EnsembleKalmanProcesses.eki_updateMethod
 eki_update(
    ekp::EnsembleKalmanProcess{FT, IT, Inversion},
    u::AbstractMatrix{FT},
    g::AbstractMatrix{FT},
    y::AbstractMatrix{FT},
    obs_noise_cov::Union{AbstractMatrix{CT}, UniformScaling{CT}},
) where {FT <: Real, IT, CT <: Real}

Returns the updated parameter vectors given their current values and the corresponding forward model evaluations, using the inversion algorithm from eqns. (4) and (5) of Schillings and Stuart (2017).

Localization is implemented following the ekp.localizer.

EnsembleKalmanProcesses.eks_updateMethod
 eks_update(
    ekp::EnsembleKalmanProcess{FT, IT, Sampler{FT}},
    u::AbstractMatrix{FT},
    g::AbstractMatrix{FT},
) where {FT <: Real, IT}

Returns the updated parameter vectors given their current values and the corresponding forward model evaluations, using the sampler algorithm.

The current implementation assumes that rows of u and g correspond to ensemble members, so it requires passing the transpose of the u and g arrays associated with ekp.

EnsembleKalmanProcesses.etki_updateMethod
 etki_update(
    ekp::EnsembleKalmanProcess{FT, IT, TransformInversion},
    u::AbstractMatrix{FT},
    g::AbstractMatrix{FT},
    y::AbstractVector{FT},
    obs_noise_cov::Union{AbstractMatrix{CT}, UniformScaling{CT}},
) where {FT <: Real, IT, CT <: Real}

Returns the updated parameter vectors given their current values and the corresponding forward model evaluations.

EnsembleKalmanProcesses.get_cov_blocksMethod
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.

EnsembleKalmanProcesses.get_gMethod
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_gMethod
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.

EnsembleKalmanProcesses.get_g_finalMethod
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_g_meanMethod
get_g_mean(ekp::EnsembleKalmanProcess, iteration::IT) where {IT <: Integer}

Returns the mean forward map evaluation at the given iteration.

EnsembleKalmanProcesses.get_uMethod
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_uMethod
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.

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

Returns the unconstrained parameter sample covariance at the given iteration.

EnsembleKalmanProcesses.get_u_covMethod
get_u_cov(uki::EnsembleKalmanProcess{FT, IT, Unscented}, iteration::IT)

Returns the unconstrained parameter covariance at the requested iteration.

EnsembleKalmanProcesses.get_u_finalMethod
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_u_meanMethod
get_u_mean(ekp::EnsembleKalmanProcess, iteration::IT) where {IT <: Integer}

Returns the mean unconstrained parameter at the given iteration.

EnsembleKalmanProcesses.get_u_meanMethod
get_u_mean(uki::EnsembleKalmanProcess{FT, IT, Unscented}, iteration::IT)

Returns the mean unconstrained parameter at the requested iteration.

EnsembleKalmanProcesses.get_u_priorMethod
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_ϕMethod
get_ϕ(prior::ParameterDistribution, ekp::EnsembleKalmanProcess)

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

EnsembleKalmanProcesses.get_ϕMethod
get_ϕ(prior::ParameterDistribution, ekp::EnsembleKalmanProcess, iteration::IT)

Returns the constrained parameters at the given iteration.

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

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

EnsembleKalmanProcesses.get_ϕ_mean_finalMethod
get_ϕ_mean_final(prior::ParameterDistribution, ekp::EnsembleKalmanProcess)

Get the constrained transform of the mean unconstrained parameter at the last iteration.

EnsembleKalmanProcesses.multiplicative_inflation!Method
multiplicative_inflation!(
    ekp::EnsembleKalmanProcess;
    s,
) where {FT, IT}

Applies multiplicative noise to particles. Inputs: - ekp :: The EnsembleKalmanProcess to update. - s :: Scaling factor for time step in multiplicative perturbation.

EnsembleKalmanProcesses.posdef_correctMethod
posdef_correct(mat::AbstractMatrix; tol) -> Any

Makes square matrix mat positive definite, by symmetrizing and bounding the minimum eigenvalue below by tol

EnsembleKalmanProcesses.sample_empirical_gaussianMethod
sample_empirical_gaussian(
    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.sparse_eki_updateMethod
 sparse_eki_update(
    ekp::EnsembleKalmanProcess{FT, IT, SparseInversion{FT}},
    u::AbstractMatrix{FT},
    g::AbstractMatrix{FT},
    y::AbstractMatrix{FT},
    obs_noise_cov::Union{AbstractMatrix{CT}, UniformScaling{CT}},
) where {FT <: Real, CT <: Real, IT}

Returns the sparse updated parameter vectors given their current values and the corresponding forward model evaluations, using the inversion algorithm from eqns. (3.7) and (3.14) of Schneider et al. (2021).

Localization is applied following Tong and Morzfeld (2022).

EnsembleKalmanProcesses.sparse_qpMethod
sparse_qp(
    ekp::EnsembleKalmanProcess{FT, IT, SparseInversion{FT}},
    v_j::Vector{FT},
    cov_vv_inv::AbstractMatrix{FT},
    H_u::SparseArrays.SparseMatrixCSC{FT},
    H_g::SparseArrays.SparseMatrixCSC{FT},
    y_j::Vector{FT};
    H_uc::SparseArrays.SparseMatrixCSC{FT} = H_u,
) where {FT, IT}

Solving quadratic programming problem with sparsity constraint.

EnsembleKalmanProcesses.split_indices_by_successMethod
 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.update_ensemble!Method
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.
EnsembleKalmanProcesses.update_ensemble!Method
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.
EnsembleKalmanProcesses.update_ensemble!Method
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.
EnsembleKalmanProcesses.update_ensemble!Method
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.
EnsembleKalmanProcesses.update_ensemble!Method
update_ensemble!(
    ekp::EnsembleKalmanProcess,
    g::AbstractMatrix{FT};
    multiplicative_inflation::Bool = false,
    additive_inflation::Bool = false,
    use_prior_cov::Bool = false,
    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.
  • usepriorcov :: Bool specifying whether to use prior covariance estimate for additive inflation. 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!.
EnsembleKalmanProcesses.update_ensemble!Method
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.
EnsembleKalmanProcesses.update_ensemble_analysis!Method
update_ensemble_analysis!(
    uki::EnsembleKalmanProcess{FT, IT, Unscented},
    u_p::AbstractMatrix{FT},
    g::AbstractMatrix{FT},
) where {FT <: AbstractFloat, IT <: Int}

UKI analysis step : g is the predicted observations Ny x N_ens matrix

EnsembleKalmanProcesses.Localizers.BernoulliDropoutType
BernoulliDropout{FT <: Real} <: LocalizationMethod

Localization method that drops cross-covariance terms with probability $1-p$, retaining a Hermitian structure.

Fields

  • prob::Real: Probability of keeping a given cross-covariance term
EnsembleKalmanProcesses.Localizers.LocalizerType
Localizer{LM <: LocalizationMethod, T}

Structure that defines a localize function, based on a localization method.

Fields

  • localize::Function: Localizing function of the form: cov -> kernel .* cov

Constructors

Localizer(localization, p, d, J)
Localizer(localization, p, d, J, T)
Localizer(localization, p, d, J)
Localizer(localization, p, d, J, T)
Localizer(localization, p, d, J)
Localizer(localization, p, d, J, T)
Localizer(localization, p, d, J)
Localizer(localization, p, d, J, T)
Localizer(localization, p, d, J)
Localizer(localization, p, d, J, T)
Localizer(localization, p, d, J)
Localizer(localization, p, d, J, T)
EnsembleKalmanProcesses.Localizers.RBFType
RBF{FT <: Real} <: LocalizationMethod

Radial basis function localization method. Covariance terms $C_{i,j}$ are damped through multiplication with a centered Gaussian with standardized deviation $d(i,j)= \vert i-j \vert / l$.

Fields

  • lengthscale::Real: Length scale defining the RBF kernel
EnsembleKalmanProcesses.Localizers.SECType
SEC{FT <: Real} <: LocalizationMethod

Sampling error correction that shrinks correlations by a factor of $\vert r \vert ^\alpha$, as per Lee (2021). Sparsity of the resulting correlations can be imposed through the parameter r_0.

Lee, Y. (2021). Sampling error correction in ensemble Kalman inversion. arXiv:2105.11341 [cs, math]. http://arxiv.org/abs/2105.11341

Fields

  • α::Real: Controls degree of sampling error correction

  • r_0::Real: Cutoff correlation

EnsembleKalmanProcesses.Localizers.SECFisherType
SECFisher <: LocalizationMethod

Sampling error correction for EKI, as per Lee (2021), but using the method from Flowerdew (2015) based on the Fisher transformation. Correlations are shrinked by a factor determined by the sample correlation and the ensemble size.

Flowerdew, J. (2015). Towards a theory of optimal localisation. Tellus A: Dynamic Meteorology and Oceanography, 67(1), 25257. https://doi.org/10.3402/tellusa.v67.25257

Lee, Y. (2021). Sampling error correction in ensemble Kalman inversion. arXiv:2105.11341 [cs, math]. http://arxiv.org/abs/2105.11341