Sparse Grids

SpectralToolbox.SparseGrids.CC(l)[source]
SpectralToolbox.SparseGrids.FEJ(l)[source]
SpectralToolbox.SparseGrids.GQN(): function for generating 1D Gaussian quadrature for integral with Gaussian weight (Gauss-Hermite)[source]
Parameters:l (int) – level of accuracy of the quadrature rule
Returns:(tuple [2] of ndarray) – nodes and weights
SpectralToolbox.SparseGrids.GQU(): function for generating 1D Gaussian quadrature rules for unweighted integral over [0, 1] (Gauss-Legendre)[source]
Parameters:l (int) – level of accuracy of the quadrature rule
Returns:(tuple [2] of ndarray) – nodes and weights
SpectralToolbox.SparseGrids.KPN(l)[source]

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

Parameters:l (int) – level of accuracy of the quadrature rule
Returns:(tuple [2] of ndarray) – nodes and weights
SpectralToolbox.SparseGrids.KPU(l)[source]

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

Parameters:l (int) – level of accuracy of the quadrature rule
Returns:(tuple [2] of ndarray) – nodes and weights
SpectralToolbox.SparseGrids.NESTEDGAUSS(l)[source]

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

Parameters:l (int) – level of accuracy of the quadrature rule
Returns:(tuple [2] of ndarray) – nodes and weights
SpectralToolbox.SparseGrids.NESTEDLOBATTO(l)[source]

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

Parameters:l (int) – level of accuracy of the quadrature rule
Returns:(tuple [2] of ndarray) – nodes and weights
class SpectralToolbox.SparseGrids.SparseGrid(qrule, dim, k, sym)[source]

Bases: object

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)
  • CC = Clenshaw-Curtis quadrature for unweighted integral over [-1,1]
  • FEJ = Fejer’s quadrature for unweighted integral over [-1,1]
  • 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)
DIM = 0
K = 0
QRULE = ''
SYM = 1
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