EnsembleKalmanProcesses.ParameterDistributions.Constraint
— TypeConstraint{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 domainbounds::Union{Nothing, Dict}
: Dictionary of values used to build the Constraint (e.g. "lowerbound" or "upperbound")
EnsembleKalmanProcesses.ParameterDistributions.GaussianRandomFieldInterface
— Typestruct 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 functionpackage::EnsembleKalmanProcesses.ParameterDistributions.GaussianRandomFieldsPackage
: the choice of GRF packagedistribution::EnsembleKalmanProcesses.ParameterDistributions.ParameterDistribution
: the distribution of the coefficients that we shall compute with
EnsembleKalmanProcesses.ParameterDistributions.GaussianRandomFieldInterface
— MethodGaussianRandomFieldInterface(gaussian_random_field::Any, package::GRFJL)
Constructor of the interface with GRFJL package. Internally this constructs a prior for the degrees of freedom of the function distribution
EnsembleKalmanProcesses.ParameterDistributions.GaussianRandomFieldsPackage
— Typeabstract type GaussianRandomFieldsPackage
Type to dispatch which Gaussian Random Field package to use:
GRFJL
uses the Julia PackageGaussianRandomFields.jl
EnsembleKalmanProcesses.ParameterDistributions.ParameterDistribution
— TypeParameterDistribution
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 spaceconstraint::AbstractVector{CType} where CType<:EnsembleKalmanProcesses.ParameterDistributions.ConstraintType
: Vector of constraints defining transformations between constrained and unconstrained spacename::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.ParameterDistribution
— MethodParameterDistribution(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.ParameterDistribution
— MethodParameterDistribution(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.ParameterDistribution
— MethodParameterDistribution(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 ofParameterDistributionType
"constraint"
- constraint(s) given as aConstraintType
or array ofConstraintType
s with length equal to the dims of the distribution"name"
- a name of the distribution as a String.
EnsembleKalmanProcesses.ParameterDistributions.Parameterized
— TypeParameterized <: ParameterDistributionType
A distribution constructed from a parameterized formula (e.g Julia Distributions.jl)
Fields
distribution::Distributions.Distribution
: A parameterized distribution
EnsembleKalmanProcesses.ParameterDistributions.Samples
— TypeSamples{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.VectorOfParameterized
— TypeVectorOfParameterized <: 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.length
— Methodlength(c<:ConstraintType)
A constraint has length 1.
Base.ndims
— Methodndims(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.ndims
— Methodndims(d<:ParametrizedDistributionType)
The number of dimensions of the parameter space
Base.size
— Methodsize(c<:ConstraintType)
A constraint has size 1.
Distributions.logpdf
— Methodlogpdf(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.batch
— Methodbatch(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 FunctionParameterDistributionType
s.
EnsembleKalmanProcesses.ParameterDistributions.bounded
— Methodbounded(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.bounded_above
— Methodbounded_above(upper_bound::FT) where {FT <: Real}
Constructs a Constraint with provided upper bound, enforced by maps x -> log(upper_bound - x)
and x -> upper_bound - exp(x)
.
EnsembleKalmanProcesses.ParameterDistributions.bounded_below
— Methodbounded_below(lower_bound::FT) where {FT <: Real}
Constructs a Constraint with provided lower bound, enforced by maps x -> log(x - lower_bound)
and x -> exp(x) + lower_bound
.
EnsembleKalmanProcesses.ParameterDistributions.build_function_sample
— Methodbuild_function_sample(grfi::GaussianRandomFieldInterface, coeff_vecormat::AbstractVecOrMat, n_draws::Int)
build function n_draw
times on the discrete grid, given the coefficients coeff_vecormat
.
Defaults: n_draw = 1
.
EnsembleKalmanProcesses.ParameterDistributions.build_function_sample
— Methodbuild_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.combine_distributions
— Methodcombine_distributions(pd_vec::AbstractVector{PD})
Form a ParameterDistribution by concatenating a vector of single ParameterDistributions.
EnsembleKalmanProcesses.ParameterDistributions.constrained_gaussian
— Methodconstrained_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.
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.
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.eval_pts
— Methodeval_pts(grfi::GaussianRandomFieldInterface)
the discrete evaluation point grid, stored as a range in each dimension
EnsembleKalmanProcesses.ParameterDistributions.get_all_constraints
— Methodget_all_constraints(grfi::GaussianRandomFieldInterface) = get_all_constraints(get_distribution(grfi))
gets all the constraints of the internally stored coefficient prior distribution of the GRFI
EnsembleKalmanProcesses.ParameterDistributions.get_all_constraints
— Methodget_all_constraints(pd::ParameterDistribution; return_dict = false)
Returns the (flattened) array of constraints of the parameter distribution. or as a dictionary ("param_name" => constraints)
EnsembleKalmanProcesses.ParameterDistributions.get_bounds
— Methodget_bounds(c::Constraint)
Gets the bounds field from the Constraint.
EnsembleKalmanProcesses.ParameterDistributions.get_constraint_type
— Methodget_bounds(c::Constraint{T})
Gets the parametric type T.
EnsembleKalmanProcesses.ParameterDistributions.get_dimensions
— Methodget_dimensions(pd::ParameterDistribution; function_parameter_opt = "dof")
The number of dimensions of the parameter space. (Also represents other dimensions of interest for FunctionParameterDistributionType
s with keyword argument)
EnsembleKalmanProcesses.ParameterDistributions.get_distribution
— Methodgets the, distribution over the coefficients
EnsembleKalmanProcesses.ParameterDistributions.get_distribution
— Methodget_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.get_grf
— Methodgets the distribution, i.e. Gaussian random field object
EnsembleKalmanProcesses.ParameterDistributions.get_n_samples
— Methodget_n_samples(pd::ParameterDistribution)
The number of samples in a Samples distribution
EnsembleKalmanProcesses.ParameterDistributions.get_name
— Methodget_name(pd::ParameterDistribution)
Returns a list of ParameterDistribution names.
EnsembleKalmanProcesses.ParameterDistributions.get_package
— Methodgets the package type used to construct the GRF
EnsembleKalmanProcesses.ParameterDistributions.input_dims
— Methodinput_dims(grfi::GaussianRandomFieldInterface
the number of input dimensions of the GRF
EnsembleKalmanProcesses.ParameterDistributions.n_dofs
— Methodn_dofs(grfi::GaussianRandomFieldInterface)
the number of degrees of freedom / coefficients (i.e. the number of parameters)
EnsembleKalmanProcesses.ParameterDistributions.n_eval_pts
— Methodn_eval_pts(grfi::GaussianRandomFieldInterface)
the number of total discrete evaluation points
EnsembleKalmanProcesses.ParameterDistributions.n_samples
— Methodn_samples(d<:Samples)
The number of samples in the array.
EnsembleKalmanProcesses.ParameterDistributions.no_constraint
— Methodno_constraint()
Constructs a Constraint with no constraints, enforced by maps x -> x and x -> x.
EnsembleKalmanProcesses.ParameterDistributions.spectrum
— Methodspectrum(grfi::GaussianRandomFieldInterface)
the spectral information of the GRF, e.g. the Karhunen-Loeve coefficients and eigenfunctions if using this decomposition
EnsembleKalmanProcesses.ParameterDistributions.transform_constrained_to_unconstrained
— Methodtransform_constrained_to_unconstrained(pd::ParameterDistribution, x::Array{Array{<:Real,2},1})
Apply the transformation to map parameter sample ensembles x
from the (possibly) constrained space into unconstrained space. Here, x
is an iterable of parameters sample ensembles for different EKP iterations.
EnsembleKalmanProcesses.ParameterDistributions.transform_constrained_to_unconstrained
— Methodtransform_constrained_to_unconstrained(d::ParameterDistribution, x::Dict)
Apply the transformation to map (possibly constrained) parameter samples x
into the unconstrained space. Here, x
contains parameter names as keys, and 1- or 2-arrays as parameter samples.
EnsembleKalmanProcesses.ParameterDistributions.transform_constrained_to_unconstrained
— Methodtransform_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_unconstrained
— Methodtransform_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_constrained_to_unconstrained
— Methodtransform_constrained_to_unconstrained(pd::ParameterDistribution, x::VecOrMat)
Apply the transformation to map (possibly constrained) parameter sample(s) x
into the unconstrained space.
Each column of x
is a sample, and each row is a parameter.
The return type is a vector if x
is a vector, and a matrix otherwise.
EnsembleKalmanProcesses.ParameterDistributions.transform_unconstrained_to_constrained
— Methodtransform_unconstrained_to_constrained(pd::ParameterDistribution, x::Array{Array{<:Real,2},1})
Apply the transformation to map parameter sample ensembles x
from the unconstrained space into (possibly constrained) space. Here, x
is an iterable of parameters sample ensembles for different EKP iterations.
EnsembleKalmanProcesses.ParameterDistributions.transform_unconstrained_to_constrained
— Methodtransform_unconstrained_to_constrained(d::ParameterDistribution, x::Dict)
Apply the transformation to map (possibly constrained) parameter samples x
into the unconstrained space. Here, x
contains parameter names as keys, and 1- or 2-arrays as parameter samples.
EnsembleKalmanProcesses.ParameterDistributions.transform_unconstrained_to_constrained
— Methodtransform_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.
- Apply the transformation to map (possibly constrained) parameter samples
x
into the unconstrained space. Using internally stored constraints (given by the coefficient prior) - Build the unconstrained (flattened) function sample at the evaluation points from these constrained coefficients.
- 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_constrained
— Methodtransform_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.
- Apply the transformation to map (possibly constrained) parameter samples
x
into the unconstrained space. Using internally stored constraints (given by the coefficient prior) - Build the unconstrained (flattened) function sample at the evaluation points from these constrained coefficients.
- 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.
EnsembleKalmanProcesses.ParameterDistributions.transform_unconstrained_to_constrained
— Methodtransform_unconstrained_to_constrained(pd::ParameterDistribution, x::VecOrMat)
Apply the transformation to map unconstrained parameter sample(s) x
into the constrained space.
Each column of x
is a sample, and each row is a parameter.
The return type is a vector if x
is a vector, and a matrix otherwise.
Statistics.cov
— Methodcov(pd::ParameterDistribution)
Returns a dense blocked (co)variance of the distributions.
Statistics.mean
— Methodmean(pd::ParameterDistribution)
Returns a concatenated mean of the parameter distributions.
Statistics.var
— Methodvar(pd::ParameterDistribution)
Returns a flattened variance of the distributions
StatsBase.sample
— Methodsample([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.sample
— Methodsample([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.sample
— Methodsample([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.sample
— Methodsample([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.
EnsembleKalmanProcesses.DataContainers.DataContainer
— TypeDataContainer{FT <: Real}
Container to store data samples as columns in an array.
EnsembleKalmanProcesses.DataContainers.PairedDataContainer
— TypePairedDataContainer{FT <: Real}
Stores input - output pairs as data containers, there must be an equal number of inputs and outputs.
Base.size
— Methodsize(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.size
— Methodsize(pdc::PairedDataContainer, idx::IT) where {IT <: Integer}
Returns the sizes of the inputs and ouputs along dimension idx
(if provided).
EnsembleKalmanProcesses.DataContainers.get_data
— Methodget_data(pdc::PairedDataContainer)
Returns both input and output data stored in pdc
.
EnsembleKalmanProcesses.DataContainers.get_inputs
— Methodget_inputs(pdc::PairedDataContainer)
Returns input data stored in pdc
.
EnsembleKalmanProcesses.DataContainers.get_outputs
— Methodget_outputs(pdc::PairedDataContainer)
Returns output data stored in pdc
.
EnsembleKalmanProcesses.Observations.Observation
— TypeObservation{FT <: AbstractFloat}
Structure that contains the observations
Fields
samples::Array{Vector{FT}, 1} where FT<:AbstractFloat
: vector of observational samples, each of length sample_dimobs_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
(wheresample_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 meandata_names::Union{String, AbstractVector{String}}
: names of the data
EnsembleKalmanProcesses.TOMLInterface.assign_values!
— Methodassign_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_expr
— Methodcollect_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_2d_array
— Methodconstruct_2d_array(arr)
Reconstructs 2d array of samples
Args: arr
- expression (has type Expr
) with head vcat
.
Returns a 2d array of samples constructed from the arguments of expr
EnsembleKalmanProcesses.TOMLInterface.construct_constraint
— Methodconstruct_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_prior
— Methodconstruct_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_names
— Methodgenerate_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 N
th name pad_zeros
- amount of digits to pad to Returns a list of directory names
EnsembleKalmanProcesses.TOMLInterface.get_admissible_parameters
— Methodget_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_distribution_from_expr
— Methodget_distribution_from_expr(d)
Parses a distribution
Args: d
- expression containing the distribution information
Returns a distribution of type Parameterized
or Samples
EnsembleKalmanProcesses.TOMLInterface.get_parameter_distribution
— Methodget_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_distribution
— Methodget_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_values
— Methodget_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_regularization
— Methodget_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_regularization
— Methodget_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.get_vector_of_parameterized
— Methodget_vector_of_parameterized(d)
Parses a distribution of type VectorOfParameterized
Args: d
- expression containing the distribution information
Returns a VectorOfParameterized
EnsembleKalmanProcesses.TOMLInterface.path_to_ensemble_member
— Methodpath_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.path_to_ensemble_member
— MethodOne can also call this without the iteration level
EnsembleKalmanProcesses.TOMLInterface.save_parameter_ensemble
— Methodsave_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.save_parameter_ensemble
— MethodOne can also call this without the iteration level
EnsembleKalmanProcesses.TOMLInterface.write_log_file
— Methodwrite_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.ConstantNesterovAccelerator
— Typemutable 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.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 iterationinv_sqrt_noise::Vector
: the inverse square-root of the noise covariance is storedterminate_at::Any
: the algorithm time for termination, default: 1.0on_terminate::Any
: the action on termination, default: "stop",
EnsembleKalmanProcesses.DefaultAccelerator
— Typestruct DefaultAccelerator <: EnsembleKalmanProcesses.Accelerator
Default accelerator provides no acceleration, runs traditional EKI
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.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.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
]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 sizeg::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 errorsscheduler::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 iterationprocess::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 toRandom.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 ensembleobs_mean
:: Vector of observationsobs_noise_cov
:: Noise covariance associated with the observationsobs_mean
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,
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.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)
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.FailureHandler
— MethodFailureHandler(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.FailureHandler
— MethodFailureHandler(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.FailureHandler
— MethodFailureHandler(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.FailureHandler
— MethodFailureHandler(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.FirstOrderNesterovAccelerator
— Typemutable 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.IgnoreFailures
— TypeFailure handling method that ignores forward model failures
EnsembleKalmanProcesses.Inversion
— TypeInversion <: Process
An ensemble Kalman Inversion process
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.NesterovAccelerator
— Typemutable 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.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.Sampler
— TypeSampler{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 spaceprior_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.SparseInversion
— TypeSparseInversion <: Process
A sparse ensemble Kalman Inversion process
Fields
γ::AbstractFloat
: upper limit of l1-normthreshold_value::AbstractFloat
: threshold below which the norm of parameters is pruned to zerouc_idx::Union{Colon, AbstractVector}
: indices of parameters included in the evaluation of l1-norm constraintreg::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.TransformInversion
— TypeTransformInversion <: 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.Unscented
— TypeUnscented{FT<:AbstractFloat, IT<:Int} <: Process
An unscented Kalman Inversion process.
Fields
u_mean::Any
: an interable of arrays of sizeN_parameters
containing the mean of the parameters (in eachuki
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 predictionuu_cov::Any
: an iterable of arrays of size (N_parameters x N_parameters
) containing the covariance of the parameters (in eachuki
iteration a new array ofcov
is added), note - this is not the same as the ensemble cov of the sigma ensemble as it is taken prior to predictionobs_pred::Any
: an iterable of arrays of sizeN_y
containing the predicted observation (in eachuki
iteration a new array of predicted observation is added)c_weights::AbstractVecOrMat{FT} where FT<:AbstractFloat
: weights in UKImean_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 parameterr::AbstractVector{FT} where FT<:AbstractFloat
: regularization vectorupdate_freq::Int64
: update frequencyimpose_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 meanprior_cov::Any
: prior covariance - defaults to initial covarianceiter::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 besymmetric
with2N_par+1
ensemble members orsimplex
withN_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!
— MethodState 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!
— MethodState 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!
— MethodState 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!
— MethodPerforms state update with modified constant Nesterov Accelerator.
EnsembleKalmanProcesses.accelerate!
— MethodPerforms traditional state update with no momentum.
EnsembleKalmanProcesses.accelerate!
— MethodPerforms state update with modified first-order Nesterov Accelerator.
EnsembleKalmanProcesses.accelerate!
— MethodPerforms 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!
— Methodadditive_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!
— Methodcalculate_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!
— Methodcompute_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_cov
— Methodconstruct_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_cov
— Methodconstruct_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_ensemble
— Methodconstruct_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_mean
— Methodconstruct_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_ensemble
— Methodconstruct_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_cov
— Methodconstruct_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_cov
— Methodconstruct_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_mean
— Methodconstruct_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_update
— Method 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_update
— Method 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_update
— Method 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_N_iterations
— Methodget_N_iterations(ekp::EnsembleKalmanProcess)
Get number of times update has been called (equals size(g)
, or size(u)-1
).
EnsembleKalmanProcesses.get_accelerator
— Methodget_accelerator(ekp::EnsembleKalmanProcess)
Return accelerator type of EnsembleKalmanProcess.
EnsembleKalmanProcesses.get_cov_blocks
— Methodget_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_error
— Methodget_error(ekp::EnsembleKalmanProcess)
Returns the mean forward model output error as a function of algorithmic time.
EnsembleKalmanProcesses.get_g
— Methodget_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_g
— Methodget_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_final
— Methodget_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_mean
— Methodget_g_mean(ekp::EnsembleKalmanProcess, iteration::IT) where {IT <: Integer}
Returns the mean forward map evaluation at the given iteration.
EnsembleKalmanProcesses.get_g_mean_final
— Methodget_g_mean_final(ekp::EnsembleKalmanProcess)
Get the mean forward model evaluation at the last iteration.
EnsembleKalmanProcesses.get_localizer
— Methodget_localizer(ekp::EnsembleKalmanProcess)
Return localizer type of EnsembleKalmanProcess.
EnsembleKalmanProcesses.get_process
— Methodget_process(ekp::EnsembleKalmanProcess)
Return process type of EnsembleKalmanProcess.
EnsembleKalmanProcesses.get_scheduler
— Methodget_scheduler(ekp::EnsembleKalmanProcess)
Return scheduler type of EnsembleKalmanProcess.
EnsembleKalmanProcesses.get_u
— Methodget_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_u
— Methodget_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_cov
— Methodget_u_cov(ekp::EnsembleKalmanProcess, iteration::IT) where {IT <: Integer}
Returns the unconstrained parameter sample covariance at the given iteration.
EnsembleKalmanProcesses.get_u_cov
— Methodget_u_cov(uki::EnsembleKalmanProcess{FT, IT, Unscented}, iteration::IT)
Returns the unconstrained parameter covariance at the requested iteration.
EnsembleKalmanProcesses.get_u_cov_final
— Methodget_u_cov_final(ekp::EnsembleKalmanProcess)
Get the mean unconstrained parameter covariance at the last iteration.
EnsembleKalmanProcesses.get_u_cov_prior
— Methodget_u_cov_prior(ekp::EnsembleKalmanProcess)
Returns the unconstrained parameter sample covariance for the initial ensemble.
EnsembleKalmanProcesses.get_u_final
— Methodget_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_mean
— Methodget_u_mean(ekp::EnsembleKalmanProcess, iteration::IT) where {IT <: Integer}
Returns the mean unconstrained parameter at the given iteration.
EnsembleKalmanProcesses.get_u_mean
— Methodget_u_mean(uki::EnsembleKalmanProcess{FT, IT, Unscented}, iteration::IT)
Returns the mean unconstrained parameter at the requested iteration.
EnsembleKalmanProcesses.get_u_mean_final
— Methodget_u_mean_final(ekp::EnsembleKalmanProcess)
Get the mean unconstrained parameter at the last iteration.
EnsembleKalmanProcesses.get_u_prior
— Methodget_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_ϕ
— Methodget_ϕ(prior::ParameterDistribution, ekp::EnsembleKalmanProcess)
Returns the constrained parameters from all iterations. The outer dimension is given by the number of iterations.
EnsembleKalmanProcesses.get_ϕ
— Methodget_ϕ(prior::ParameterDistribution, ekp::EnsembleKalmanProcess, iteration::IT)
Returns the constrained parameters at the given iteration.
EnsembleKalmanProcesses.get_ϕ_final
— Methodget_ϕ_final(ekp::EnsembleKalmanProcess)
Get the constrained parameters at the last iteration.
EnsembleKalmanProcesses.get_ϕ_mean
— Methodget_ϕ_mean(prior::ParameterDistribution, ekp::EnsembleKalmanProcess, iteration::IT)
Returns the constrained transform of the mean unconstrained parameter at the given iteration.
EnsembleKalmanProcesses.get_ϕ_mean_final
— Methodget_ϕ_mean_final(prior::ParameterDistribution, ekp::EnsembleKalmanProcess)
Get the constrained transform of the mean unconstrained parameter at the last iteration.
EnsembleKalmanProcesses.multiplicative_inflation!
— Methodmultiplicative_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_correct
— Methodposdef_correct(mat::AbstractMatrix; tol) -> Any
Makes square matrix mat
positive definite, by symmetrizing and bounding the minimum eigenvalue below by tol
EnsembleKalmanProcesses.sample_empirical_gaussian
— Methodsample_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.set_ICs!
— MethodSets u_prev to the initial parameter values
EnsembleKalmanProcesses.set_ICs!
— MethodSets u_prev to the initial parameter values
EnsembleKalmanProcesses.set_ICs!
— MethodSets u_prev to the initial parameter values
EnsembleKalmanProcesses.sparse_eki_update
— Method 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_qp
— Methodsparse_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_success
— Method 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!
— Methodupdate_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!
— Methodupdate_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!
— Methodupdate_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.
EnsembleKalmanProcesses.update_ensemble!
— Methodupdate_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!
— Methodupdate_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!
— Methodupdate_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.update_ensemble_analysis!
— Methodupdate_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.update_ensemble_prediction!
— Methodupdate_ensemble_prediction!(process::Unscented, Δt::FT) where {FT <: AbstractFloat}
UKI prediction step : generate sigma points.
EnsembleKalmanProcesses.Localizers.BernoulliDropout
— TypeBernoulliDropout{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.Delta
— TypeDirac delta localization method, with an identity matrix as the kernel.
EnsembleKalmanProcesses.Localizers.Localizer
— TypeLocalizer{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.Localizer
— MethodRandomized Bernoulli dropout kernel localizer constructor
EnsembleKalmanProcesses.Localizers.Localizer
— MethodDelta kernel localizer constructor
EnsembleKalmanProcesses.Localizers.Localizer
— MethodUniform kernel constructor
EnsembleKalmanProcesses.Localizers.Localizer
— MethodRBF kernel localizer constructor
EnsembleKalmanProcesses.Localizers.Localizer
— MethodSampling error correction (Lee, 2021) constructor
EnsembleKalmanProcesses.Localizers.Localizer
— MethodSampling error correction (Flowerdew, 2015) constructor
EnsembleKalmanProcesses.Localizers.NoLocalization
— TypeIdempotent localization method.
EnsembleKalmanProcesses.Localizers.RBF
— TypeRBF{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.SEC
— TypeSEC{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 correctionr_0::Real
: Cutoff correlation
EnsembleKalmanProcesses.Localizers.SECFisher
— TypeSECFisher <: 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
EnsembleKalmanProcesses.Localizers.bernoulli_kernel
— MethodLocalization kernel with Bernoulli trials as off-diagonal terms (symmetric)
EnsembleKalmanProcesses.Localizers.bernoulli_kernel_function
— MethodLocalize using a Schur product with a random draw of a Bernoulli kernel matrix. Only the u–G(u) block is localized.
EnsembleKalmanProcesses.Localizers.get_localizer
— Methodget_localizer(loc::Localizer)
Return localizer type.
EnsembleKalmanProcesses.Localizers.kernel_function
— MethodLocalize using a Schur product with a kernel matrix. Only the u–G(u) block is localized.
EnsembleKalmanProcesses.Localizers.sec
— MethodFunction that performs sampling error correction as per Lee (2021). The input is assumed to be a covariance matrix, hence square.
EnsembleKalmanProcesses.Localizers.sec_fisher
— MethodFunction that performs sampling error correction as per Flowerdew (2015). The input is assumed to be a covariance matrix, hence square.