Welcome to DABISpectralToolbox’s documentation!

Contents:

Created on Wed Jul 11 10:40:30 2012

@author: Daniele Bigoni (dabi@imm.dtu.dk)

@version: 0.1

DABISpectralToolbox is the collection of tools for spectral methods and numerical integration.

@copyright: 2012 Daniele Bigoni

DABI Spectral 1D

Created on Thu Mar 29 11:33:32 2012

@author: Daniele Bigoni (dabi@imm.dtu.dk)

Description

Implementation of Spectral Methods in 1 dimension.

Available polynomials:
  • Jacobi Polynomials
  • Hermite Physicist
  • Hermite Function
  • Hermite Probabilistic
  • Laguerre Polynomial
  • Laguerre Function
  • ORTHPOL package (generation of recursion coefficients using [4])

Jacobi Polynomials

Jacobi polynomials are defined on the domain Ω=[-1,1] by the recurrence relation

xPn (α,β)(x)= 2(n+1)(n+α+β+1) (2n+α+β+1)(2n+α+β+2)P n+1 (α,β)(x) + β2-α2 (2n+α+β)(2n+α+β+2)P n (α,β)(x) + 2(n+α)(n+β) (2n+α+β)(2n+α+β+1)P n-1 (α,β)(x)

with weight function

w(x;α,β)= Γ(α+β+2) 2 α+β+1Γ(α+1)Γ(β+1)(1-x)α(1+x)β

Note

In probability theory, the Beta distribution is defined on Ψ=[0,1] and its the Probability Distribution Function is

ρB(x;α,β)= Γ(α+β) Γ(α)Γ(β)x α-1(1-x)**(β-1)

The relation betwen w(x;α,β) and ρB(x;α,β) for xΨ is

ρB(x;α,β)=2*w(2*x-1;β-1,α-1)

For example:

>>> from scipy import stats
>>> plot(xp,stats.beta(3,5).pdf(xp))
>>> plot(xp,2.*Bx(xx,4,2),'--')
>>> plot(xp,stats.beta(3,8).pdf(xp))
>>> plot(xp,2.*Bx(xx,7,2),'--')

References

[1](1, 2, 3, 4, 5, 6, 7, 8, 9) “Implemenenting Spectral Methods for Partial Differential Equations” by David A. Kopriva, Springer, 2009
[2]
  1. Shen and L.L. Wang, “Some recent advances on spectral methods for unbounded domains”. Communications in Computational Physics, vol. 5, no. 2–4, pp. 195–241, 2009
[3]W.H. Press, Numerical Recipes: “The Art of Scientific Computing”. Cambridge University Press, 2007
[4](1, 2, 3, 4)
  1. Gautschi, “Algorithm 726: ORTHPOL – a package of routines for generating orthogonal polynomials and Gauss-type quadrature rules”. ACM Trans. Math. Softw., vol. 20, issue 1, pp. 21-62, 1994
[5]Y.Shi-jun, “Gauss-Radau and Gauss-Lobatto formulae for the Jacobi weight and Gori-Micchelli weight functions”. Journal of Zhejiang University Science, vol. 3, issue 4, pp. 455-460

Examples

class DABISpectralToolbox.DABISpectral1D.Poly1D(poly, params, sdout=<open file '<stderr>', mode 'w' at 0x2ab42e489270>)[source]

Initialization of the Polynomial instance.

Syntax:
p = Poly1D(poly,params)
Input:
  • poly = The orthogonal polynomial type desired
  • params = The parameters needed by the selected polynomial
  • sdout = (optional,default=sys.stderr) output stream for logging
Description:

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’ or DABISpectral1D.JACOBI
  • ‘HermiteP’ or DABISpectral1D.HERMITEP
  • ‘HermiteF’ or DABISpectral1D.HERMITEF
  • ‘HermitePprob’ or DABISpectral1D.HERMITEP_PROB
  • ‘LaguerreP’ or DABISpectral1D.LAGUERREP
  • ‘LaguerreF’ or DABISpectral1D.LAGUERREF
  • ‘ORTHPOL’ or DABISpectral1D.ORTHPOL

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

Note

