Spectral 1D

class SpectralToolbox.Spectral1D.Basis[source]

Bases: object

This is an abstract class for 1-d basis

Raises:NotImplementedError – if not overridden.
AssemblyDerivativeMatrix(x, N, k)[source]

Assemble the k-th derivative matrix using polynomials of order N.

Parameters:
  • x (ndarray [m]) – evaluate the polynomials
  • N (int) – maximum polynomial order
  • k (int) – order of the derivative
Returns:

(ndarray [m,``m``]) – derivative matrix

DerivativeMatrix(x, N, k)[source]
Evaluate(r, N)[source]

Evaluate the N-th order polynomial at points r

Parameters:
  • r (ndarray [m]) – evaluate the polynomial
  • N (int) – order of the polynomial
Returns:

(ndarray [m]) – polynomial evaluated on

the r points.

Raises:

NotImplementedError – if not overridden.

GradEvaluate(r, N, k=0)[source]

Evaluate the k-th derivative of the N-th order polynomial at points r

Parameters:
  • r (ndarray [m]) – evaluate the polynomial
  • N (int) – order of the polynomial
  • k (int) – order of the derivative
Returns:

(ndarray [m]) – k-th derivative of the

polynomial evaluated on the r points.

Raises:

NotImplementedError – if not overridden.

GradVandermonde(r, N, k=0)[source]

Generate the k-th derivative of the N-th order Vandermoned matrix.

Parameters:
  • r (ndarray [m]) – evaluate the polynomials
  • N (int) – maximum polynomial order
  • k (int) – order of the derivative
Returns:

(ndarray [m,``N+1``]) – polynomials evaluated

at the r points.

interpolate(x, f, xi, order)[source]

Interpolate function values f from points x to points xi.

Parameters:
  • x (ndarray [m]) –
  • f (ndarray [m]) –
  • xi (ndarray [n]) – points
  • order (int) – polynomial order
Returns:

(ndarray [n]) – interpolated function values at

points xi

project(r, f, N)[source]

Project f onto the span of the polynomial up to order N.

Parameters:
  • r (ndarray [m]) –
  • f (ndarray [m]) – evaluation points
  • N (int) – maximum polynomial order
Returns:

(ndarray [N]) – projection coefficients

class SpectralToolbox.Spectral1D.OrthogonalBasis(normalized=None)[source]

Bases: SpectralToolbox.Spectral1D.AbstractClasses.Basis

This is an abstract class for 1-d orthogonal basis

Parameters:normalized (bool) – whether to normalize the polynomials. Default=``None`` which leaves the choice at evaluation time.
Evaluate(r, N, norm=True)[source]

Evaluate the N-th order polynomial at points r

Parameters:
  • r (ndarray [m]) – evaluate the polynomial
  • N (int) – order of the polynomial
  • norm (bool) – whether to return normalized (True) or unnormalized (False) polynomial. The parameter is ignored if the normalized parameter is provided at construction time.
Returns:

(ndarray [m]) – polynomial evaluated on

the r points.

Raises:

NotImplementedError – if not overridden.

Gamma(N)[source]

Return the normalization constant for the N-th polynomial.

Parameters:N (int) – polynomial order
Returns:g – normalization constant
Return type:float
Raises:NotImplementedError – if not overridden.
GradEvaluate(r, N, k=0, norm=True)[source]

Evaluate the k-th derivative of the N-th order polynomial at points r

Parameters:
  • r (ndarray [m]) – evaluate the polynomial
  • N (int) – order of the polynomial
  • k (int) – order of the derivative
  • norm (bool) – whether to return normalized (True) or unnormalized (False) polynomial. The parameter is ignored if the normalized parameter is provided at construction time.
Returns:

(ndarray [m]) – k-th derivative of the

polynomial evaluated on the r points.

Raises:

NotImplementedError – if not overridden.

GradVandermonde(r, N, k=0, norm=True)[source]

Generate the k-th derivative of the N-th order Vandermoned matrix.

Parameters:
  • r (ndarray [m]) – evaluate the polynomials
  • N (int) – maximum polynomial order
  • k (int) – order of the derivative
  • norm (bool) – whether to return normalized (True) or unnormalized (False) polynomial. The parameter is ignored if the normalized parameter is provided at construction time.
Returns:

(ndarray [m,``N+1``]) – polynomials evaluated

at the r points.

Quadrature(N, quadType=None, norm=False)[source]

Generate quadrature rules of the selected type.

Parameters:
  • N (int) – polynomial order
  • quadType (string) – quadrature type
  • norm (bool) – whether to normalize the weights
Returns:

(tuple of ndarray) – format

(x,w) where x and w are N+1 nodes and weights

Raises:

NotImplementedError – if not overridden.

class SpectralToolbox.Spectral1D.OrthogonalPolynomial(normalized=None)[source]

Bases: SpectralToolbox.Spectral1D.AbstractClasses.OrthogonalBasis

This is an abstract class for 1-d polynomials

GaussLobattoQuadrature(N, norm=False)[source]

Generate Gauss Lobatto quadrature nodes and weights.

Parameters:
  • N (int) – polynomial order
  • norm (bool) – whether to normalize the weights
Returns:

(tuple of ndarray) – format

(x,w) where x and w are N+1 nodes and weights

Raises:

NotImplementedError – if not overridden.

GaussQuadrature(N, norm=False)[source]

Generate Gauss quadrature nodes and weights.

Parameters:
  • N (int) – polynomial order
  • norm (bool) – whether to normalize the weights
Returns:

(tuple of ndarray) – format

(x,w) where x and w are N+1 nodes and weights

Raises:

NotImplementedError – if not overridden.

GaussRadauQuadrature(N, norm=False)[source]

Generate Gauss Radau quadrature nodes and weights.

Parameters:
  • N (int) – polynomial order
  • norm (bool) – whether to normalize the weights
Returns:

(tuple of ndarray) – format

(x,w) where x and w are N+1 nodes and weights

Raises:

NotImplementedError – if not overridden.

Quadrature(N, quadType=None, norm=False)[source]

Generate quadrature rules of the selected type.

Parameters:
  • N (int) – polynomial order
  • quadType (string) – quadrature type (Guass, Gauss-Lobatto, Gauss-Radau)
  • norm (bool) – whether to normalize the weights
Returns:

(tuple of ndarray) – format

(x,w) where x and w are N+1 nodes and weights

Raises:

ValueError – if quadrature type is not available

RecursionCoeffs(N)[source]

