ONBasis
An ONBasis (orthonormal basis) encapsulates information about the orthogonal polynomials associated with a given probability distribution. This includes norms, quadrature rules, and cached evaluations at quadrature points. The ONBasis serves as a fundamental building block for constructing the tensorized basis linked to the multi-indices in the stochastic discretization.
ExtendableASGFEM.ONBasis
— TypeONBasis(
OBT::Type{<:OrthogonalPolynomialType},
maxorder;
...
) -> ONBasis{Float64}
ONBasis(
OBT::Type{<:OrthogonalPolynomialType},
maxorder,
maxquadorder;
T
) -> ONBasis{Float64}
Constructs an ONBasis
for the given orthogonal polynomial type and associated weight function, up to the specified maximum polynomial order.
Arguments
OBT
: The type of orthogonal polynomial (subtype ofOrthogonalPolynomialType
).maxorder
: The highest polynomial degree to include in the basis.maxquadorder
: The quadrature order to use for integration (default:2 * maxorder
).T
: The floating-point type for computations (default:Float64
).
Returns
An ONBasis
object containing:
- The norms of all basis functions.
- A workspace for storing evaluations.
- The quadrature rule (points and weights).
- Precomputed values of all basis functions at the quadrature points.
ExtendableASGFEM.ONBasis
— Typestruct ONBasis{T<:Real, OBT<:OrthogonalPolynomialType, npoly, nquad}
A structure that stores all relevant information for an orthonormal polynomial basis (ONBasis), including:
- The norms of all basis functions.
- The quadrature rule (points and weights) used for integration.
- Precomputed values of all basis functions at the quadrature points.
- A workspace for storing the result of the last evaluation.
This structure enables efficient evaluation, integration, and manipulation of orthogonal polynomial bases for stochastic Galerkin methods and related applications.
ExtendableASGFEM.OrthogonalPolynomialType
— MethodOrthogonalPolynomialType(ONB::ONBasis{T, OBT}) -> Any
returns the OrthogonalPolynomialType
ExtendableASGFEM.distribution
— Methoddistribution(ONB::ONBasis{T, OBT}) -> Any
returns the distribution associated to the orthogonal polynomials
ExtendableASGFEM.evaluate
— Methodevaluate(
ONB::ONBasis{T, OBT, maxorder},
x;
normalize
) -> Vector
Evaluates all basis polynomials at point x
ExtendableASGFEM.integral
— Methodintegral(ONB::ONBasis, j; normalize) -> Any
Evaluates the integral of a single basis function.
ExtendableASGFEM.norm4poly
— Methodnorm4poly(ONB::ONBasis, p) -> Any
returns the norm of the p-th polynomial
ExtendableASGFEM.qp
— Methodqp(
ONB::ONBasis
) -> StaticArraysCore.SVector{nquad, T} where {T<:Real, nquad}
returns the quadrature points
ExtendableASGFEM.qw
— Methodqw(
ONB::ONBasis
) -> StaticArraysCore.SVector{nquad, T} where {T<:Real, nquad}
returns the quadrature weights
ExtendableASGFEM.scalar_product
— Methodscalar_product(ONB::ONBasis, j, k; normalize) -> Any
Evaluates the scalar product of two basis functions.
ExtendableASGFEM.triple_product
— Methodtriple_product(ONB::ONBasis, j, k, l; normalize) -> Any
Evaluates the triple product of three basis functions.
ExtendableASGFEM.triple_product_y
— Methodtriple_product_y(ONB::ONBasis, j, k; normalize) -> Any
Evaluates the triple product of two basis functions and y.
ExtendableASGFEM.vals4poly
— Methodvals4poly(ONB::ONBasis, p) -> Any
returns the values of the p-th polynomial at all quadrature points
ExtendableASGFEM.vals4qp
— Methodvals4qp(ONB::ONBasis, k) -> Any
returns the values of all polynomials at the k-th quadrature point
ExtendableASGFEM.vals4xref
— Methodvals4xref(
ONB::ONBasis
) -> StaticArraysCore.SMatrix{nquad, npoly, T} where {T<:Real, npoly, nquad}
returns the values of all polynomials at the quadrature points