Stochastic Coefficients
Stochastic coefficients play a central role in uncertainty quantification and stochastic finite element methods. In this package, coefficients are represented using a Karhunen-Loève expansion (KLE), which allows for the efficient representation of random fields with prescribed covariance structure.
Overview
The Karhunen-Loève expansion expresses a stochastic process as a series of orthogonal functions weighted by uncorrelated random variables:
$
a(x, \omega) = a_0(x) + \sum_{n=1}^N \sqrt{\lambda_n} \phi_n(x) \xi_n(\omega)$
where $a_0(x)$ is the mean, $\lambda_n$ and $\phi_n(x)$ are the eigenvalues and eigenfunctions of the covariance operator, and $\xi_n$ are independent standard normal random variables.
API
ExtendableASGFEM.AbstractStochasticCoefficient — Typeabstract type AbstractStochasticCoefficient{T}A stochastic coefficient is assumed to have the Karhunen-Loeve expansion form
$a(y,x) = a_0(x) + \sum_{m=1}^\infty y_m a_m(x)$ with (centered independent) random variables $y_m$ and basis functions $a_m(x)$ that need to be specified (together with their gradients) and expectation value $a_0$ (in general they steam from a spectral analysis of the covariance operator of $a$).
ExtendableASGFEM.expa_PCE_mop — Methodexpa_PCE_mop(
TB::TensorizedBasis{T, ONBType},
SC::AbstractStochasticCoefficient{T};
N_truncate,
factor
) -> Union{Tuple{ExtendableASGFEM.var"#closure#48"{_A, TensorizedBasis{T, ONBType, MIType}, ExtendableASGFEM.var"#lambda_mu#47"{T1, N_truncate, factor, SC, mu_fac_storage, multi_indices, eval_amu_storage, eval_am}, Int64} where {_A, T, ONBType, MIType, T1, N_truncate, factor, SC, mu_fac_storage, multi_indices, eval_amu_storage, eval_am}, ExtendableASGFEM.var"#lambda_mu#47"{_A, Int64, Int64, SingleStochasticCoefficient{T, k}} where {_A, T, k}}, Tuple{ExtendableASGFEM.var"#closure#48"{_A, TensorizedBasis{T, ONBType, MIType}, ExtendableASGFEM.var"#lambda_mu#47"{T1, N_truncate, factor, SC, mu_fac_storage, multi_indices, eval_amu_storage, eval_am}, Int64} where {_A, T, ONBType, MIType, T1, N_truncate, factor, SC, mu_fac_storage, multi_indices, eval_amu_storage, eval_am}, ExtendableASGFEM.var"#lambda_mu#47"{_A, Int64, Int64, StochasticCoefficientCosinus{T}} where {_A, T}}}
prepares two functions
lambda_mu!(result, input, qpinfo)
expa!(result, input, qpinfo)that calculate λ_μ, which are orthogonal decomposition coefficient functions of exp(a) ≈ ∑ λ_μ H_μ w.r.t. to the multi-indices μ and their associated orthogonal basis functions H_μ in TB, as well as (an approximation based on this decomposition of) exp(a).
ExtendableASGFEM.get_a! — Methodget_a!(
SC::AbstractStochasticCoefficient{T};
factor
) -> ExtendableASGFEM.var"#closure#34"{Int64, <:AbstractStochasticCoefficient{T}} where T
prepares a function of interface
a!(result, x, y)that evaluates the coefficient (times a factor) at space coordinates x and random variables y into result (a vector of length 1).
ExtendableASGFEM.get_am_x — Methodget_am_x(
m,
SC::AbstractStochasticCoefficient{T}
) -> ExtendableASGFEM.var"#closure#get_am_x##0"{_A, <:AbstractStochasticCoefficient{T}} where {_A, T}
prepares a function of interface
get_am_x!(result, input, qpinfo)that evaluates a_m(qpinfo.x) input (used to define ExtendableFEM operator kernels).
ExtendableASGFEM.get_expa! — Methodget_expa!(
SC::AbstractStochasticCoefficient;
factor
) -> ExtendableASGFEM.var"#closure#38"{Int64, ExtendableASGFEM.var"#closure#34"{Int64, var"#s179", _A}} where {T, var"#s179"<:AbstractStochasticCoefficient{T}, _A}
prepares a function of interface
get_expa!(result, x, y)that evaluates the exponential of the coefficient (times a factor) at space coordinates x and random variables y into result (a vector of length 1).
ExtendableASGFEM.get_grada! — Methodget_grada!(
SC::AbstractStochasticCoefficient{T};
factor
) -> ExtendableASGFEM.var"#closure#36"{Int64, <:AbstractStochasticCoefficient{T}} where T
prepares a function of interface
get_grada!(result, x, y)that evaluates the spatial gradient of the coefficient (times a factor) at space coordinates x and random variables y into result (a vector of same length as x).
ExtendableASGFEM.get_grada_x_sigma — Methodget_grada_x_sigma(
dim,
SC::AbstractStochasticCoefficient{T},
y
) -> ExtendableASGFEM.var"#closure#get_grada_x_sigma##0"{_A, <:AbstractStochasticCoefficient{T}} where {_A, T}
prepares a function of interface
get_grada_x_sigma!(result, input, qpinfo)that evaluates ∇a(qpinfo.x) ⋅ σ(qpinfo.x) (used to define ExtendableFEM operator kernels and σ is expected to be some vector-valued quantity of same length).
ExtendableASGFEM.get_gradam_x_u — Methodget_gradam_x_u(
m,
SC::AbstractStochasticCoefficient{T}
) -> ExtendableASGFEM.var"#closure#get_gradam_x_u##0"{_A, <:AbstractStochasticCoefficient{T}} where {_A, T}
prepares a function of interface
get_gradam_x_u!(result, input, qpinfo)that evaluates ∇a_m(qpinfo.x) u(qpinfo.x) (used to define ExtendableFEM operator kernels and input is expected to be some scalar quantity).
ExtendableASGFEM.plot_am — Methodplot_am(
xgrid::ExtendableGrids.ExtendableGrid,
m,
SC::AbstractStochasticCoefficient;
Plotter,
kwargs...
)
returns a scalarplot of the m-th coefficient functiona a_m interpolated on the given grid. The Plotter backend can be changed with the Plotter argument (e.g. GLMakie, CairoMakie, PyPlot, Plots).
StochasticCoefficientCosinus
An example of a particular KLE with cosinus type basis function is given by the following subtype.
ExtendableASGFEM.StochasticCoefficientCosinus — Typestruct StochasticCoefficientCosinus{T} <: AbstractStochasticCoefficient{T}Expansion of the form
a(x,y) = a_0(x) + ∑_m y_m a_m(x)where the a_m are of cosinus type.
ExtendableASGFEM.StochasticCoefficientCosinus — MethodStochasticCoefficientCosinus(
;
T,
τ,
start,
decay,
mean,
maxm
) -> StochasticCoefficientCosinus{Float64}
constructor for StochasticCoefficientCosinus of type T (default = Float64), where decay (default = 2) steers the decay of the coefficient basis functions (the larger the faster), mean (default = 0) is the mean value of the coefficient, maxm (default = 100) is the maximal number of stochastic random variables, and τ (default = 1) is a uniform scaling factors