The ORTHPOL polynomials are built up using the “Multiple-Component Discretization Procedure” described in [4]. 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:

  • ncapm: (int) maximum integer N0 (default = 500)
  • 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, yp: (Numpy 1d-array) of dimension mp, containing the abscissas and the jumps of the point spectrum
  • mu: (function) measure function that returns the mass (float) with arguments: x (float) absissa, i (int) interval number in the continuous part
  • 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
  • finl, finr: (bool) specify whether the extreme left/right interval is finite (false for infinite)
  • endl, 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 [4] are suppressed. Instead the routine qgp of ORTHPOL [4] is used by default (iq=0 and idelta=2)

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
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 [1]
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 [1]
DiscretePolynomialTransform(r, f, N)[source]

DiscretePolynomialTransform(): computes the Discrete Polynomial Transform of function values f

Syntax:
fhat = DiscretePolynomialTransform(r, f, N)
Input:
  • r = (1d-array,float) set of points on which to the polynomials are evaluated
  • f = (1d-array,float) function values
  • N = (int) maximum order in the generalized vanermonde matrix
Output: one of the two following output is given, depending on the length of r
  • (fhat, residues, rank,s) = list of Polynomial coefficients and additional ouputs from function numpy.linalg.lstsq (if len(r) > N+1)
Description:
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.

See also

FastChebyshevTransform

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 [1]

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 [1]
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, norm=True)[source]

GradEvaluate(): evaluate the k-th derivative of the N-th order polynomial at points r

Syntax:
P = GradEvaluate(r,N,k[,norm=True])
Input:
  • r = (1d-array,float) set of points on which to evaluate the polynomial
  • N = (int) order of the polynomial
  • k = (int) order of the derivative
  • norm = (bool) whether to return normalized (True) or non normalized (False) polynomials
Output:
  • P = Polynomial evaluated on r
GradVandermonde1D(r, N, k, 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
  • 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]

InverseDiscretePolynomialTransform(): computes the nodal values from the modal form fhat.

Syntax:
f = InverseDiscretePolynomialTransform(r, fhat, alpha, beta, N)
Input:
  • x = (1d-array,float) set of points on which to the polynomials are evaluated
  • fhat = (1d-array,float) list of Polynomial coefficients
  • N = (int) maximum order in the generalized vanermonde matrix
Output:
  • f = (1d-array,float) function values
Description:
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 [1]

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 [1]
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 [1]
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 [1]
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:

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 [1]
Quadrature(N, quadType=None, normed=False, left=None, right=None, end=None)

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_POLYS) 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

DABI Spectral ND

Created on Mon Jul 9 11:35:12 2012

@author: Daniele Bigoni (dabi@imm.dtu.dk)

Implementation of Spectral Methods in n dimension.

It uses the package DABISpectral1D as basic polynomials in order to construct higher dimensional rules by tensor product.

class DABISpectralToolbox.DABISpectralND.PolyND(polys)[source]

Initialization of the N-dimensional Polynomial instance

Syntax:
p = PolyND(polys)
Input:
  • polys = (list,DABISpectral1D.Poly1D) list of polynomial instances of the class DABISpectral1D.Poly1D

See also

DABISpectral1D.Poly1D

GaussLobattoQuadrature(Ns, normed=False, warnings=True)[source]

GaussLobattoQuadrature(): computes the tensor product of the Guass Lobatto Points and weights

Syntax:
(x,w) = GaussLobattoQuadrature(Ns,[normed=False],[warnings=True])
Input:
  • Ns = (list,int) n-dimensional list with the order of approximation of each polynomial
  • normed = (optional,boolean) whether the weights will be normalized or not
  • warnings = (optional,boolean) set whether to ask for confirmation when it is required to allocate more then 100Mb of memory
Output:
  • x = tensor product of the collocation points
  • w = tensor product of the weights

Warning

The lengths of Ns has to be conform to the number of polynomials with which you have instantiated PolyND

GaussQuadrature(Ns, normed=False, warnings=True)[source]

GaussQuadrature(): computes the tensor product of the Guass Points and weights

