ParameterDistributions

ParameterDistributionTypes

EnsembleKalmanProcesses.ParameterDistributions.SamplesType
Samples{FT <: Real} <: ParameterDistributionType

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

Fields

  • distribution_samples::AbstractMatrix{FT} where FT<:Real

    Samples defining an empirical distribution, stored as columns

source
EnsembleKalmanProcesses.ParameterDistributions.VectorOfParameterizedType
VectorOfParameterized <: ParameterDistributionType

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

Fields

  • distribution::AbstractVector{DT} where DT<:Distributions.Distribution

    A vector of parameterized distributions

source

Constraints

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

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

  • NoConstraint
  • BoundedBelow
  • BoundedAbove
  • Bounded

Fields

  • constrained_to_unconstrained::Function

    A map from constrained domain -> (-Inf,Inf)

  • c_to_u_jacobian::Function

    The jacobian of the map from constrained domain -> (-Inf,Inf)

  • unconstrained_to_constrained::Function

    Map from (-Inf,Inf) -> constrained domain

  • bounds::Union{Nothing, Dict}

    Dictionary of values used to build the Constraint (e.g. "lowerbound" or "upperbound")

source
EnsembleKalmanProcesses.ParameterDistributions.boundedFunction
bounded(lower_bound::Real, upper_bound::Real)

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

source

ParameterDistributions

EnsembleKalmanProcesses.ParameterDistributions.ParameterDistributionType
ParameterDistribution

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

Fields

  • distribution::AbstractVector{PDType} where PDType<:EnsembleKalmanProcesses.ParameterDistributions.ParameterDistributionType

    Vector of parameter distributions, defined in unconstrained space

  • constraint::AbstractVector{CType} where CType<:EnsembleKalmanProcesses.ParameterDistributions.ConstraintType

    Vector of constraints defining transformations between constrained and unconstrained space

  • name::AbstractVector{ST} where ST<:AbstractString

    Vector of parameter names

Constructors

ParameterDistribution(distribution, constraint, name)

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.

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

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

Note

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

Note

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

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

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

source
EnsembleKalmanProcesses.ParameterDistributions.get_distributionFunction
get_distribution(pd::ParameterDistribution)

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

source

gets the, distribution over the coefficients

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

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

source
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.

source
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.

source
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.

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

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

source
Statistics.meanFunction
mean(pd::ParameterDistribution)

Returns a concatenated mean of the parameter distributions.

source
Statistics.varFunction
var(pd::ParameterDistribution)

Returns a flattened variance of the distributions

source
Statistics.covFunction
cov(pd::ParameterDistribution)

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

source

FunctionParameterDistributions

EnsembleKalmanProcesses.ParameterDistributions.GaussianRandomFieldInterfaceType
struct GaussianRandomFieldInterface <: FunctionParameterDistributionType

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

Fields

  • gaussian_random_field::Any

    GRF object, containing the mapping from the field of coefficients to the discrete function

  • package::EnsembleKalmanProcesses.ParameterDistributions.GaussianRandomFieldsPackage

    the choice of GRF package

  • distribution::EnsembleKalmanProcesses.ParameterDistributions.ParameterDistribution

    the distribution of the coefficients that we shall compute with

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

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

  • "dof" : returns n_dofs(grfi), the degrees of freedom in the function
  • "eval" : returns n_eval_pts(grfi), the number of discrete evaluation points of the function
  • "constraint": returns 1, the number of constraints in the evaluation space
source
EnsembleKalmanProcesses.ParameterDistributions.transform_constrained_to_unconstrainedMethod
transform_constrained_to_unconstrained(d::GaussianRandomFieldInterface, constraint::AbstractVector, x::AbstractMatrix)

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

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

Optional Args build_flag::Bool = true

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

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

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

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

Optional args: build_flag::Bool = true

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

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

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

source