Generate the first N recursion coefficients.

These coefficients define the polynomials: \(P_{n+1}(x) = (x - a_n) P_n(x) - b_n P_{n-1}(x)\)

Parameters:N (int) – number of recursion coefficients
Returns:
(tuple [2] of ndarray [N]) –
recursion coefficients a and b
Raises:NotImplementedError – if not overridden.
class SpectralToolbox.Spectral1D.JacobiPolynomial(alpha, beta, span=None, normalized=None)[source]

Bases: SpectralToolbox.Spectral1D.AbstractClasses.OrthogonalPolynomial

Construction of Jacobi polynomials

Parameters:
  • alpha (float) – parameter \(\alpha > -1\)
  • beta (float) – parameter \(\beta > -1\).
  • span (list) – span of the domain
  • normalized (bool) – whether to normalize the polynomials. Default=``None`` which leaves the choice at evaluation time.
Evaluate(x, N, alpha=None, beta=None, norm=True)[source]

Evaluate the N-th order Jacobi polynomial

Parameters:
  • x ((ndarray [m]) – evaluate the polynomial
  • N (int) – order of the polynomial
Returns:

(ndarray [m]) – polynomial evaluated on the

x points.

FastChebyshevTransform(f)[source]

FastChebyshevTransform(): Returns the coefficients of the Fast Chebyshev Transform.

Syntax:
fhat = FastChebyshevTransform(f)
Input:
  • f = (1d-array,float) function values
Output:
  • fhat = (1d-array,float) list of Polynomial coefficients

Warning

It is assumed that the values f are computed at Chebyshev-Gauss-Lobatto points.

Note

If f is odd, the vector is interpolated to even Chebyshev-Gauss-Lobatto points.

Note

Modification of algorithm (29) from [Kop09]

Gamma(N, alpha=None, beta=None)[source]

Return the normalization constant for the N-th Jacobi polynomial.

See also

OrthogonalPolynomial.Gamma()

GaussLobattoQuadrature(N, norm=False)[source]

Jacobi Gauss Lobatto quadrature points

GaussQuadrature(N, norm=False)[source]

Jacobi Gauss quadrature points

GaussRadauQuadrature(N, norm=False)[source]

Jacobi Gauss Radau quadrature points

GradEvaluate(r, N, k=0, norm=True)[source]

k-th derivative of the N-th order Jacobi polynomial.

See also

OrthogonalPolynomial.GradEvaluate()

InverseFastChebyshevTransform(fhat)[source]

InverseFastChebyshevTransform(): Returns the coefficients of the Inverse Fast Chebyshev Transform.

Syntax:
f = InverseFastChebyshevTransform(fhat)
Input:
  • fhat = (1d-array,float) list of Polynomial coefficients
Output:
  • f = (1d-array,float) function values

Note

If f is odd, the vector is padded with a zero value (highest freq.)

Note

Modification of algorithm (29) from [Kop09]

RecursionCoeffs(N, alpha=None, beta=None)[source]

Get coefficients for the Jacobi polynomials up to order N

eps = 1.1102230246251565e-16
class SpectralToolbox.Spectral1D.LaguerrePolynomial(alpha, normalized=None)[source]

Bases: SpectralToolbox.Spectral1D.AbstractClasses.OrthogonalPolynomial

Construction of Laguerre polynomials

Parameters:
  • alpha (float) – parameter \(\alpha\).
  • normalized (bool) – whether to normalize the polynomials. Default=``None`` which leaves the choice at evaluation time.
Evaluate(x, N, alpha=None, norm=True)[source]

Evaluate the N-th order Laguerre polynomial

Parameters:
  • x ((ndarray [m]) – evaluate the polynomial
  • N (int) – order of the polynomial
Returns:

(ndarray [m]) – polynomial evaluated on the

x points.

Gamma(N, alpha=None, beta=None)[source]

Return the normalization constant for the N-th Laguerre polynomial.

See also

OrthogonalPolynomial.Gamma()

GaussQuadrature(N, norm=False)[source]

Laguerre Gauss quadrature points

GaussRadauQuadrature(N, norm=False)[source]

Laguerre Gauss Radau quadrature points

GradEvaluate(r, N, k=0, norm=True)[source]

k-th derivative of the N-th order Laguerre polynomial.

See also

OrthogonalPolynomial.GradEvaluate()

RecursionCoeffs(N, alpha=None)[source]

Get coefficients for the Laguerre polynomials up to order N

eps = 1.1102230246251565e-16
class SpectralToolbox.Spectral1D.HermitePhysicistsPolynomial(normalized=None)[source]

Bases: SpectralToolbox.Spectral1D.AbstractClasses.OrthogonalPolynomial

Construction of the Hermite Physicists’ polynomials

Parameters:normalized (bool) – whether to normalize the polynomials. Default=``None`` which leaves the choice at evaluation time.
Evaluate(x, N, norm=True)[source]

Evaluate the N-th order Hermite Physicists’ polynomial

Parameters:
  • x ((ndarray [m]) – evaluate the polynomial
  • N (int) – order of the polynomial
Returns:

(ndarray [m]) – polynomial evaluated on the

x points.

Gamma(N, alpha=None, beta=None)[source]

Return the normalization constant for the N-th Hermite Physicists’ polynomial.

See also

OrthogonalPolynomial.Gamma()

GaussQuadrature(N, norm=False)[source]

Hermite Physicists’ Gauss quadrature points

GradEvaluate(r, N, k=0, norm=True)[source]

k-th derivative of the N-th order Hermite Physicists’ polynomial.

See also

OrthogonalPolynomial.GradEvaluate()

RecursionCoeffs(N)[source]

Get coefficients for the Hermite Physicists’ polynomials up to order N

class SpectralToolbox.Spectral1D.HermiteProbabilistsPolynomial(normalized=None)[source]

Bases: SpectralToolbox.Spectral1D.AbstractClasses.OrthogonalPolynomial

Construction of the Hermite Probabilists polynomials

Parameters:normalized (bool) – whether to normalize the polynomials. Default=``None`` which leaves the choice at evaluation time.
Evaluate(x, N, norm=True)[source]

Evaluate the N-th order Hermite Probabilists’ polynomial

Parameters:
  • x ((ndarray [m]) – evaluate the polynomial
  • N (int) – order of the polynomial
Returns:

(ndarray [m]) – polynomial evaluated on the

x points.

Gamma(N, alpha=None, beta=None)[source]

Return the normalization constant for the N-th Hermite Probabilists’ polynomial.

See also

OrthogonalPolynomial.Gamma()

GaussQuadrature(N, norm=False)[source]

Hermite Probabilists’ Gauss quadrature points

GradEvaluate(r, N, k=0, norm=True)[source]

k-th derivative of the N-th order Hermite Probabilists’ polynomial.

See also

OrthogonalPolynomial.GradEvaluate()

RecursionCoeffs(N)[source]

Get coefficients for the Hermite Probabilists’ polynomials up to order N

class SpectralToolbox.Spectral1D.GenericOrthogonalPolynomial(mu, endl, endr, mc=1, mp=0, xp=None, yp=None, ncapm=500, irout=1, normalized=None)[source]

Bases: SpectralToolbox.Spectral1D.AbstractClasses.OrthogonalPolynomial

Construction of polynomials orthogonal with respect to a generic measure

Parameters:
  • mu (float (float x, int i) – function returning the mass value at the point x of interval i of the continuous spectrum.
  • endl (ndarray [mc]) – left endpoints of the intervals in the continuous spectrum.
  • endr (ndarray [mc]) – right endpoints of the intervals in the continuous spectrum.
  • mc (int) – number of component intervals in the continuous part of the spectrum.
  • mp (int) – number of points in the discrete part of the spectrum. If the measure has no discrete part, set mp=0.
  • xp (ndarray [mp]) – abscissas of the points in the discrete spectrum.
  • yp (ndarray [mp]) – jumps of the points in the discrete spectrum.
  • ncapm (int) – maximum integer N0
  • irout (int) – selects the routine for generating the recursion coefficients from the discrete inner product; irout=1 selects the routine sti, irout!=1 selects the routine lancz.
  • normalized (bool) – whether to normalize the polynomials. Default=``None`` which leaves the choice at evaluation time.

Note

Parameters iq, quad, idelta in [Gau94] are suppressed. Instead the routine qgp of ORTHPOL [Gau94] is used by default (iq=0 and idelta=2)

Evaluate(x, N, norm=True)[source]

Evaluate the N-th order polynomial

See also

OrthogonalPolynomial.Evaluate()

Gamma(N, alpha=None, beta=None)[source]

Return the normalization constant for the N-th polynomial.

See also

OrthogonalPolynomial.Gamma()

GaussLobattoQuadrature(N, norm=False)[source]

Gauss Lobatto quadrature points.

GaussQuadrature(N, norm=False)[source]

Gauss quadrature points.

GaussRadauQuadrature(N, norm=False)[source]

Gauss Radau quadrature points.

GradEvaluate(r, N, k=0, norm=True)[source]

k-th derivative of the N-th order polynomial.

Warning

no derivatives are implemented for this type of polynomials. Works only with k==0.

See also

OrthogonalPolynomial.GradEvaluate()

GradVandermonde(r, N, k=0, norm=True)[source]

See also

OrthogonalPolynomial.GradVandermonde()

RecursionCoeffs(N)[source]

Get the recursion coefficients up to order N

eps = 1.1102230246251565e-16
class SpectralToolbox.Spectral1D.HermitePhysicistsFunction(normalized=None)[source]

Bases: SpectralToolbox.Spectral1D.AbstractClasses.OrthogonalPolynomial

Construction of the Hermite Physiticists’ functions

Parameters:normalized (bool) – whether to normalize the polynomials. Default=``None`` which leaves the choice at evaluation time.
Evaluate(x, N, norm=True)[source]

Evaluate the N-th order Hermite Physicists’ function

See also

OrthogonalPolynomial.Evaluate()

Gamma(N)[source]

Return the normalization constant for the N-th Hermite Physiticsts’ function.

See also

OrthogonalPolynomial.Gamma()

GaussQuadrature(N, norm=False)[source]

Hermite Physicists’ function Gauss quadratures

GradEvaluate(x, N, k=0, norm=True)[source]

Evaluate the k-th derivative of the N-th order Hermite Physicists’ function

One can write it as

\[\psi_n^{(k)} = \sum_{i=0}^k (-1)^i F_{k,i} \psi_i H_n^{(k-i)}\]

where \(F_{k,i}\) is the \(i\)-th component of the \(k\)-th row of the Pascal’s triangle, \(psi_i\) is the \(i\)-th Hermite Physicists’ function and \(H_n^{(n-i)}\) is the \((k-1)\)-th derivative of the \(n\)-th Hermite Physicists’ polynomial.

See also

OrthogonalPolynomial.GradEvaluate()

class SpectralToolbox.Spectral1D.HermiteProbabilistsFunction(normalized=None)[source]

Bases: SpectralToolbox.Spectral1D.AbstractClasses.OrthogonalPolynomial

Construction of the Hermite Probabilists’ functions

This are rescaling of \(\sqrt{2}\) of the Hermite Physicists’ functions.

Parameters:normalized (bool) – whether to normalize the polynomials. Default=``None`` which leaves the choice at evaluation time.
Evaluate(x, N, norm=True)[source]

Evaluate the N-th order Hermite Probabilists’ function

See also

OrthogonalPolynomial.Evaluate()

Gamma(N)[source]

Return the normalization constant for the N-th Hermite Probabilists’ function.

See also

OrthogonalPolynomial.Gamma()

GaussQuadrature(N, norm=False)[source]

Hermite Probabilists’ function Gauss quadratures

GradEvaluate(x, N, k=0, norm=True)[source]

Evaluate the k-th derivative of the N-th order Hermite Probabilists’ function

class SpectralToolbox.Spectral1D.LaguerreFunction(alpha, normalized=None)[source]

Bases: SpectralToolbox.Spectral1D.AbstractClasses.OrthogonalPolynomial

Construction of the Laguerre functions

Parameters:
  • alpha (float) – parameter \(\alpha\).
  • normalized (bool) – whether to normalize the polynomials. Default=``None`` which leaves the choice at evaluation time.
Evaluate(x, N, norm=True)[source]

Evaluate the N-th order Laguerre function

See also

OrthogonalPolynomial.Evaluate()

Gamma(N)[source]
GaussQuadrature(N, norm=False)[source]

Laguerre function Gauss quadratures

GaussRadauQuadrature(N, norm=False)[source]
GradEvaluate(x, N, k=0, norm=True)[source]

Evaluate the k-th derivative of the N-th order Laguerre function

See also

OrthogonalPolynomial.GradEvaluate()

class SpectralToolbox.Spectral1D.Fourier[source]

Bases: SpectralToolbox.Spectral1D.AbstractClasses.OrthogonalBasis

Evaluate(x, N)[source]

Evaluate the N-th Fourier basis

See also

OrthogonalBase.Evaluate()

Gamma(N)[source]

Evaluate the N-th normalization constant for the Fourier basis

See also

OrthogonalBase.Gamma()

GradEvaluate(x, N, k=0, norm=True)[source]

Evaluate the k-th derivative of the N-th Fourier basis

See also

OrthogonalBase.GradEvaluate()

Quadrature(N, norm=False)[source]

Generate quadrature points for the Fourier series

See also

OrthogonalBase.Quadrature()

interpolate(x, f, xi, order)[source]

Interpolate function values f from points x to points xi.

Parameters:
  • x (ndarray [m]) –
  • f (ndarray [m]) –
  • xi (ndarray [n]) – points
  • order (int) – polynomial order
Returns:

(ndarray [n]) – interpolated function values at

points xi

project(r, f, N)[source]

Project f onto the span of Fourier basis up to order N

See also

OrthogonalBase.project()

SpectralToolbox.Spectral1D.FirstPolynomialDerivativeMatrix(x)[source]

PolynomialDerivativeMatrix(): Assemble the first Polynomial Derivative Matrix using matrix multiplication.

Syntax:
D = FirstPolynomialDerivativeMatrix(x)
Input:
  • x = (1d-array,float) set of points on which to evaluate the derivative matrix
Output:
  • D = derivative matrix

Notes

Algorithm (37) from [Kop09]

SpectralToolbox.Spectral1D.PolynomialDerivativeMatrix(x, k)[source]

PolynomialDerivativeMatrix(): Assemble the Polynomial k-th Derivative Matrix using the matrix recursion. This algorithm is generic for every types of polynomials.

Syntax:
D = PolynomialDerivativeMatrix(x,k)
Input:
  • x = (1d-array,float) set of points on which to evaluate the derivative matrix
  • k = derivative order
Output:
  • D = k-th derivative matrix

Notes

Algorithm (38) taken from [Kop09]

SpectralToolbox.Spectral1D.BarycentricWeights(x)[source]

BarycentricWeights(): Returns a 1-d array of weights for Lagrange Interpolation

Syntax:
w = BarycentricWeights(x)
Input:
  • x = (1d-array,float) set of points
Output:
  • w = (1d-array,float) set of barycentric weights

Notes

Algorithm (30) from [Kop09]

SpectralToolbox.Spectral1D.LagrangeInterpolationMatrix(x, w, xi)[source]

LagrangeInterpolationMatrix(): constructs the Lagrange Interpolation Matrix from points x to points xi

Syntax:
T = LagrangeInterpolationMatrix(x, w, xi)
Input:
  • x = (1d-array,float) set of N original points
  • w = (1d-array,float) set of N barycentric weights
  • xi = (1d-array,float) set of M interpolating points
Output:
  • T = (2d-array(MxN),float) Lagrange Interpolation Matrix

Notes

Algorithm (32) from [Kop09]

SpectralToolbox.Spectral1D.LagrangeInterpolate(x, f, xi)[source]

LagrangeInterpolate(): Interpolate function values f from points x to points xi using Lagrange weights

Syntax:
fi = LagrangeInterpolate(x, f, xi)
Input:
  • x = (1d-array,float) set of N original points where f is evaluated
  • f = (1d-array/2d-array,float) set of N function values (if K functions are passed, the values are stored in a NxK matrix)
  • xi = (1d-array,float) set of M points where the function is interpolated
Output:
  • fi = (1d-array,float) set of M function values (if K functions are passed, the values are stored in a MxK matrix)

Notes

Modification of Algorithm (33) from [Kop09]

SpectralToolbox.Spectral1D.LinearShapeFunction(x, xm, xp, xi)[source]

Hat function used for linear interpolation

Parameters:
  • x (array) – 1d original points
  • xm,xp (float) – bounding points of the support of the shape function
  • xi (array) – 1d interpolation points
Returns array N:
 

evaluation of the shape function on xi

SpectralToolbox.Spectral1D.SparseLinearShapeFunction(x, xm, xp, xi)[source]

Hat function used for linear interpolation. Returns sparse indices for construction of scipy.sparse.coo_matrix.

Parameters:
  • x (array) – 1d original points
  • xm,xp (float) – bounding points of the support of the shape function
  • xi (array) – 1d interpolation points
Returns tuple (idxs,vals):
 

List of indexes and evaluation of the shape function on xi

SpectralToolbox.Spectral1D.LinearInterpolationMatrix(x, xi)[source]

LinearInterpolationMatrix(): constructs the Linear Interpolation Matrix from points x to points xi

Syntax:
T = LagrangeInterpolationMatrix(x, xi)
Input:
  • x = (1d-array,float) set of N original points
  • xi = (1d-array,float) set of M interpolating points
Output:
  • T = (2d-array(MxN),float) Linear Interpolation Matrix
SpectralToolbox.Spectral1D.SparseLinearInterpolationMatrix(x, xi)[source]

LinearInterpolationMatrix(): constructs the Linear Interpolation Matrix from points x to points xi. Returns a scipy.sparse.coo_matrix

Syntax:
T = LagrangeInterpolationMatrix(x, xi)
Input:
  • x = (1d-array,float) set of N original points
  • xi = (1d-array,float) set of M interpolating points
Output:
  • T = (scipy.sparse.coo_matrix(MxN),float) Linear Interpolation Matrix
class SpectralToolbox.Spectral1D.ConstantExtendedHermiteProbabilistsFunction(normalized=None)[source]

Bases: SpectralToolbox.Spectral1D.AbstractClasses.Basis

Construction of the Hermite Probabilists’ functions extended with the constant basis

The basis is defined by:

\[\phi_0(x) = 1 \qquad \phi_i(x) = \psi_{i-1}(x) \quad \text{for } i=1\ldots\]

where \(\psi_j\) are the Hermite Probabilists’ functions.

Parameters:normalized (bool) – whether to normalize the underlying polynomials. Default=``None`` which leaves the choice at evaluation time.
Evaluate(x, N, norm=True)[source]

Evaluate the N-th order constant extended Hermite Probabilists’ function

See also

OrthogonalPolynomial.Evaluate()

GaussQuadrature(N, norm=False)[source]

Hermite Probabilists’ function Gauss quadratures

GradEvaluate(x, N, k=0, norm=True)[source]

Evaluate the k-th derivative of the N-th order constant extended Hermite Probabilists’ function

GradVandermonde(r, N, k=0, norm=True)[source]

Generate the k-th derivative of the N-th order Vandermoned matrix.

Parameters:
  • r (ndarray [m]) – evaluate the polynomials
  • N (int) – maximum polynomial order
  • k (int) – order of the derivative
  • norm (bool) – whether to return normalized (True) or unnormalized (False) polynomial. The parameter is ignored if the normalized parameter is provided at construction time.
Returns:

(ndarray [m,``N+1``]) – polynomials evaluated

at the r points.

Quadrature(N, quadType=None, norm=False)[source]

Generate quadrature rules of the selected type.

class SpectralToolbox.Spectral1D.HermiteProbabilistsRadialBasisFunction(order, scale=1.0)[source]

Bases: SpectralToolbox.Spectral1D.AbstractClasses.Basis

Construction of the Hermite Probabilists’ Radial Basis Functions

For the set \(\left\{x_i\right\}_{i=1}^N\) of Gauss-Hermite points, the basis are defined by:

\[\phi_i(x) = \exp\left( -\frac{(x-x_i)^2}{2\sigma^2_{i}} \right)\]

where \(\sigma_i=(x_{i+1} - x_{i-1})/2.\), \(\sigma_0=\sigma_1\) and \(\sigma_N=\sigma_{N-1}\)

Parameters:
  • nknots (int) – number of knots points \(x_i\)
  • scale (float) – scaling for the badwidth \(\sigma\).
Evaluate(x, N, norm=True, extended_output=False)[source]

Evaluate the N-th Hermite Probabilists’ Radial Basis Function

See also

OrthogonalPolynomial.Evaluate()

GaussQuadrature(N, norm=False)[source]

Hermite Probabilists’ function Gauss quadratures

GradEvaluate(x, N, k=0, norm=True)[source]

Evaluate the k-th derivative of the N-th Hermite Probabilists’ Radial Basis Function

See also

OrthogonalPolynomial.GradEvaluate()

GradVandermonde(r, N, k=0, norm=True)[source]

Generate the k-th derivative of the N-th order Vandermoned matrix.

Parameters:
  • r (ndarray [m]) – evaluate the polynomials
  • N (int) – maximum polynomial order
  • k (int) – order of the derivative
  • norm (bool) – whether to return normalized (True) or unnormalized (False) polynomial. The parameter is ignored if the normalized parameter is provided at construction time.
Returns:

(ndarray [m,``N+1``]) – polynomials evaluated

at the r points.

Quadrature(N, quadType=None, norm=False)[source]

Generate quadrature rules of the selected type.

class SpectralToolbox.Spectral1D.ConstantExtendedHermiteProbabilistsRadialBasisFunction(nbasis, scale=1.0)[source]

Bases: SpectralToolbox.Spectral1D.AbstractClasses.Basis

Construction of the Hermite Probabilists’ Radial Basis Functions

For the set \(\left\{x_i\right\}_{i=1}^N\) of Gauss-Hermite points, the basis \(\{\phi_i\}_{i=0}^M\) are defined by:

\[\begin{split}\phi_0(x) = 1 \\ \phi_i(x) = \exp\left( -\frac{(x-x_i)^2}{2\sigma^2_{i-1}} \right)\end{split}\]

where \(\sigma_i=x_{i+1} - x_{i}\), \(\sigma_0=\sigma_1\) and \(\sigma_N=\sigma_{N-1}\)

Parameters:
  • order (int) – maximum order \(M\)
  • scale (float) – scaling for the badwidth \(\sigma\).
Evaluate(x, N, norm=True)[source]

Evaluate the N-th Hermite Probabilists’ Radial Basis Function

See also

OrthogonalPolynomial.Evaluate()

GaussQuadrature(N, norm=False)[source]

Hermite Probabilists’ function Gauss quadratures

GradEvaluate(x, N, k=0, norm=True)[source]

Evaluate the k-th derivative of the N-th Hermite Probabilists’ Radial Basis Function

See also

OrthogonalPolynomial.GradEvaluate()

GradVandermonde(r, N, k=0, norm=True)[source]

Generate the k-th derivative of the N-th order Vandermoned matrix.

Parameters:
  • r (ndarray [m]) – evaluate the polynomials
  • N (int) – maximum polynomial order
  • k (int) – order of the derivative
  • norm (bool) – whether to return normalized (True) or unnormalized (False) polynomial. The parameter is ignored if the normalized parameter is provided at construction time.
Returns:

(ndarray [m,``N+1``]) – polynomials evaluated

at the r points.

Quadrature(N, quadType=None, norm=False)[source]

Generate quadrature rules of the selected type.

class SpectralToolbox.Spectral1D.LinearExtendedHermiteProbabilistsRadialBasisFunction(nbasis, scale=1.0)[source]

Bases: SpectralToolbox.Spectral1D.AbstractClasses.Basis

Construction of the Hermite Probabilists’ Radial Basis Functions

For the set \(\left\{x_i\right\}_{i=1}^N\) of Gauss-Hermite points, the basis \(\{\phi_i\}_{i=0}^M\) are defined by:

\[\begin{split}\phi_0(x) = 1 \\ \phi_1(x) = x \\ \phi_i(x) = \exp\left( -\frac{(x-x_i)^2}{2\sigma^2_{i-1}} \right)\end{split}\]

where \(\sigma_i=x_{i+1} - x_{i}\), \(\sigma_0=\sigma_1\) and \(\sigma_N=\sigma_{N-1}\)

Parameters:
  • order (int) – maximum order \(M\)
  • scale (float) – scaling for the badwidth \(\sigma\).
Evaluate(x, N, norm=True, extended_output=False)[source]

Evaluate the N-th Hermite Probabilists’ Radial Basis Function

See also

OrthogonalPolynomial.Evaluate()

GaussQuadrature(N, norm=False)[source]

Hermite Probabilists’ function Gauss quadratures

GradEvaluate(x, N, k=0, norm=True)[source]

Evaluate the k-th derivative of the N-th Hermite Probabilists’ Radial Basis Function

See also

OrthogonalPolynomial.GradEvaluate()

GradVandermonde(r, N, k=0, norm=True)[source]

Generate the k-th derivative of the N-th order Vandermoned matrix.

Parameters:
  • r (ndarray [m]) – evaluate the polynomials
  • N (int) – maximum polynomial order
  • k (int) – order of the derivative
  • norm (bool) – whether to return normalized (True) or unnormalized (False) polynomial. The parameter is ignored if the normalized parameter is provided at construction time.
Returns:

(ndarray [m,``N+1``]) – polynomials evaluated

at the r points.

Quadrature(N, quadType=None, norm=False)[source]

Generate quadrature rules of the selected type.

class SpectralToolbox.Spectral1D.Poly1D(poly, params, sdout=<_io.TextIOWrapper name='<stderr>' mode='w' encoding='UTF-8'>)[source]

Bases: SpectralToolbox.Spectral1D.AbstractClasses.OrthogonalPolynomial

Initialization of the Polynomial instance.

This method generates an instance of the Poly1D class, to be used in order to generate orthogonal basis of the polynomial type selected. Avaliable polynomial types can be selected using their string name or by predefined attributes

  • Jacobi Polynomials: ‘Jacobi’ or Spectral1D.JACOBI
  • Hermite Physicist: ‘HermiteP’ or Spectral1D.HERMITEP
  • Hermite Function: ‘HermiteF’ or Spectral1D.HERMITEF
  • Hermite Probabilistic: ‘HermitePprob’ or Spectral1D.HERMITEP_PROB
  • Hermite Function Probabilistic: ‘HermiteFprob’ or Spectral1D.HERMITEF_PROB
  • Laguerre Polynomial: ‘LaguerreP’ or Spectral1D.LAGUERREP
  • Laguerre Function: ‘LaguerreF’ or Spectral1D.LAGUERREF
  • General orthogonal polynomial: ‘ORTHPOL’ or Spectral1D.ORTHPOL
  • Fourier: ‘Fourier’ or Spectral1D.FOURIER

Additional parameters are required for some polynomials.

Polynomial Parameters
Jacobi (alpha,beta)
HermiteP None
HermiteF None
HermitePprob None
LaguerreP alpha
LaguerreF alpha
ORTHPOL see notes
Fourier None
Available quadrature rules (related to selected polynomials):
  • Gauss or Spectral1D.GAUSS
  • Gauss-Lobatto or Spectral1D.GAUSSLOBATTO
  • Gauss-Radau or Spectral1D.GAUSSRADAU
Available quadrature rules (without polynomial selection):
  • Kronrod-Patterson on the real line or Spectral1D.KPN (function Spectral1D.kpn(n))
  • Kronrod-Patterson uniform or Spectral1D.KPU (function Spectral1D.kpu(n))
  • Clenshaw-Curtis or Spectral1D.CC (function Spectral1D.cc(n))
  • Fejer’s or Spectral1D.FEJ (function Spectral1D.fej(n))
Parameters:
  • poly (string) – The orthogonal polynomial type desired
  • params (list) – The parameters needed by the selected polynomial
  • sdout (stream) – output stream for logging

Note

The ORTHPOL polynomials are built up using the “Multiple-Component Discretization Procedure” described in [Gau94]. The following parameters describing the measure function are required in order to use the procedure for finding the recursion coefficients (alpha,beta) and have to be provided at construction time:

  • 0 ncapm: (int) maximum integer N0 (default = 500)
  • 1 mc: (int) number of component intervals in the continuous part of the spectrum
  • 2 mp: (int) number of points in the discrete part of the spectrum. If the measure has no discrete part, set mp=0
  • 3 xp, 4 yp: (Numpy 1d-array) of dimension mp, containing the abscissas and the jumps of the point spectrum
  • 5 mu: (function) measure function that returns the mass (float) with arguments: x (float) absissa, i (int) interval number in the continuous part
  • 6 irout: (int) selects the routine for generating the recursion coefficients from the discrete inner product; irout=1 selects the routine sti, irout!=1 selects the routine lancz
  • 7 finl, 8 finr: (bool) specify whether the extreme left/right interval is finite (false for infinite)
  • 9 endl, 10 endr: (Numpy 1d-array) of dimension mc containing the left and right endpoints of the component intervals. If the first of these extends to -infinity, endl[0] is not being used by the routine.

Parameters iq, quad, idelta in [Gau94] are suppressed. Instead the routine qgp of ORTHPOL [Gau94] is used by default (iq=0 and idelta=2)

Deprecated since version 0.2.0: Use OrthogonalPolynomial and its sub-classes instead.

Raises:DeprecationWarnin – use OrthogonalPolynomial and its sub-classes instead.
AssemblyDerivativeMatrix(x, N, k)[source]

AssemblyDerivativeMatrix(): Assemble the k-th derivative matrix using polynomials of order N.

Syntax:
Dk = AssemblyDerivativeMatrix(x,N,k)
Input:
  • x = (1d-array,float) Set of points on which to evaluate the polynomials
  • N = (int) maximum order in the vanermonde matrix
  • k = (int) derivative order
Output:
  • Dk = Derivative matrix
Description:
This function performs D = linalg.solve(V.T, Vx.T) where V and Vx are a Generalized Vandermonde Matrix and its derivative respectively.

Notes

For Chebyshev Polynomial, this function refers to the recursion form implemented in PolynomialDerivativeMatrix

ChebyshevDerivativeCoefficients(fhat)[source]

ChebyshevDerivativeCoefficients(): computes the Chebyshev coefficients of the derivative of a function

Syntax:
dfhat = ChebyshevDerivativeCoefficients(fhat)
Input:
  • fhat = (1d-array,float) list of Chebyshev coefficients of the original function
Output:
  • dfhat = (1d-array,float) list of Chebyshev coefficients of the derivative of the original function

Notes

Algorithm (5) from [Kop09]

DiscretePolynomialTransform(r, f, N)[source]

Computes the Discrete Polynomial Transform of function values f, i.e. project on the space spanned by the selected polynomials up to order N.

param r:set of points on which to the polynomials are evaluated
type r:1d-array or float
param f:function values
type f:1d-array or float
param int N:maximum order in the generalized vanermonde matrix
return:projection coefficients
rtype:1d-ndarray containing the projection coefficients. If len(r) == N+1 exact projection is used. If len(r) != N+1 the least square solution is obtained by numpy.linalg.lstsq.

Note

If the Chebyshev polynomials are chosen and r contains Chebyshev-Gauss-Lobatto points, the Fast Chebyshev Transform is used. Otherwise uses the Generalized Vandermonde Matrix in order to transform from physical space to transform space.

Note

If the Fourier expansion is chosen, no check on r and N is done. Though the periodic function is assumed to be evaluated over N+1 == len(f) equidistant points over [0, 2pi-2pi/(N+2)]. The returned coefficients a are ordered such that the following expansion is obtained: \(f(x)=a_0 + \sum_{i=1}^M a_{2i-1} cos( (2i-1) x) - \sum_{i=1}^M a_{2i} sin( 2i x )\), where :math:`M=lfloor N/2

floor + 1` if N is even, \(M=(N+1)/2\) if N is odd.

Note

this is the same of calling :py:method:`DiscretePolynomialTransform`

See also

FastChebyshevTransform

See also

numpy.linalg.lstsq

FastChebyshevTransform(f)[source]

FastChebyshevTransform(): Returns the coefficients of the Fast Chebyshev Transform.

Syntax:
fhat = FastChebyshevTransform(f)
Input:
  • f = (1d-array,float) function values
Output:
  • fhat = (1d-array,float) list of Polynomial coefficients

Warning

It is assumed that the values f are computed at Chebyshev-Gauss-Lobatto points.

Note

If f is odd, the vector is interpolated to even Chebyshev-Gauss-Lobatto points.

Note

Modification of algorithm (29) from [Kop09]

Gamma(N)[source]

Gamma(): returns the normalization constant for the N-th polynomial

Syntax:
g = Gamma(N)
Input:
  • N = polynomial order
Output:
  • g = normalization constant
GaussLobattoQuadrature(N, normed=False, left=None, right=None)[source]

GaussLobattoQuadrature(): Generates list of nodes for the Gauss-Lobatto quadrature rule using selected Polynomial basis

Syntax:
x = GaussLobattoQuadrature(N,[normed=False],[left=None],[right=None])
Input:
  • N = (int) accuracy level required
  • normed = (optional,bool) whether the weights will be normalized or not
  • left = (optional,float) containing the left endpoint (used by ORTHPOL)
  • right = (optional,float) containing the right endpoint (used by ORTHPOL)
Output:
  • x = (1d-array,float) containing the nodes
  • w = (1d-array,float) containing the weights

Note

Available only for Jacobi Polynomials and ORTHPOL

GaussQuadrature(N, normed=False)[source]

GaussQuadrature(): Generates list of nodes and weights for the Gauss quadrature rule using the selected Polynomial basis

Syntax:
(x,w) = GaussQuadrature(N,[normed=False])
Input:
  • N = (int) accuracy level required
  • normed = (optional,bool) whether the weights will be normalized or not
Output:
  • x = (1d-array,float) containing the nodes
  • w = (1d-array,float) containing the weights
GaussRadauQuadrature(N, normed=False, end=None)[source]

GaussRadauQuadrature(): Generates list of nodes for the Gauss-Radau quadrature rule using selected Polynomial basis

Syntax:
``x = GaussRadauQuadrature(N,[normed=False],[end=None])’‘
Input:
  • ``N’’ = (int) accuracy level required
  • normed = (optional,bool) whether the weights will be normalized or not
  • end = (optional,float) containing the endpoint (used by ORTHPOL)
Output:
  • ``x’’ = (1d-array,float) containing the nodes
  • ``w’’ = (1d-array,float) weights

Note

Available only for Laguerre Polynomials/Functions and ORTHPOL

GradEvaluate(r, N, k=0, norm=True, isConstructingVandermonde=False, recursion_storage=None)[source]

Evaluate the k-th derivative of the N-th order polynomial at points r

Parameters:
  • r (1d-array or float) – set of points on which to evaluate the polynomial
  • N (int) – order of the polynomial
  • k (int) – order of the derivative [default=0]
  • norm (bool) – whether to return normalized (True) or non normalized (False) polynomials
  • isConstructingVandermonde (bool) – whether to store values in the recursion relation for the construciton of the Vandermonde matrix
  • recursion_storage (ndarray len(r) x 2) – storage needed for the recursion whether one is computing the Vandermonde matrix (isConstructingVandermonde==True)
Returns:

polynomial evaluated on r

Return type:

1d-ndarray

GradVandermonde(r, N, k=0, norm=True)[source]
GradVandermonde1D(r, N, k=0, norm=True)[source]

GradVandermonde1D(): Initialize the k-th gradient of the modal basis N at r

Syntax:
V = GradVandermonde1D(r,N,k,[norm])
Input:
  • r = (1d-array,float) set of M points on which to evaluate the polynomials
  • N = (int) maximum order in the vanermonde matrix
  • k = (int) derivative order [default=0]
  • norm = (optional,boolean) True -> orthonormal polynomials, False -> non orthonormal polynomials
Output:
  • V = (2d-array(MxN),float) Generalized Vandermonde matrix
GramSchmidt(p, N, w)[source]

GramSchmidt(): creates a Generalized Vandermonde Matrix of orthonormal polynomials with respect to the weights w

Syntax:
V = GramSchmidt(p, N, w)
Input:
  • p = (1d-array,float) points at which to evaluate the new polynomials
  • N = (int) the maximum order of the polynomials
  • w = (1d-array,float) weights to be used for the orthogonoalization
Output:
  • V = Generalized Vandermonde Matrix containing the new orthogonalized polynomials
Description:
Takes the points where the polynomials have to be evaluated and computes a Generalized Gram Schmidth procedure, where a weighted projection is used. If w==1 then the usual inner product is used for the orthogonal projection.
InverseDiscretePolynomialTransform(r, fhat, N)[source]

Computes the nodal values from the modal form fhat.

Parameters:
  • x (1d-array or float) – set of points on which to the polynomials are evaluated
  • fhat (1d-array or float) – list of Polynomial coefficients
  • N (int) – maximum order in the generalized vanermonde matrix
Returns:

function values

Return type:

1d-array or float

Note

If the Chebyshev polynomials are chosen and r contains Chebyshev-Gauss-Lobatto points, the Inverse Fast Chebyshev Transform is used. Otherwise uses the Generalized Vandermonde Matrix in order to transform from transform space to physical space.

See also

InverseFastChebyshevTransform

InverseFastChebyshevTransform(fhat)[source]

InverseFastChebyshevTransform(): Returns the coefficients of the Inverse Fast Chebyshev Transform.

Syntax:
f = InverseFastChebyshevTransform(fhat)
Input:
  • fhat = (1d-array,float) list of Polynomial coefficients
Output:
  • f = (1d-array,float) function values

Note

If f is odd, the vector is padded with a zero value (highest freq.)

Note

Modification of algorithm (29) from [Kop09]

LegendreDerivativeCoefficients(fhat)[source]

LegendreDerivativeCoefficients(): computes the Legendre coefficients of the derivative of a function

Syntax:
dfhat = LegendreDerivativeCoefficients(fhat)
Input:
  • fhat = (1d-array,float) list of Legendre coefficients of the original function
Output:
  • dfhat = (1d-array,float) list of Legendre coefficients of the derivative of the original function

Notes

Algorithm (4) from [Kop09]

PolyInterp(x, f, xi, order)[source]

PolyInterp(): Interpolate function values f from points x to points xi using Forward and Backward Polynomial Transform

Syntax:
fi = PolyInterp(x, f, xi)
Input:
  • x = (1d-array,float) set of N original points where f is evaluated
  • f = (1d-array,float) set of N function values
  • xi = (1d-array,float) set of M points where the function is interpolated
  • order = (integer) order of polynomial interpolation
Output:
  • fi = (1d-array,float) set of M function values

Notes:

Quadrature(N, quadType=None, normed=False, left=None, right=None, end=None)[source]

Quadrature(): Generates list of nodes and weights for the quadType quadrature rule using the selected Polynomial basis

Syntax:
(x,w) = Quadrature(N, [quadType=None], [normed=False], [left=None], [right=None], [end=None])
Input:
  • N = (int) accuracy level required
  • quadType = (AVAIL_QUADPOINTS) type of quadrature to be used. Default is Gauss quadrature rule.
  • normed = (optional,bool) whether the weights will be normalized or not
  • left = (optional,float) containing the left endpoint (used by ORTHPOL Gauss-Lobatto rules)
  • right = (optional,float) containing the right endpoint (used by ORTHPOL Gauss-Lobatto rules)
  • end = (optional,float) containing the endpoint (used by ORTHPOL Gauss-Radau rules)
Output:
  • x = (1d-array,float) containing the nodes
  • w = (1d-array,float) containing the weights
interpolate(x, f, xi, order)[source]

Interpolates function values f from points x to points xi using Forward and Backward Polynomial Transform

Parameters:
  • x (1d-array,float) – set of N original points where f is evaluated
  • f (1d-array,float) – set of N function values
  • xi (1d-array,float) – set of M points where the function is interpolated
  • order (int) – order of polynomial interpolation
Returns:

fi (1d-array,float) set of M function values

..note:: this is the same of calling Ploy1D.PolyInterp() .. seealso:: PolyInterp

orthpol_rec_coeffs(N)[source]

Compute the recursion coefficients for polynomials up to the N-th order

project(r, f, N)[source]

Computes the Discrete Polynomial Transform of function values f, i.e. project on the space spanned by the selected polynomials up to order N.

Parameters:
  • r (1d-array or float) – set of points on which to the polynomials are evaluated
  • f (1d-array or float) – function values
  • N (int) – maximum order in the generalized vanermonde matrix
Returns:

projection coefficients

Return type:

1d-ndarray containing the projection coefficients. If len(r) == N+1 exact projection is used. If len(r) != N+1 the least square solution is obtained by numpy.linalg.lstsq.

Note

If the Chebyshev polynomials are chosen and r contains Chebyshev-Gauss-Lobatto points, the Fast Chebyshev Transform is used. Otherwise uses the Generalized Vandermonde Matrix in order to transform from physical space to transform space.

Note

this is the same of calling DiscretePolynomialTransform()

See also

DiscretePolynomialTransform

See also

numpy.linalg.lstsq

SpectralToolbox.Spectral1D.nestedlobatto(N, norm=True)[source]

nestedlobatto(): function for generating 1D Nested rule for integral with Uniform weight with 2**l scaling

Syntax:
(n,w) = nestedlobatto(N)
Input:
N = level of accuracy of the quadrature rule
Output:
n = nodes w = weights
SpectralToolbox.Spectral1D.nestedgauss(N, norm=True)[source]

nestedgauss(): function for generating 1D Nested rule for integral with Uniform weight with 2**l scaling

Syntax:
(n,w) = nestedgauss(N)
Input:
N = level of accuracy of the quadrature rule
Output:
n = nodes w = weights
SpectralToolbox.Spectral1D.fej(N, norm=True)[source]

fej(): function for generating 1D Nested Fejer’s rule [-1,1]

Syntax:
(n,w) = fej(N)
Input:
N = level of accuracy of the quadrature rule
Output:
n = nodes w = weights
SpectralToolbox.Spectral1D.cc(N, norm=True)[source]

cc(): function for generating 1D Nested Clenshaw-Curtis [-1,1]

Syntax:
(n,w) = cc(N)
Input:
N = level of accuracy of the quadrature rule
Output:
n = nodes w = weights
SpectralToolbox.Spectral1D.kpn(N)[source]

KPN(): function for generating 1D Nested rule for integral with Gaussian weight

Note:: Max ``N’’ is 25

Syntax:
(n,w) = GQU(l)
Input:
l = level of accuracy of the quadrature rule
Output:
n = nodes w = weights
SpectralToolbox.Spectral1D.kpu(N, norm=True)[source]

KPU(): function for generating 1D Nested rule for unweighted integral over [-1,1]

Note:: Max ``N’’ is 25

Syntax:
(n,w) = GQU(l)
Input:
l = level of accuracy of the quadrature rule
Output:
n = nodes w = weights
SpectralToolbox.Spectral1D.gqn(N)[source]

GQN(): function for generating 1D Gaussian quadrature for integral with Gaussian weight (Gauss-Hermite)

Note:: Max ``N’’ is 25

Syntax:
(n,w) = GQU(l)
Input:
l = level of accuracy of the quadrature rule
Output:
n = nodes w = weights
SpectralToolbox.Spectral1D.gqu(N, norm=True)[source]

GQU(): function for generating 1D Gaussian quadrature rules for unweighted integral over [-1,1] (Gauss-Legendre)

Note:: Max ``N’’ is 25

Syntax:
(n,w) = GQU(l)
Input:
l = level of accuracy of the quadrature rule
Output:
n = nodes w = weights
SpectralToolbox.Spectral1D.generate(ptype, params)[source]

Generate orthogonal basis objects from Spectral1D.AVAIL_POLY.

Parameters:
  • ptype (string) – one of the available polynomial types as listed in Spectral1D.AVAIL_POLY
  • params (list) – list of parameters need.
Returns:

(OrthogonalBasis) – the orthogonal basis required