Syntax:
(x,w) = GaussQuadrature(Ns, [normed=False],[warnings=True])
Input:
  • Ns = (list,int) n-dimensional list with the order of approximation of each polynomial
  • normed = (optional,boolean) whether the weights will be normalized or not
  • warnings = (optional,boolean) set whether to ask for confirmation when it is required to allocate more then 100Mb of memory
Output:
  • x = tensor product of the collocation points
  • w = tensor product of the weights

Warning

The lengths of Ns has to be conform to the number of polynomials with which you have instantiated PolyND

GradVandermonde(rs, Ns, ks, norms=None, usekron=True, output=True, warnings=True)[source]

GradVandermonde(): initialize the tensor product of the k-th gradient of the modal basis.

Syntax:
V = GradVandermonde(r,N,k,[norms=None],[usekron=True],[output=True],[warnings=True])
Input:
  • rs = (list of 1d-array,float) n-dimensional list of set of points on which to evaluate the polynomials (by default they are not the kron product of the points. See usekron option)
  • Ns = (list,int) n-dimensional list with the maximum orders of approximation of each polynomial
  • ks = (list,int) n-dimensional list with derivative orders
  • norms = (default=None,list,boolean) n-dimensional list of boolean, True -> orthonormal, False -> orthogonal, None -> all orthonormal
  • usekron = (optional,boolean) set whether to apply the kron product of the single dimensional Vandermonde matrices or to multiply column-wise. kron(rs) and usekron==False is equal to rs and usekron==True
  • output = (optional,boolean) set whether to print out information about memory allocation
  • warnings = (optional,boolean) set whether to ask for confirmation when it is required to allocate more then 100Mb of memory
OUTPUT:
  • V = Tensor product of the Generalized Vandermonde matrices

Warning

The lengths of Ns , rs , ks , norms has to be conform to the number of polynomials with which you have instantiated PolyND

GradVandermondePascalSimplex(rs, N, ks, norms=None, output=True, warnings=True)[source]

GradVandermondePascalSimplex(): initialize k-th gradient of the modal basis up to the total order N

Syntax:
V = GradVandermonde(r,N,k,[norms=None],[output=True],[warnings=True])
Input:
  • rs = (list of 1d-array,float) n-dimensional list of set of points on which to evaluate the polynomials
  • N = (int) the maximum orders of the polynomial basis
  • ks = (list,int) n-dimensional list with derivative orders
  • norms = (default=None,list,boolean) n-dimensional list of boolean, True -> orthonormal, False -> orthogonal, None -> all orthonormal
  • output = (optional,boolean) set whether to print out information about memory allocation
  • warnings = (optional,boolean) set whether to ask for confirmation when it is required to allocate more then 100Mb of memory
OUTPUT:
  • V = Generalized Vandermonde matrix up to the N-th order

Warning

The lengths of rs , ks , norms has to be conform to the number of polynomials with which you have instantiated PolyND

Quadrature(Ns, quadTypes=None, left=None, right=None, end=None, normed=False, warnings=True)

GaussQuadrature(): computes the tensor product of the Guass Points and weights

Syntax:
(x,w) = GaussQuadrature(Ns, [quadTypes=None], [normed=False],[warnings=True])
Input:
  • Ns = (list,int) n-dimensional list with the order of approximation of each polynomial
  • quadTypes = (list,``DABISpectral1D.AVAIL_QUADPOINTS``) n-dimensional list of quadrature point types chosen among Gauss, Gauss-Radau, Gauss-Lobatto (using the definition in DABISpectral1D). If None, Gauss points will be generated by default
  • left: (list,float) list of left values used by ORTHPOL for Gauss-Lobatto rules (the dimensions where the value is not used can be set to anything)
  • right: (list,float) list of left values used by ORTHPOL for Gauss-Lobatto rules (the dimensions where the value is not used can be set to anything)
  • end: (list,float) list of left values used by ORTHPOL for Gauss-Radau rules (the dimensions where the value is not used can be set to anything)
  • normed = (optional,boolean) whether the weights will be normalized or not
  • warnings = (optional,boolean) set whether to ask for confirmation when it is required to allocate more then 100Mb of memory
Output:
  • x = tensor product of the collocation points
  • w = tensor product of the weights

