Orthogonal Polynomials

In the stochastic discretization of random variables, global polynomials that are orthogonal with respect to the probability distribution of each random variable $y_m$ are used. These orthogonal polynomials can be generated via recurrence relations, with coefficients determined by the underlying distribution.

Recurrence Relations

Orthogonal polynomials $H_n$ with respect to a weight function $\omega$ satisfy

\[\int_\Gamma \omega(y) H_{n}(y) H_m(y) \,dy = N^2_{nm}\delta_{nm}\]

where

\[N_{nn} := \| H_n \|_{\omega}^2 := \int_\Gamma \omega(y) H_{n}(y) H_n(y) \,dy\]

The polynomials satisfy the three-term recurrence relation:

\[\begin{aligned} H_{n+2}(y) & = (a_n y - b_n) H_{n+1}(y) - c_n H_{n}(y) \end{aligned}\]

with initial values $H_0 = 0$ and $H_1 = 1$.

ExtendableASGFEM.evaluate!Method
evaluate!(
    y,
    basis::Type{<:OrthogonalPolynomialType},
    n::Integer,
    x::Real
)

Evaluates the first n+1 orthogonal polynomials of the specified type at a single value x and writes the results into the vector y.

Arguments

  • y: Preallocated vector to store the polynomial values. Must have length at least n+1.
  • basis: The type of orthogonal polynomial (subtype of OrthogonalPolynomialType).
  • n: The highest polynomial degree to evaluate (computes degrees 0 to n).
  • x: The point at which to evaluate the polynomials.

Example

y = zeros(Float64, 6)
evaluate!(y, HermitePolynomials, 5, 0.0)
# y now contains H_0(0), H_1(0), ..., H_5(0)
source
ExtendableASGFEM.evaluateMethod
evaluate(
    basis::Type{<:OrthogonalPolynomialType},
    n::Integer,
    x::Real
) -> Any

Evaluates the first n+1 orthogonal polynomials of the specified type at a single value x and returns a vector containing the results.

Arguments

  • basis: The type of orthogonal polynomial (subtype of OrthogonalPolynomialType).
  • n: The highest polynomial degree to evaluate (computes degrees 0 to n).
  • x: The point at which to evaluate the polynomials.

Returns

A vector of length n+1 containing the values [p_0(x), p_1(x), ..., p_n(x)].

Example

y = evaluate(HermitePolynomials, 5, 0.0)
# y now contains H_0(0), H_1(0), ..., H_5(0)
source
ExtendableASGFEM.evaluateMethod
evaluate(
    basis::Type{<:OrthogonalPolynomialType},
    n::Integer,
    x::AbstractArray{T, 1}
) -> Any

Evaluates the first n+1 orthogonal polynomials at a vector of x

source
LinearAlgebra.normMethod
norm(basis::Type{<:OrthogonalPolynomialType}, n; gr) -> Any

Computes the norms of the first n+1 orthogonal polynomials of the specified type using Gauss quadrature.

Arguments

  • basis: The type of orthogonal polynomial (subtype of OrthogonalPolynomialType).
  • n: The highest polynomial degree for which to compute the norm (computes degrees 0 to n).
  • gr: Optional. A Gauss quadrature rule as a tuple (nodes, weights). By default, uses gauss_rule(basis, 2 * n).

Returns

A vector of length n+1 containing the L2 norms [||p_0||, ||p_1||, ..., ||p_n||] of the orthogonal polynomials with respect to the weight function of the basis.

Example

norms = norm(HermitePolynomials, 3)
# norms contains the L^2 norms of H_0, H_1, H_2, H_3
source

Legendre Polynomials (Uniform Distribution)

For the weight function $\omega(y) = 1/2$ on the interval $[-1,1]$ (uniform distribution), the recurrence coefficients are $a_n = (2n+1)/(n+1)$, $b_n = 0$, and $c_n = n/(n+1)$. The norms of the resulting Legendre polynomials are

\[ \| H_n \|^2_\omega = \frac{2}{2n+1}\]

ExtendableASGFEM.LegendrePolynomialsType
abstract type LegendrePolynomials <: OrthogonalPolynomialType

Type for dispatching Legendre polynomials, which are orthogonal with respect to the uniform distribution on [-1, 1].

source
ExtendableASGFEM.distributionMethod
distribution(
    _::Type{LegendrePolynomials}
) -> Distributions.Uniform{Float64}

Returns the probability distribution associated with the Legendre polynomials, which is the uniform distribution on [-1, 1].

source
ExtendableASGFEM.normsMethod
norms(_::Type{LegendrePolynomials}, k) -> Any

Returns the norm of the k-th Legendre polynomial, i.e.,

||P_k|| = sqrt(2 / (2k + 1))
source
ExtendableASGFEM.recurrence_coefficientsMethod
recurrence_coefficients(
    _::Type{LegendrePolynomials},
    k::Integer
) -> Tuple{Int64, Any, Any}

Returns the recurrence coefficients (a, b, c) for the k-th Legendre polynomial, corresponding to the three-term recurrence relation:

P_{k+1}(x) = (a_k x - b_k) P_k(x) - c_k P_{k-1}(x)

For Legendre polynomials:

  • a = 0
  • b = (2k + 1) / (k + 1)
  • c = k / (k + 1)
source

Hermite Polynomials (Normal Distribution)

For the weight function $\omega(y) = \exp(-y^2/2)/(2\pi)$ (normal distribution), the recurrence coefficients are $a_n = 1$, $b_n = 0$, and $c_n = n$. The first six Hermite polynomials are:

\[\begin{aligned} H_0 & = 0\\ H_1 & = 1\\ H_2 & = y\\ H_3 & = y^2 - 1\\ H_4 & = y^3 - 3y\\ H_5 & = y^4 - 6y^2 + 3\\ \end{aligned}\]

Their norms are given by

\[ \| H_n \|^2_\omega = n!\]

ExtendableASGFEM.HermitePolynomialsType
abstract type HermitePolynomials <: OrthogonalPolynomialType

Type for dispatching Hermite polynomials, which are orthogonal with respect to the standard normal distribution.

source
ExtendableASGFEM.distributionMethod
distribution(
    _::Type{HermitePolynomials}
) -> Distributions.Normal{Float64}

Returns the probability distribution associated with the Hermite polynomials, which is the standard normal distribution Normal(0, 1).

source
ExtendableASGFEM.normsMethod
norms(_::Type{HermitePolynomials}, k) -> Any

Returns the norm of the k-th Hermite polynomial, i.e.,

||H_k|| = sqrt(k!)
source
ExtendableASGFEM.recurrence_coefficientsMethod
recurrence_coefficients(
    _::Type{HermitePolynomials},
    k::Integer
) -> Tuple{Int64, Int64, Integer}

Returns the recurrence coefficients (a, b, c) for the k-th Hermite polynomial, corresponding to the three-term recurrence relation:

H_{k+1}(x) = (a_k x - b_k) H_k(x) - c_k H_{k-1}(x)

For Hermite polynomials:

  • a = 0
  • b = 1
  • c = k
source