ParameterDistributions
ParameterDistributionTypes
EnsembleKalmanProcesses.ParameterDistributions.Parameterized — Typestruct Parameterized <: EnsembleKalmanProcesses.ParameterDistributions.ParameterDistributionTypeA distribution constructed from a parameterized formula (e.g Julia Distributions.jl)
Fields
distribution::Distributions.Distribution: A parameterized distribution
EnsembleKalmanProcesses.ParameterDistributions.Samples — Typestruct Samples{FT<:Real} <: EnsembleKalmanProcesses.ParameterDistributions.ParameterDistributionTypeA 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 — Typestruct VectorOfParameterized{DT<:Distributions.Distribution} <: EnsembleKalmanProcesses.ParameterDistributions.ParameterDistributionTypeA 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
Constraints
EnsembleKalmanProcesses.ParameterDistributions.Constraint — Typestruct Constraint{T} <: EnsembleKalmanProcesses.ParameterDistributions.ConstraintTypeClass 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.no_constraint — Functionno_constraint(
) -> EnsembleKalmanProcesses.ParameterDistributions.Constraint{EnsembleKalmanProcesses.ParameterDistributions.NoConstraint}
Constructs a Constraint with no constraints, enforced by maps x -> x and x -> x.
EnsembleKalmanProcesses.ParameterDistributions.bounded_below — Functionbounded_below(
lower_bound::Real
) -> Union{EnsembleKalmanProcesses.ParameterDistributions.Constraint{EnsembleKalmanProcesses.ParameterDistributions.BoundedBelow}, EnsembleKalmanProcesses.ParameterDistributions.Constraint{EnsembleKalmanProcesses.ParameterDistributions.NoConstraint}}
Constructs a Constraint with provided lower bound, enforced by maps x -> log(x - lower_bound) and x -> exp(x) + lower_bound.
EnsembleKalmanProcesses.ParameterDistributions.bounded_above — Functionbounded_above(
upper_bound::Real
) -> Union{EnsembleKalmanProcesses.ParameterDistributions.Constraint{EnsembleKalmanProcesses.ParameterDistributions.BoundedAbove}, EnsembleKalmanProcesses.ParameterDistributions.Constraint{EnsembleKalmanProcesses.ParameterDistributions.NoConstraint}}
Constructs a Constraint with provided upper bound, enforced by maps x -> log(upper_bound - x) and x -> upper_bound - exp(x).
EnsembleKalmanProcesses.ParameterDistributions.bounded — Functionbounded(
lower_bound::Real,
upper_bound::Real
) -> EnsembleKalmanProcesses.ParameterDistributions.Constraint
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).
Base.length — Methodlength(
c::EnsembleKalmanProcesses.ParameterDistributions.ConstraintType
) -> Int64
A constraint has length 1.
Base.size — Methodsize(
c::EnsembleKalmanProcesses.ParameterDistributions.ConstraintType
) -> Tuple{Int64}
A constraint has size 1.
ParameterDistributions
EnsembleKalmanProcesses.ParameterDistributions.ParameterDistribution — Typestruct ParameterDistribution{PDType<:EnsembleKalmanProcesses.ParameterDistributions.ParameterDistributionType, CType<:EnsembleKalmanProcesses.ParameterDistributions.ConstraintType, ST<:AbstractString}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
Recommended construction (for most problems) is via the constrained_gaussian() utility.
ParameterDistribution(distribution, constraint, name)defined at /home/runner/work/EnsembleKalmanProcesses.jl/EnsembleKalmanProcesses.jl/src/ParameterDistributions.jl:291.
ParameterDistribution(param_dist_dict)defined at /home/runner/work/EnsembleKalmanProcesses.jl/EnsembleKalmanProcesses.jl/src/ParameterDistributions.jl:309.
ParameterDistribution(distribution, constraint, name)defined at /home/runner/work/EnsembleKalmanProcesses.jl/EnsembleKalmanProcesses.jl/src/ParameterDistributions.jl:391.
ParameterDistribution(
distribution_samples,
constraint,
name;
params_are_columns
)defined at /home/runner/work/EnsembleKalmanProcesses.jl/EnsembleKalmanProcesses.jl/src/ParameterDistributions.jl:430.
EnsembleKalmanProcesses.ParameterDistributions.constrained_gaussian — Functionconstrained_gaussian(
name::AbstractString,
μ_c::Real,
σ_c::Real,
lower_bound::Real,
upper_bound::Real;
repeats,
optim_algorithm,
optim_kwargs...
) -> Union{EnsembleKalmanProcesses.ParameterDistributions.ParameterDistribution{EnsembleKalmanProcesses.ParameterDistributions.Parameterized, T} where T<:EnsembleKalmanProcesses.ParameterDistributions.Constraint, EnsembleKalmanProcesses.ParameterDistributions.ParameterDistribution{T, T1} where {T<:(EnsembleKalmanProcesses.ParameterDistributions.VectorOfParameterized{T} where T<:Distributions.Normal), T1<:EnsembleKalmanProcesses.ParameterDistributions.Constraint}}
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.
kwargs:
repeats=1: K-dimensional distribution by stackingKindependent copies of the defined univariate distribution by settingrepeats = K.
Example usage:
d1 = constrained_gaussian("mean2-sd1-positive", 2.0, 1.0, 0, Inf)
d2 = constrained_gaussian("4-dim-mean0-sd10", 0.0, 4.0, -Inf, Inf, repeats=4)
# combine with:
d = combine_distributions([d1,d2])EnsembleKalmanProcesses.ParameterDistributions.n_samples — Functionn_samples(
d::EnsembleKalmanProcesses.ParameterDistributions.Samples
) -> Any
The number of samples in the array.
EnsembleKalmanProcesses.ParameterDistributions.get_name — Functionget_name(
pd::EnsembleKalmanProcesses.ParameterDistributions.ParameterDistribution
) -> AbstractVector{ST} where ST<:AbstractString
Returns a list of ParameterDistribution names.
EnsembleKalmanProcesses.ParameterDistributions.get_dimensions — Functionget_dimensions(
pd::EnsembleKalmanProcesses.ParameterDistributions.ParameterDistribution;
function_parameter_opt
) -> Any
The number of dimensions of the parameter space. (Also represents other dimensions of interest for FunctionParameterDistributionTypes with keyword argument)
EnsembleKalmanProcesses.ParameterDistributions.get_n_samples — Functionget_n_samples(
pd::EnsembleKalmanProcesses.ParameterDistributions.ParameterDistribution
) -> Dict{String, Any}
The number of samples in a Samples distribution
EnsembleKalmanProcesses.ParameterDistributions.get_all_constraints — Methodget_all_constraints(
pd::EnsembleKalmanProcesses.ParameterDistributions.ParameterDistribution;
return_dict
) -> Union{Dict{Any, Any}, AbstractVector{CType} where CType<:EnsembleKalmanProcesses.ParameterDistributions.ConstraintType}
Returns the (flattened) array of constraints of the parameter distribution. or as a dictionary ("param_name" => constraints)
EnsembleKalmanProcesses.ParameterDistributions.get_constraint_type — Functionget_constraint_type(
c::EnsembleKalmanProcesses.ParameterDistributions.Constraint{T}
) -> Any
Gets the parametric type T.
EnsembleKalmanProcesses.ParameterDistributions.get_bounds — Functionget_bounds(
c::EnsembleKalmanProcesses.ParameterDistributions.Constraint
) -> Union{Nothing, Dict}
Gets the bounds field from the Constraint.
EnsembleKalmanProcesses.ParameterDistributions.batch — Functionbatch(
pd::Union{EnsembleKalmanProcesses.ParameterDistributions.ParameterDistribution, EnsembleKalmanProcesses.ParameterDistributions.VectorOfParameterized};
function_parameter_opt
) -> Any
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.get_distribution — Functionget_distribution(
pd::EnsembleKalmanProcesses.ParameterDistributions.ParameterDistribution
) -> Dict{String, Any}
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.
get_distribution(
grfi::EnsembleKalmanProcesses.ParameterDistributions.GaussianRandomFieldInterface
) -> EnsembleKalmanProcesses.ParameterDistributions.ParameterDistribution
gets the, distribution over the coefficients
StatsBase.sample — Functionsample(
rng::Random.AbstractRNG,
pd::EnsembleKalmanProcesses.ParameterDistributions.ParameterDistribution,
n_draws::Integer
) -> Any
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.
sample(
rng::Random.AbstractRNG,
d::EnsembleKalmanProcesses.ParameterDistributions.Samples,
n_draws::Integer
) -> Any
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.
sample(
rng::Random.AbstractRNG,
d::EnsembleKalmanProcesses.ParameterDistributions.Parameterized,
n_draws::Integer
) -> Any
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.
sample(
rng::Random.AbstractRNG,
d::EnsembleKalmanProcesses.ParameterDistributions.VectorOfParameterized,
n_draws::Integer
) -> Any
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.
Distributions.logpdf — Functionlogpdf(
d::EnsembleKalmanProcesses.ParameterDistributions.Parameterized,
x::Real
) -> Any
Obtains the independent logpdfs of the parameter distributions at xarray (non-Samples Distributions only), and returns their sum.
Statistics.mean — Functionmean(
d::EnsembleKalmanProcesses.ParameterDistributions.Parameterized
) -> Any
Returns a concatenated mean of the parameter distributions.
Statistics.var — Functionvar(
d::EnsembleKalmanProcesses.ParameterDistributions.Parameterized
) -> Any
Returns a flattened variance of the distributions
Statistics.cov — Functioncov(
d::EnsembleKalmanProcesses.ParameterDistributions.Parameterized
) -> Any
Returns a dense blocked (co)variance of the distributions.
EnsembleKalmanProcesses.ParameterDistributions.transform_constrained_to_unconstrained — Methodtransform_constrained_to_unconstrained(
pd::EnsembleKalmanProcesses.ParameterDistributions.ParameterDistribution,
x
) -> Union{Matrix{T} where T<:Real, Vector{T} where T<:Real}
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(
pd::EnsembleKalmanProcesses.ParameterDistributions.ParameterDistribution,
x
) -> Union{Matrix{T} where T<:Real, Vector{T} where T<:Real}
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(
pd::EnsembleKalmanProcesses.ParameterDistributions.ParameterDistribution,
x::Dict
) -> Dict{Any, Any}
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(
pd::EnsembleKalmanProcesses.ParameterDistributions.ParameterDistribution,
x
) -> Union{Matrix{T} where T<:Real, Vector{T} where T<:Real}
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. The build_flag will construct any function parameters
EnsembleKalmanProcesses.ParameterDistributions.transform_unconstrained_to_constrained — Methodtransform_unconstrained_to_constrained(
pd::EnsembleKalmanProcesses.ParameterDistributions.ParameterDistribution,
x
) -> Union{Matrix{T} where T<:Real, Vector{T} where T<:Real}
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. The build_flag will construct any function parameters
EnsembleKalmanProcesses.ParameterDistributions.transform_unconstrained_to_constrained — Methodtransform_unconstrained_to_constrained(
pd::EnsembleKalmanProcesses.ParameterDistributions.ParameterDistribution,
x::Dict;
build_flag
) -> Dict{Any, Any}
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. The build_flag will reconstruct any function parameters onto their flattened, discretized, domains.
FunctionParameterDistributions
EnsembleKalmanProcesses.ParameterDistributions.GaussianRandomFieldsPackage — Typeabstract type GaussianRandomFieldsPackageType to dispatch which Gaussian Random Field package to use:
GRFJLuses the Julia PackageGaussianRandomFields.jl
EnsembleKalmanProcesses.ParameterDistributions.GaussianRandomFieldInterface — Typestruct GaussianRandomFieldInterface <: EnsembleKalmanProcesses.ParameterDistributions.FunctionParameterDistributionTypeGaussianRandomFieldInterface 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
Base.ndims — Methodndims(
grfi::EnsembleKalmanProcesses.ParameterDistributions.GaussianRandomFieldInterface;
function_parameter_opt
) -> Any
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
EnsembleKalmanProcesses.ParameterDistributions.get_all_constraints — Methodget_all_constraints(
grfi::EnsembleKalmanProcesses.ParameterDistributions.GaussianRandomFieldInterface
) -> Union{Dict{Any, Any}, AbstractVector{CType} where CType<:EnsembleKalmanProcesses.ParameterDistributions.ConstraintType}
gets all the constraints of the internally stored coefficient prior distribution of the GRFI
EnsembleKalmanProcesses.ParameterDistributions.transform_constrained_to_unconstrained — Methodtransform_constrained_to_unconstrained(
d::EnsembleKalmanProcesses.ParameterDistributions.GaussianRandomFieldInterface,
constraint::AbstractVector,
x::AbstractArray{FT<:Real, 1}
) -> Any
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(
d::EnsembleKalmanProcesses.ParameterDistributions.GaussianRandomFieldInterface,
constraint::AbstractVector,
x::AbstractArray{FT<:Real, 2}
) -> Any
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_unconstrained_to_constrained — Methodtransform_unconstrained_to_constrained(
d::EnsembleKalmanProcesses.ParameterDistributions.GaussianRandomFieldInterface,
constraint::AbstractVector,
x::AbstractArray{FT<:Real, 1};
build_flag
) -> Any
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
xinto 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
constraintto 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(
d::EnsembleKalmanProcesses.ParameterDistributions.GaussianRandomFieldInterface,
constraint::AbstractVector,
x::AbstractArray{FT<:Real, 2};
build_flag
) -> Any
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
xinto 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
constraintto 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.get_grf — Functionget_grf(
grfi::EnsembleKalmanProcesses.ParameterDistributions.GaussianRandomFieldInterface
) -> Any
gets the distribution, i.e. Gaussian random field object
EnsembleKalmanProcesses.ParameterDistributions.build_function_sample — Functionbuild_function_sample(
grfi::EnsembleKalmanProcesses.ParameterDistributions.GaussianRandomFieldInterface,
coeff_vecormat::AbstractVecOrMat,
n_draws::Int64
) -> Any
build function n_draw times on the discrete grid, given the coefficients coeff_vecormat.
Defaults: n_draw = 1.
build_function_sample(
rng::Random.AbstractRNG,
grfi::EnsembleKalmanProcesses.ParameterDistributions.GaussianRandomFieldInterface,
n_draws::Int64
) -> Any
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.get_package — Functionget_package(
grfi::EnsembleKalmanProcesses.ParameterDistributions.GaussianRandomFieldInterface
) -> EnsembleKalmanProcesses.ParameterDistributions.GaussianRandomFieldsPackage
gets the package type used to construct the GRF
EnsembleKalmanProcesses.ParameterDistributions.spectrum — Functionspectrum(
grfi::EnsembleKalmanProcesses.ParameterDistributions.GaussianRandomFieldInterface
) -> Any
the spectral information of the GRF, e.g. the Karhunen-Loeve coefficients and eigenfunctions if using this decomposition
EnsembleKalmanProcesses.ParameterDistributions.n_dofs — Functionn_dofs(
grfi::EnsembleKalmanProcesses.ParameterDistributions.GaussianRandomFieldInterface
) -> Any
the number of degrees of freedom / coefficients (i.e. the number of parameters)
EnsembleKalmanProcesses.ParameterDistributions.eval_pts — Functioneval_pts(
grfi::EnsembleKalmanProcesses.ParameterDistributions.GaussianRandomFieldInterface
) -> Any
the discrete evaluation point grid, stored as a range in each dimension
EnsembleKalmanProcesses.ParameterDistributions.n_eval_pts — Functionn_eval_pts(
grfi::EnsembleKalmanProcesses.ParameterDistributions.GaussianRandomFieldInterface
) -> Any
the number of total discrete evaluation points
EnsembleKalmanProcesses.ParameterDistributions.input_dims — Functioninput_dims(
grfi::EnsembleKalmanProcesses.ParameterDistributions.GaussianRandomFieldInterface
) -> Any
the number of input dimensions of the GRF