Warning

The lengths of Ns has to be conform to the number of polynomials with which you have instantiated PolyND

DABI Sparse Grids

Created on Mon Jul 9 14:08:30 2012

@author: Daniele Bigoni (dabi@imm.dtu.dk)

Implementation of Sparse Grid (Smolyak) methods for numerical integration.

This implementation is a porting of the Sparse-Grid package developped in MatLab(R) at http://www.sparse-grids.de

@copyright: 2007 Florian Heiss, Viktor Winschel

@copyright: 2012 Daniele Bigoni

class DABISpectralToolbox.DABISparseGrids.SparseGrid(qrule, dim, k, sym)[source]

Initialization of the Sparse Grid instance.

Syntax:
sg = SparseGrid(qrule, dim, k, sym)
Input:
  • qrule = Function of 1D integration rule
  • dim = dimension of the integration problem
  • k = Accuracy level. The rule will be exact for polynomial up to total order 2k-1
  • sym = (optional) only used for own 1D quadrature rule (type not “KPU”,...). If sym is supplied and not=0, the code will run faster but will produce incorrect results if 1D quadrature rule is asymmetric.
Description:

Several 1D integration rules are available to be chosen for the qrule input parameter

  • KPU = Nested rule for unweighted integral over [0,1]
  • KPN = Nested rule for integral with Gaussian weight
  • GQU = Gaussian quadrature for unweighted integral over [0,1] (Gauss-Legendre)
  • GQN = Gaussian quadrature for integral with Gaussian weight (Gauss-Hermite)
  • func = any function provided by the user that accept level l and returns nodes n and weights w for univariate quadrature rule with polynomial exactness 2l-1 as [n w] = feval(func,level)
sparseGrid()[source]

sparseGrid(): main function for generating nodes & weights for sparse grids intergration

Syntax:
(n,w) = sparseGrid()
Output:
  • n = matrix of nodes with dim columns
  • w = row vector of corresponding weights
DABISpectralToolbox.DABISparseGrids.GQU(): function for generating 1D Gaussian quadrature rules for unweighted integral over [0, 1] (Gauss-Legendre)[source]
Syntax:
(n,w) = GQU(l)
Input:
l = level of accuracy of the quadrature rule
Output:
n = nodes w = weights
DABISpectralToolbox.DABISparseGrids.GQN(): function for generating 1D Gaussian quadrature for integral with Gaussian weight (Gauss-Hermite)[source]
Syntax:
(n,w) = GQU(l)
Input:
l = level of accuracy of the quadrature rule
Output:
n = nodes w = weights
DABISpectralToolbox.DABISparseGrids.KPU(l)[source]

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

(n,w) = GQU(l)
Input:
l = level of accuracy of the quadrature rule
Output:
n = nodes w = weights
DABISpectralToolbox.DABISparseGrids.KPN(l)[source]

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

(n,w) = GQU(l)
Input:
l = level of accuracy of the quadrature rule
Output:
n = nodes w = weights

DABI Miscellaneous

Created on Wed Feb 27

@author: Daniele Bigoni (dabi@dtu.dk)

Collection of miscellaneous functions used in the DABISpectralToolbox

class DABISpectralToolbox.DABIMisc.ExpandingArray(initData, allocInitDim=None, dtype=<type 'numpy.float64'>, maxIncrement=None)[source]

ExpandingArray is used for the dynamic allocation of memory in applications where the total allocated memory needed cannot be predicted. Memory is preallocated with increases of 50% all the time data exceed the allocated memory.

Initialization of the Expanding Array.

>>> EA = ExpandingArray(initData,[allocInitDim=None,[dtype=np.float64,[maxIncrement=None]]])
Parameters:
  • initData (ndarray) – InitialData with which to be initially filled. This must provide the number of dimensions of the array
  • allocInitDim (1darray-integer) – Initial allocated dimension (optional)
  • dtype (dtype) – type for the data that will be contained in EA (optional,default=np.float64)
  • maxIncrement (integer) – upper limit for the allocation increment
concatenate(X, axis=0)[source]

Concatenate data to the existing array. If needed the array is resized in the axis direction by a factor of 50%.

