TensorizedBasis
Each multi-index $\mu = [\mu_1, \mu_2, \ldots, \mu_M]$ defines a tensorized basis function for the parameter space of the form $H_\mu = \prod_{k=1}^M H_k$, where each $H_k$ is an orthogonal polynomial.
The TensorizedBasis
object collects all information required to evaluate these basis functions, including the set of multi-indices and the triple products of the form $(y_m H_\mu, H_\lambda)$ for each $m$ and $\mu, \lambda$ in the set of multi-indices, stored as a sparse matrix. While analytic formulas exist to compute these triple products using recurrence coefficients, storing them can significantly speed up evaluations.
ExtendableASGFEM.TensorizedBasis
— TypeTensorizedBasis(
OBT::Type{<:OrthogonalPolynomialType},
M,
order,
maxorder;
...
) -> TensorizedBasis{Float64, ONBasis{T, OBT, npoly, nquad}} where {T<:Real, OBT<:OrthogonalPolynomialType, npoly, nquad}
TensorizedBasis(
OBT::Type{<:OrthogonalPolynomialType},
M,
order,
maxorder,
maxquadorder;
T,
multi_indices
) -> TensorizedBasis{Float64, ONBasis{T, OBT, npoly, nquad}} where {T<:Real, OBT<:OrthogonalPolynomialType, npoly, nquad}
constructor for a tensorized basis for the given OrthogonalPolynomialType. If no multi-indices are provided it automatically generates all multi-indices up to support length M and polynomial order maxorder.
ExtendableASGFEM.TensorizedBasis
— Typestruct TensorizedBasis{T<:Real, ONBType<:ONBasis, MIType}
A structure representing a tensorized orthogonal polynomial basis for stochastic Galerkin methods, associated with a given set of multi-indices.
Fields
ONB::ONBType
: The underlying univariate orthogonal basis (e.g., Hermite, Legendre).vals::Vector{Vector{T}}
: Cached evaluations of the univariate basis functions at the most recent sample(s).nmodes::Int
: Number of multi-indices (i.e., the number of stochastic modes).G::ExtendableSparseMatrix{T, Int64}
: Coupling matrix encoding triple products or recurrence relations between basis functions.multi_indices::MIType
: Collection of multi-indices (typically an array of integer arrays), each representing a multi-dimensional polynomial degree.
Description
The TensorizedBasis
type encapsulates all information required to evaluate and manipulate a tensor-product polynomial basis in multiple stochastic dimensions. It supports efficient evaluation, storage of basis values, and access to coupling coefficients for use in stochastic Galerkin finite element methods (SGFEM).
Example
# Construct a tensorized Hermite basis with 3 variables and total degree 2
TB = TensorizedBasis(HermitePolynomials, 3, 2, 2)
# Evaluate all basis functions at a sample point
set_sample!(TB, [0.1, -0.2, 0.3])
# Get the value of the 5th basis function at the current sample
val = evaluate(TB, 5)
Base.show
— Methodshow(io::IO, TB::TensorizedBasis)
shows information on the tensorized basis.
ExtendableASGFEM.distribution
— Methoddistribution(TB::TensorizedBasis{T, ONBT}) -> Any
returns distribution associated to the orthogonal basis
ExtendableASGFEM.evaluate
— Methodevaluate(TB::TensorizedBasis{T}, j) -> Any
evaluates the basis function for the j-th multi-index at the sample that was set with set_sample.
ExtendableASGFEM.get_coupling_coefficient
— Methodget_coupling_coefficient(
TB::TensorizedBasis{T},
m,
j,
k
) -> Any
returns the triple products between $y_m$ and $H_j$ and $H_k$ for two multi-indices j
and k
.
ExtendableASGFEM.get_multiindex
— Methodget_multiindex(TB::TensorizedBasis, j) -> Any
returns the j-th multi-index.
ExtendableASGFEM.maxlength_multiindices
— Methodmaxlength_multiindices(TB::TensorizedBasis{T, ONBT}) -> Any
returns the maximal length of the stored multi-indices
ExtendableASGFEM.num_multiindices
— Methodnum_multiindices(TB::TensorizedBasis) -> Int64
returns the number of multi-indices
ExtendableASGFEM.sample_distribution
— Methodsample_distribution(
TB::TensorizedBasis,
nsamples;
M,
Mweights
) -> Tuple{Matrix{Float64}, Any}
Generate samples and weights for the distribution of the ON basis that the tensorized basis is based upon (that can be used for a Monte-Carlo estimator).
ExtendableASGFEM.set_sample!
— Methodset_sample!(
TB::TensorizedBasis,
x;
normalize
) -> Array{Vector{T}, 1} where T<:Real
evaluates all basis functions at sample vector x; call this before using evaluate!