Polynomial Families

The computational building blocks of polynomial Chaos expansions (PCE, see Polynomial Chaos Expansions) are methods for univariate orthogonal polynomial families.

Contains classes/methods for general univariate orthogonal polynomial families. - evaluation - gauss quadrature - ratio evaluations - linear/quadratic measure modifications

class UncertainSCI.opoly1d.OrthogonalPolynomialBasis1D(recurrence=[], probability_measure=True)
apply_jacobi_matrix(v)

Premultiplies the input array by the Jacobi matrix of the appropriate size for the polynomial family. Applies the Jacobi matrix across the first dimension of v.

Parameters:

v (ndarray) – Input vector or array

Returns:

Jv – J*v, where J is the Jacobi matrix of size v.shape[0].

Return type:

ndarray

canonical_connection(N)

Returns the N x N matrix C, where row n of C contains expansion coefficients for p_n in the monomial basis. I.e.,

p_n(x) = sum_{j=0}^{n} C[n,j] x**j,

for n = 0, …, N-1.

canonical_connection_inverse(N)

Returns the N x N matrix C, where row n of C contains expansion coefficients for x^n in the orthonormal basis . I.e.,

x^n = sum_{j=0}^{n} C[n,j] p_j(x)

for n = 0, …, N-1.

christoffel_function(x, k)

Computes the normalized (inverse) Christoffel function lambda, defined as

lambda**2 = k / sum(p**2, axi=1),

where p is a matrix containing evaluations of an orthonormal polynomial family up to degree k-1, defined by the recurrence coefficients ab.

derivative_expansion(s, N, K=None)

Computes the coefficients

\[\sigma^{(s)}_{n,k} = \left\langle p_n^{(s)}, p_k \right\rangle,\]

where \(p_n^{(s)}\) is the \(s`th derivative of the degree-:math:`k\) orthonormal polynomial \(p_k\). These are, equivalently, expansion coefficients of \(p_n^{(s)}\) in the basis \(\{p_k\}_k\).

Parameters:
  • s – The integer order of the derivative. Must be non-negative.

  • N – Computes coefficients for \(n \leq N\)

  • K – Computes coefficients for \(k \leq K\) (optional). If not given, is set to N.

Returns:

(N+1) x (K+1) numpy array containing coefficients.

Return type:

C

eval(x, n, d=0)
gauss_quadrature(N)

Computes the N-point Gauss quadrature rule associated to the recurrence coefficients ab.

gauss_radau_quadrature(N, anchor=0.0)

Computes the N-point Gauss quadrature rule associated to the polynomial family, with a node at the specified anchor.

jacobi_matrix(N)

Returns the N x N jacobi matrix associated to the polynomial family.

jacobi_matrix_driver(N)

Returns the N x N jacobi matrix associated to the input recurrence coefficients ab. (Requires ab.shape[0] >= N+1.)

leading_coefficient(N)

Returns the leading coefficients for the first N polynomial basis elements.

qpoly1d_eval(x, n, d=0)

Evalutes Christoffel-function normalized polynomials. These are given by

q_k(x) = p_k(x) / sqrt( sum_{j=0}^{n-1} p_j^2 ), k = 0, …, n-1

The output is a x.size x n array

r_eval(x, n, d=0)

Evalutes ratios of orthonormal polynomials. These are given by

r_n(x) = p_n(x) / p_{n-1}(x), n >= 1

The output is a x.size x n.size array.

recurrence(N)

Returns the first N+1 orthogonal polynomial recurrence pairs. The orthonormal polynomial family satisfies the recurrence

p_{-1}(x) = 0 p_0(x) = 1/ab[0,1]

x p_n(x) = ab[n+1,1] p_{n+1}(x) + ab[n+1,0] p_n(x) + ab[n,1] p_{n-1}(x)

(n >= 0)

The value ab[0,0] is ignored and never used.

Recurrence coefficients ab, once computed, are stored as an instance variable. On subsequent calls to this function, the stored values are returned if the instance variable already contains the desired coefficients. If the instance variable does not contain enough coefficients, then a call to recurrence_driver is performed to compute the desired coefficients, and the output is stored in the instance variable.