Parameters:
  • X (ndarray) – data to be concatenated to the array. Note that X.shape[i]==EA.shape()[i] is required for all i!=axis
  • axis (integer) – axis along which to concatenate the additional data (optional)
>>> import numpy as np
>>> EA = DABIMisc.ExpandingArray(np.random.rand(25,4))
>>> EA.shape()
(25, 4)
>>> EA.getAllocArray().shape
(25, 4)
>>> EA.concatenate(np.random.rand(13,4))
>>> EA.shape()
(38, 4)
getAllocArray()[source]

Return the allocated array.

Returns:allocated array
Return type:ndarray
>>> import numpy as np
>>> EA = DABIMisc.ExpandingArray(np.random.rand(25,4))
>>> EA.shape()
(25, 4)
>>> EA.getAllocArray().shape
(25, 4)
>>> EA.concatenate(np.random.rand(13,4))
>>> EA.shape()
(38, 4)
>>> EA.getAllocArray().shape
(45, 4)
getDataArray()[source]

Get the view of the array with data.

Returns:allocated array
Return type:ndarray
>>> EA.shape()
(38, 4)
>>> EA.getAllocArray().shape
(45, 4)
>>> EA.getDataArray().shape
(38, 4)
shape()[source]

Returns the shape of the data inside the array. Note that the allocated memory is always bigger or equal to shape().

Returns:shape of the data
Return type:tuple of integer
>>> import numpy as np
>>> EA = DABIMisc.ExpandingArray(np.random.rand(25,4))
>>> EA.shape()
(25, 4)
>>> EA.getAllocArray().shape
(25, 4)
>>> EA.concatenate(np.random.rand(13,4))
>>> EA.shape()
(38, 4)
>>> EA.getAllocArray().shape
(45, 4)
trim(N, axis=0)[source]

Trim the axis dimension of N elements. The allocated data is not reinitialized or deallocated. Only the dimensions of the view are redefined.

Parameters:
  • N (integer) – number of elements to be removed along the axis dimension
  • axis (integer) – axis along which to remove elements (optional)
>>> EA = DABIMisc.ExpandingArray(np.random.rand(4,2))
>>> EA.getDataArray()
array([[ 0.42129746,  0.76220921],
       [ 0.9238783 ,  0.11256142],
       [ 0.42031437,  0.87349243],
       [ 0.83187297,  0.555708  ]])
>>> EA.trim(2,axis=0)
>>> EA.getDataArray()
array([[ 0.42129746,  0.76220921],
       [ 0.9238783 ,  0.11256142]])
>>> EA.getAllocArray()
array([[ 0.42129746,  0.76220921],
       [ 0.9238783 ,  0.11256142],
       [ 0.42031437,  0.87349243],
       [ 0.83187297,  0.555708  ]])
DABISpectralToolbox.DABIMisc.MultiIndex(d, N)[source]

Generates the multi index ordering for the construction of multidimensional Generalized Vandermonde matrices

Parameters:
  • d (integer) – dimension of the simplex
  • N (integer) – maximum value of the sum of the indices
Returns:

array containing the ordered multi-indices

Return type:

2d-array of integer

>>> DABIMisc.MultiIndex(2,3)
array([[0, 0],
       [1, 0],
       [0, 1],
       [2, 0],
       [1, 1],
       [0, 2],
       [3, 0],
       [2, 1],
       [1, 2],
       [0, 3]])
DABISpectralToolbox.DABIMisc.almostEqual(x, y, tol)[source]

Check equality of two iterable objects up to a certain tolerance

Parameters:
  • x,y (iterable objects of floats) – values to be compared
  • tol (float) – tolerance to be used
Returns:

true if equal, false otherwise

Return type:

bool

>>> eps2 = 2.*DABIMisc.machineEpsilon(np.float64)
>>> DABIMisc.almostEqual(np.array([2.]),np.array([2.+0.5*eps2]),eps2)
True
>>> DABIMisc.almostEqual(np.array([2.]),np.array([2.+2.*eps2]),eps2)
False
DABISpectralToolbox.DABIMisc.almostEqualList(xArray, y, tol)[source]

