ParameterDistributions
ParameterDistributionTypes
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
Constraints
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 domain
bounds::Union{Nothing, Dict}
Dictionary of values used to build the Constraint (e.g. "lowerbound" or "upperbound")
EnsembleKalmanProcesses.ParameterDistributions.no_constraint
— Functionno_constraint()
Constructs a Constraint with no constraints, enforced by maps x -> x and x -> x.
EnsembleKalmanProcesses.ParameterDistributions.bounded_below
— Functionbounded_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.bounded_above
— Functionbounded_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
— Functionbounded(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)
.
Base.length
— Methodlength(c<:ConstraintType)
A constraint has length 1.
Base.size
— Methodsize(c<:ConstraintType)
A constraint has size 1.
ParameterDistributions
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 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)
defined at /home/runner/work/EnsembleKalmanProcesses.jl/EnsembleKalmanProcesses.jl/src/ParameterDistributions.jl:289
.
ParameterDistribution(param_dist_dict)
defined at /home/runner/work/EnsembleKalmanProcesses.jl/EnsembleKalmanProcesses.jl/src/ParameterDistributions.jl:305
.
ParameterDistribution(distribution, constraint, name)
defined at /home/runner/work/EnsembleKalmanProcesses.jl/EnsembleKalmanProcesses.jl/src/ParameterDistributions.jl:389
.
ParameterDistribution(distribution_samples, constraint, name)
defined at /home/runner/work/EnsembleKalmanProcesses.jl/EnsembleKalmanProcesses.jl/src/ParameterDistributions.jl:431
.
EnsembleKalmanProcesses.ParameterDistributions.constrained_gaussian
— Functionconstrained_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.n_samples
— Functionn_samples(d<:Samples)
The number of samples in the array.
EnsembleKalmanProcesses.ParameterDistributions.get_name
— Functionget_name(pd::ParameterDistribution)
Returns a list of ParameterDistribution names.
EnsembleKalmanProcesses.ParameterDistributions.get_dimensions
— Functionget_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_n_samples
— Functionget_n_samples(pd::ParameterDistribution)
The number of samples in a Samples distribution
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_constraint_type
— Functionget_bounds(c::Constraint{T})
Gets the parametric type T.
EnsembleKalmanProcesses.ParameterDistributions.get_bounds
— Functionget_bounds(c::Constraint)
Gets the bounds field from the Constraint.
EnsembleKalmanProcesses.ParameterDistributions.batch
— Functionbatch(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.get_distribution
— Functionget_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.
gets the, distribution over the coefficients
StatsBase.sample
— Functionsample([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.
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.
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.
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.
Distributions.logpdf
— Functionlogpdf(pd::ParameterDistribution, xarray::Array{<:Real,1})
Obtains the independent logpdfs of the parameter distributions at xarray
(non-Samples Distributions only), and returns their sum.
Statistics.mean
— Functionmean(pd::ParameterDistribution)
Returns a concatenated mean of the parameter distributions.
Statistics.var
— Functionvar(pd::ParameterDistribution)
Returns a flattened variance of the distributions
Statistics.cov
— Functioncov(pd::ParameterDistribution)
Returns a dense blocked (co)variance of the distributions.
EnsembleKalmanProcesses.ParameterDistributions.transform_constrained_to_unconstrained
— Methodtransform_constrained_to_unconstrained(pd::ParameterDistribution, x::iterable{AbstractVecOrMat})
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::ParameterDistribution, x::iterable{AbstractVecOrMat})
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_unconstrained_to_constrained
— Methodtransform_unconstrained_to_constrained(pd::ParameterDistribution, x::iterable{AbstractVecOrMat})
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::ParameterDistribution, x::iterable{AbstractVecOrMat})
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(d::ParameterDistribution, x::Dict; build_flag = true)
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 GaussianRandomFieldsPackage
Type to dispatch which Gaussian Random Field package to use:
GRFJL
uses the Julia PackageGaussianRandomFields.jl
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 function
package::EnsembleKalmanProcesses.ParameterDistributions.GaussianRandomFieldsPackage
the choice of GRF package
distribution::EnsembleKalmanProcesses.ParameterDistributions.ParameterDistribution
the distribution of the coefficients that we shall compute with
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
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.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(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_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(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.