Parameters:

N (positive integer) – Maximum polynomial degree for desired recurrence coefficients

Returns:

ab – (N+1) x 2 array of recurrence coefficients.

Return type:

ndarray

recurrence_driver(N)
s_eval(x, n)

The output is a x.size x (n+1) array.

s_n(x) = p_n(x) / sqrt(sum_{j=0}^{n-1} p_j^2(x)), n >= 0

s_0(x) = p_0(x)

s_1(x) = 1 / b_1 * (x - a_1)

s_2(x) = 1 / (b_2 * sqrt(1+s_1^2)) * ((x - a_2)*s_1 - b_1)

Need {a_k, b_k} k up to n

tuple_product(N, alpha)

Computes integrals of polynomial products. Returns an N x N matrix C with entries

C[n,m] = < p_n p_m, p_alpha >,

where alpha is a vector of integers and p_alpha is defined

p_alpha = prod_{j=1}^{alpha.size} p_[alpha[j]](x),

The notation <., .> denotes the inner product under which the polynomial family is orthogonal.

Parameters:
  • N (integer) – Size of matrix to return

  • alpha (ndarray (1d)) – Multi-index defining a polynomial product

Returns:

C – Output N x N matrix containing integral values

Return type:

ndarray

tuple_product_generator(IC, ab=None)

Helper function that increments indices for a polynomial product expansion.

IC is a vector with entries

IC[j] = < p_j, p_alpha >, j in range(N),

where N = IC.size. The notation < ., .> is the inner product under which the polynomial family is orthonormal. alpha is a multi-index of arbitrary shape with a polynomial defined by

p_alpha = prod_{j in range(alpha.size)} p_{alpha[j]}(x).

The value of alpha is not needed by this function.

This function returns an N x N matrix C with entries

C[n,j] = < p_n p_j, p_alpha >, j, n in range(N)

Parameters:
  • IC (vector (1d array)) – Values of input inner products

  • ab (ndarray, optional) – Recurrence coefficients

Returns:

C – Output coefficient expansion vector for beta

Return type:

ndarray

Contains routines that specialize opoly1d things for classical orthogonal polynomial families - jacobi polys - hermite poly - laguerre polys

class UncertainSCI.families.JacobiPolynomials(alpha=0.0, beta=0.0, domain=[-1.0, 1.0], **options)

Jacobi Polynomial family.

eval_1d(x, n)

Evaluates univariate orthonormal polynomials given their three-term recurrence coefficients ab

eval_nd(x, lambdas)

Evaluates tensorial orthonormal polynomials associated with the univariate recurrence coefficients ab

fidistinv(u, n)
fidistinv_jacobi_setup(n, data)
idist(x, n, M=50)

Computes the order-n induced distribution at the locations x using M=10 quadrature nodes.

idist_medapprox(n)

Computes an approximation to the median of the degree-n induced distribution.

idistinv(u, n)
recurrence_driver(N)
class UncertainSCI.families.HermitePolynomials(alpha=2, rho=0.0, **options)
idist(x, n)
idistinv(u, n)
recurrence_driver(N)
class UncertainSCI.families.LaguerrePolynomials(alpha=1.0, rho=0.0, **options)
idist(x, n)
idist_medapprox(n)
idistinv(u, n)
recurrence_driver(N)
class UncertainSCI.families.DiscreteChebyshevPolynomials(M=2, domain=[0.0, 1.0])

Class for polynomials orthonormal on [0,1] with respect to an M-point discrete uniform measure with support equidistributed on the interval.

eval(x, n, **options)
fidistinv(u, n)

Fast routine for idistinv. (In this case, the “slow” routine is already very fast, so this is just an alias for idistinv.)

idist(x, n, nugget=False)

Evalutes the order-n induced distribution at the locations x.

Optionally, add a nugget to ensure correct computation on the support points.

idistinv(u, n)

Computes the inverse order-n induced distribution at the locations u.

recurrence_driver(N)