Check equality of a list of floats against an iterable value up to certain tolerance

Parameters:
  • xArray (2d-array of floats) – values to be compared to y
  • y (iterable objects of floats) – values to be compared to
  • tol (float) – tolerance to be used
Returns:

array of booleans containing true where equal, false elsewhere.

Return type:

1d-array of bool

Syntax:
b = almostEqualList(xArray,y,tol)
>>> eps2 = 2.*DABIMisc.machineEpsilon(np.float64)
>>> X = np.random.rand(4,2)
>>> DABIMisc.almostEqualList(X,X[1,:],eps2)
array([False,  True, False, False], dtype=bool)
DABISpectralToolbox.DABIMisc.argsort_insertion(X, tol, start_idx=1, end_idx=None)[source]

Implements the insertion sort with binary_search. Returns permutation indices.

Parameters:
  • X (2d-array of floats) – values ordered by row according to the compare function
  • tol (float) – tolerance to be used
  • start_idx,end_idx (int) – starting and ending indices for the ordering (optional)
Returns:

permutation indices

Return type:

1d-array of integers

>>> X = np.random.rand(5,2)
>>> X
array([[ 0.56865133,  0.18490129],
       [ 0.01411459,  0.46076606],
       [ 0.64384365,  0.24998971],
       [ 0.47840414,  0.32554137],
       [ 0.12961966,  0.43712056]])
>>> perm = DABIMisc.argsort_insertion(X,eps2)
>>> X[perm,:]
array([[ 0.01411459,  0.46076606],
       [ 0.12961966,  0.43712056],
       [ 0.47840414,  0.32554137],
       [ 0.56865133,  0.18490129],
       [ 0.64384365,  0.24998971]])

Search for the minimum X bigger than val

Parameters:
  • X (2d-array of floats) – values ordered by row according to the compare function
  • val (1d-array of floats) – value to be compared to
  • lo,hi (integer) – staring and ending indices
  • tol (float) – tolerance to be used
  • perm (1d-array of integers) – possible permutation to be used prior to the search (optional)
Returns:

index pointing to the maximum X smaller than val. If perm is provided, perm[idx] points to the maximum X smaller than val

Return type:

integer

>>> X = np.arange(1,5).reshape((4,1))
>>> X
array([[1],
       [2],
       [3],
       [4]])
>>> DABIMisc.binary_search(X,np.array([2.5]),0,4,eps2)
>>> idx = DABIMisc.binary_search(X,np.array([2.5]),0,4,eps2)
>>> idx
2
>>> X[idx,:]
array([3])
DABISpectralToolbox.DABIMisc.compare(x, y, tol)[source]

Compares two iterable objects up to a certain tolerance

Parameters:
  • x,y (iterable objects of floats) – values to be compared
  • tol (float) – tolerance to be used
Returns:

-1 if (x-y) < tol, 1 if (x-y) > tol, 0 otherwise

Return type:

integer

>>> eps2 = 2.*DABIMisc.machineEpsilon(np.float64)
>>> DABIMisc.compare(np.array([2.]),np.array([2.+0.5*eps2]),eps2)
0
>>> DABIMisc.compare(np.array([2.]),np.array([2.+2.*eps2]),eps2)
-1
DABISpectralToolbox.DABIMisc.findOverlapping(XF, X, tol)

Finds overlapping points of XF on X grids of points. The two grids are ordered with respect to DABIMisc.compare().

Parameters:
  • XF,X (2d-array of floats) – values ordered by row according to the DABIMisc.compare().
  • tol (float) – tolerance to be used
Returns:

true values for overlapping points of XF on X, false for not overlapping points. Note: the overlapping return argument is a true-false indexing for XF.

Return type:

1d-array of bool

Example

>>> XF
array([[ -1.73205081e+00,   0.00000000e+00],
       [ -1.00000000e+00,  -1.00000000e+00],
       [ -1.00000000e+00,   0.00000000e+00],
       [ -1.00000000e+00,   1.00000000e+00],
       [  0.00000000e+00,  -1.73205081e+00],
       [  0.00000000e+00,  -1.00000000e+00],
       [  0.00000000e+00,   2.16406754e-16],
       [  0.00000000e+00,   1.00000000e+00],
       [  0.00000000e+00,   1.73205081e+00],
       [  1.00000000e+00,  -1.00000000e+00],
       [  1.00000000e+00,   0.00000000e+00],
       [  1.00000000e+00,   1.00000000e+00],
       [  1.73205081e+00,   0.00000000e+00]])
>>> X
array([[ -1.73205081e+00,   0.00000000e+00],
       [ -1.00000000e+00,  -1.00000000e+00],
       [ -1.00000000e+00,   0.00000000e+00],
       [ -1.00000000e+00,   1.00000000e+00],
       [  0.00000000e+00,  -1.00000000e+00],
       [  2.16406754e-16,   0.00000000e+00],
       [  0.00000000e+00,   1.00000000e+00],
       [  1.00000000e+00,  -1.00000000e+00],
       [  1.00000000e+00,   0.00000000e+00],
       [  1.00000000e+00,   1.00000000e+00],
       [  1.73205081e+00,   0.00000000e+00]])
>>> tol = 2. * DABIMisc.machineEpsilon()
>>> bool_idx_over = DABIMisc.findOverlapping(XF,X,tol)
>>> XF[np.logical_not(bool_idx_over),:]
array([[ 0.        , -1.73205081],
       [ 0.        ,  1.73205081]])
DABISpectralToolbox.DABIMisc.machineEpsilon(func=<type 'float'>)[source]

Returns the abolute machine precision for the type passed as argument

Parameters:func (dtype) – type
Returns:absolute machine precision
Return type:float
>>> DABIMisc.machineEpsilon(np.float64)
2.2204460492503131e-16
>>> DABIMisc.machineEpsilon(np.float128)
1.084202172485504434e-19
DABISpectralToolbox.DABIMisc.powerset(iterable)

Compute the power set of an iterable object.

DABISpectralToolbox.DABIMisc.unique_cuts(X, tol, retIdxs=False)

Returns the unique values and a list of arrays of boolean indicating the positions of the unique values. If retIdx is true, then it returns the group of indices with the same values as a indicator function (true-false array)

DABI Heterogeneous Sparse Grids

Created on Wed Feb 27

@author: Daniele Bigoni (dabi@dtu.dk)

Description

This module is used to construct Heterogeneous Sparse Grids. These are integration rules based on Smolyak’s Sparse Grids [6] built up in order to provide flexibility for both:
  1. the types of polynomials to be used per direction
  2. the accuracy to be used per direction

The types of polynomials available are all the ones included in the module DABISpectral1D. The rules don’t need to be symmetric and the accuracy per each direction can vary.

For rules with Heterogeneous accuracy, two sparse grids will be construced: one partial sparse grid and one full sparse grid (up to the maximum accuracy). The values computed for the partial sparse grid can then be used to interpolate on the points of the full sparse grid. This latter rule can then be used to compute the integral.

References

[6]
    1. Smolyak, “Quadrature and interpolation formulas for tensor products of certain classes of functions”. Soviet Math. Dokl., 4:240–243, 1963

Examples

Let’s consider the following space Ω=[-,]imes[0,1] with the associated measures μ1(x)= 1 2πe - x2 2 and μ2(x)=1 for each dimension. The product measure is given by μ( x¯)= i=1nμi(xi). We will consider the function

f( x¯)=x1 p1x2 p2
>>> f_mult = lambda x,y,xp,yp: x**xp * y**yp

with exact value of the integral given by

Ωf( x¯)μ( x¯)= 2 -1+p2 1+(-1) p2Γ 1+p2 2 (1+p1) π
>>> def I_mult(q,p):
>>>     return 2.**(-.5+.5*(-1.+q)) * (1.+(-1.)**q) * scipy.special.gamma((1.+q)/2.) / ((1.+p)*np.sqrt(np.pi))

Let p1=2 and p2=4. We can obtain a sparse grid composed using Hermite basis (DABISpectral1D.HERMITEP_PROB) and Legendre basis (DABISpectral1D.JACOBI) with orders 2 and 4 respectively.

>>> from DABISpectralToolbox import DABIHeterogeneousSparseGrids as HSG
>>> pH = DABISpectral1D.Poly1D(DABISpectral1D.HERMITEP_PROB,None)
>>> pL = DABISpectral1D.Poly1D(DABISpectral1D.JACOBI,[0.0,0.0])
>>> polys = [pH,pL]
>>> Ns = [2,4]
>>> sg = HSG.HSparseGrid(polys,Ns)
>>> (XF,W,X) = sg.sparseGrid()
[SG] Sparse Grid Generation [============================================] 100%
[SG] Sparse Grid Generation: 0.01s
>>> XF[:,1] = (XF[:,1]+1.)/2.
>>> X[:,1] = (X[:,1]+1.)/2.
>>> plt.figure()
>>> plt.plot(XF[:,0],XF[:,1],'o')
>>> plt.plot(X[:,0],X[:,1],'or')

The resulting partial and full sparse grids are shown in the following figure.

_images/DABIHeterSparseGrids-GridExample.png

Partial (red) and full (blue) sparse grid. The full sparse grid is overlapping over the partial sparse grid.

The values on the partial grid can be computed and then the interpolation is taken over the full sparse grid.

>>> fX = f_mult(X[:,0],X[:,1],Q,P)
>>> fXF = sg.sparseGridInterp(X,fX,XF)
[SG] Sparse Grid Interpolation [=========================================] 100%
[SG] Sparse Grid Interpolation: 0.00s

Finally the error of the quadrature rule is

>>> IErr = np.abs(0.5*np.dot(fXF,W)-I_mult(Q,P))
>>> print IErr
3.33066907388e-16
class DABISpectralToolbox.DABIHeterogeneousSparseGrids.HSparseGrid(polys, Ns, tol=None, sdout=<open file '<stderr>', mode 'w' at 0x2ab42e489270>)[source]

Heterogeneous Sparse Grid class

Constructor of Heterogeneous Sparse Grid object (this does not allocate the sparse grid)

Parameters:
  • polys (list of DABISpectral1D.Poly1D) – orthogonal polynomials to be used as basis functions
  • Ns (list of integers) – accuracy for each dimension. It can be a list of accuracies or a single accuracy, in which case uniform accuracy is assumed
  • tol (float) – tolerance to be used when comparing points of the grid (optional, default=:py:func:DABIMisc.machineEpsilon())
  • sdout (stream) – default output stream for the class (optional,default=``sys.stderr``)

Note

one of the following must hold: len(polys)==len(Ns) or len(Ns)==1, in which case the same order is used for all the directions.

Example

>>> from DABISpectralToolbox import DABIHeterogeneousSparseGrids as HSG
>>> pH = DABISpectral1D.Poly1D(DABISpectral1D.HERMITEP_PROB,None)
>>> pL = DABISpectral1D.Poly1D(DABISpectral1D.JACOBI,[0.0,0.0])
>>> polys = [pH,pL]
>>> Ns = [2,4]
>>> sg = HSG.HSparseGrid(polys,Ns)
sparseGrid(heter=False)

Generates the full and partial sparse grids

Parameters:heter (bool) – if Ns is homogeneous, this parameter will force the output of the partial sparse grid as well
Returns:tuple (XF,WF,X) containing:
  • XF: full grid points
  • WF: full grid weights
  • X: partial grid points

Example

>>> (XF,W,X) = sg.sparseGrid()
[SG] Sparse Grid Generation [============================================] 100%
[SG] Sparse Grid Generation: 0.01s
sparseGridInterp(X, fX, XF)

Interpolate values of the Sparse Grid using 1D Polynomial interpolation along cuts.

Parameters:
  • X (2d-array of floats) – partial grid on which a function has been evaluated
  • fX (1d-array of floats) – values for the points in the partial grid
  • XF (2d-array of floats) – full grid on which to interpolate
Returns:

fXF the interpolated values on the full grid

Return type:

1d-array of floats

..note:: The partial and full grid must be overlapping

Example

>>> fXF = sg.sparseGridInterp(X,fX,XF)
[SG] Sparse Grid Interpolation [=========================================] 100%
[SG] Sparse Grid Interpolation: 0.00s

Indices and tables