Quadrature¶

Base classes¶

class modepy.Quadrature(nodes: ndarray, weights: ndarray, exact_to: int | None = None)[source]¶

The basic interface for a quadrature rule.

nodes: ndarray¶

An array of shape (dim, nnodes), where dim is the dimension of the qudrature rule. In 1D, the shape is just (nnodes,).

weights: ndarray¶

A vector of length nnodes that contains the quadrature weights.

property dim: int¶

Dimension of the space on which the quadrature rule applies.

property exact_to: int¶

Summed polynomial degree up to which the quadrature is exact.

In higher-dimensions, the quadrature is supposed to be exact on (at least) \(P^N\), where \(N\) = exact_to. If the quadrature accuracy is not known, exact_to will not be set, and an AttributeError will be raised when attempting to access this information.

__call__(f: Callable[[ndarray], ndarray]) ndarray[source]¶

Evaluate the callable f at the quadrature nodes and return its integral.

f is assumed to accept arrays of shape (dims, npts), or of shape (npts,) for 1D quadrature.

class modepy.ZeroDimensionalQuadrature[source]¶

Bases: Quadrature

A quadrature rule that should be used for 0d domains (i.e. points).

modepy.quadrature_for_space(space: FunctionSpace, shape: Shape) Quadrature[source]¶
modepy.quadrature_for_space(space: PN, shape: Simplex)
modepy.quadrature_for_space(space: TensorProductSpace, shape: TensorProductShape)
Returns:

a Quadrature that exactly integrates the functions in space over shape.

Jacobi-Gauss quadrature in one dimension¶

class modepy.JacobiGaussQuadrature(alpha: float, beta: float, N: int, backend: str | None = None, force_dim_axis: bool = False)[source]¶

Bases: Quadrature

A Gauss quadrature of order N associated with the Jacobi polynomials.

The quadrature rule can be used for weighted integrals of the form

\[I[f] = \int_{-1}^1 f(x) (1 - x)^\alpha (1 + x)^\beta\, \mathrm{d}x,\]

where \(\alpha, \beta > -1\). The quadrature rule is exact up to degree \(2N + 1\).

__init__(alpha: float, beta: float, N: int, backend: str | None = None, force_dim_axis: bool = False) None[source]¶
Parameters:

backend – Either "builtin" or "scipy". When the "builtin" backend is in use, there is an additional requirement that \(\alpha + \beta \ne -1\), with the exception of the Chebyshev quadrature \(\alpha = \beta = -1/2\). The "scipy" backend has no such restriction.

alpha: float¶

Power of \((1 - x)\) term in Jacobi quadrature.

beta: float¶

Power of \((1 + x)\) term in Jacobi quadrature.

static compute_weights_and_nodes(N: int, alpha: float, beta: float) Tuple[ndarray, ndarray][source]¶
Parameters:
  • N – order of the Gauss quadrature (the order of exactly integrated polynomials is \(2 N + 1\)).

  • alpha – power of \(1 - x\) in the Jacobi polynomial weight.

  • beta – power of \(1 + x\) in the Jacobi polynomial weight.

Returns:

a tuple (nodes, weights) of quadrature notes and weights.

class modepy.LegendreGaussQuadrature(N: int, backend: str | None = None, force_dim_axis: bool = False)[source]¶

Bases: JacobiGaussQuadrature

A Gauss quadrature rule with weight \(1\).

Corresponds to a Gauss-Jacobi quadrature rule with \(\alpha = \beta = 0\).

class modepy.ChebyshevGaussQuadrature(N: int, kind: int = 1, backend: str | None = None, force_dim_axis: bool = False)[source]¶

Bases: JacobiGaussQuadrature

A Gauss quadrature rule with weight \(\sqrt{1-x^2}^{\mp 1}\).

The Chebyshev-Gauss quadrature rules of the first kind and second kind correspond to Gauss-Jacobi quadrature rules with \(\alpha = \beta = -0.5\) and \(\alpha = \beta = 0.5\), respectively.

Added in version 2019.1.

class modepy.GaussGegenbauerQuadrature(alpha: float, N: int, backend: str | None = None, force_dim_axis: bool = False)[source]¶

Bases: JacobiGaussQuadrature

Gauss-Gegenbauer quadrature is a special case of Gauss-Jacobi quadrature with \(\alpha = \beta\).

Added in version 2019.1.

modepy.quadrature.jacobi_gauss.jacobi_gauss_lobatto_nodes(alpha: float, beta: float, N: int, backend: str | None = None, force_dim_axis: bool = False) ndarray[source]¶

Compute the Gauss-Lobatto quadrature nodes corresponding to the JacobiGaussQuadrature with the same parameters.

Exact to degree \(2N - 3\).

modepy.quadrature.jacobi_gauss.legendre_gauss_lobatto_nodes(N: int, backend: str | None = None, force_dim_axis: bool = False) ndarray[source]¶

Compute the Legendre-Gauss-Lobatto quadrature nodes.

Exact to degree \(2N - 1\).

Clenshaw-Curtis and Fejér quadrature in one dimension¶

class modepy.ClenshawCurtisQuadrature(N: int, force_dim_axis: bool = False)[source]¶

Bases: Quadrature

Clenshaw-Curtis quadrature of order N with N + 1 points.

The quadrature rule is exact up to degree \(N\) and can be nested. Its performance for differentiable functions is comparable with the classic Gauss-Legendre quadrature, which is exact for polynomials of degree up to \(2N + 1\). Implementation is based on [Waldvogel2003].

Integrates on the interval \((-1, 1)\).

[Waldvogel2003] (1,2)

Jörg Waldvogel, Fast Construction of the Fejer and Clenshaw-Curtis Quadrature Rules, BIT Numerical Mathematics, 2003, Vol. 43, No. 1, pp. 001-018. DOI

class modepy.FejerQuadrature(N: int, kind: int = 1, force_dim_axis: bool = False)[source]¶

Bases: Quadrature

Fejér quadrature rules of order N.

  • Fejér quadrature of the first kind has N points and uses Chebyshev nodes, i.e. the roots of Chebyshev polynomials.

  • Fejér quadrature of the second kind has N - 1 points and uses only the interior extrema of the Chebyshev nodes, i.e. the true stationary points. This rule is almost identical to Clenshaw-Curtis and can be nested.

Integrates on the interval \((-1, 1)\). Implementation is based on [Waldvogel2003].

kind: int¶

Kind of the Fejér quadrature, either first-kind or second-kind.

Quadratures on the simplex¶

exception modepy.QuadratureRuleUnavailable[source]¶

Added in version 2013.3.

class modepy.GrundmannMoellerSimplexQuadrature(order: int, dimension: int)[source]¶

Bases: Quadrature

Cubature on an n-simplex from [Grundmann1978].

This cubature rule has both negative and positive weights. It is exact for polynomials up to order \(2s + 1\), where \(s\) is given as order. The integration domain is the unit simplex (see Coordinates on the triangle and Coordinates on the tetrahedron).

[Grundmann1978]

A. Grundmann and H.M. Moeller, Invariant integration formulas for the n-simplex by combinatorial methods, SIAM J. Numer. Anal. 15 (1978), 282–290. DOI

__init__(order: int, dimension: int) None[source]¶
Parameters:
  • order – a parameter correlated with the total degree of polynomials that are integrated exactly (see also exact_to).

  • dimension – the number of dimensions for the quadrature rule.

__call__(f: Callable[[ndarray], ndarray]) ndarray[source]¶

Evaluate the callable f at the quadrature nodes and return its integral.

f is assumed to accept arrays of shape (dims, npts), or of shape (npts,) for 1D quadrature.

class modepy.XiaoGimbutasSimplexQuadrature(order: int, dims: int)[source]¶

Bases: Quadrature

A (nearly) Gaussian simplicial quadrature with very few quadrature nodes from [Xiao2010].

This rule is available for low-to-moderate orders. The integration domain is the unit simplex (see Coordinates on the triangle and Coordinates on the tetrahedron).

[Xiao2010]

H. Xiao and Z. Gimbutas, A numerical algorithm for the construction of efficient quadrature rules in two and higher dimensions, Computers & Mathematics with Applications, vol. 59, no. 2, pp. 663-676, 2010. DOI

__init__(order: int, dims: int) None[source]¶
Parameters:
  • order – the total degree to which the quadrature rule is exact.

  • dims – the number of dimensions for the quadrature rule. 2 for quadrature on triangles and 3 for tetrahedra.

Raises:

modepy.QuadratureRuleUnavailable if no quadrature rule for therequested parameters is available.

__call__(f: Callable[[ndarray], ndarray]) ndarray[source]¶

Evaluate the callable f at the quadrature nodes and return its integral.

f is assumed to accept arrays of shape (dims, npts), or of shape (npts,) for 1D quadrature.

class modepy.VioreanuRokhlinSimplexQuadrature(order: int, dims: int)[source]¶

Bases: Quadrature

Simplicial quadratures with symmetric node sets and positive weights suitable for well-conditioned interpolation.

The integration domain is the unit simplex (see Coordinates on the triangle and Coordinates on the tetrahedron). When using these nodes, please acknowledge Zydrunas Gimbutas, who generated them as follows:

  • The 2D nodes are based on the interpolation node set derived in the article [Vioreanu2011].

    Note that in Vioreanu’s tables, only orders 5, 6, 9, and 12 are rotationally symmetric, which gives one extra order for integration and better interpolation conditioning. Also note that since the tables have been re-generated independently, the nodes and weights may be different.

  • The 3D nodes were derived from the modepy.warp_and_blend_nodes().

  • A tightening algorithm was then applied, as described in [Vioreanu2012].

[Vioreanu2011]

B. Vioreanu and V. Rokhlin, Spectra of Multiplication Operators as a Numerical Tool, Yale CS Tech Report 1443. PDF

[Vioreanu2012]

B. Vioreanu, Spectra of Multiplication Operators as a Numerical Tool, Yale University, 2012. PDF

Added in version 2013.3.

__init__(order: int, dims: int) None[source]¶
Parameters:
  • order – The total degree to which the quadrature rule is exact for interpolation.

  • dims – The number of dimensions for the quadrature rule. 2 for quadrature on triangles and 3 for tetrahedra.

Raises:

modepy.QuadratureRuleUnavailable if no quadrature rule for the requested parameters is available.

__call__(f: Callable[[ndarray], ndarray]) ndarray[source]¶

Evaluate the callable f at the quadrature nodes and return its integral.

f is assumed to accept arrays of shape (dims, npts), or of shape (npts,) for 1D quadrature.

Quadratures on the hypercube¶

class modepy.WitherdenVincentQuadrature(order: int, dims: int)[source]¶

Bases: Quadrature

Symmetric quadrature rules with positive weights for rectangles and hexahedra from [Witherden2015].

The integration domain is the unit hypercube \([-1, 1]^d\), where \(d\) is the dimension.

[Witherden2015]

F. D. Witherden, P. E. Vincent, On the Identification of Symmetric Quadrature Rules for Finite Element Methods, Computers & Mathematics with Applications, Vol. 69, pp. 1232–1241, 2015. DOI

__init__(order: int, dims: int) None[source]¶
Parameters:
  • nodes – an array of shape (d, nnodes), where d is the dimension of the qudrature rule.

  • weights – an array of length nnodes.

  • exact_to – an optional argument denoting the summed polynomial degree to which the quadrature is exact. By default, exact_to is None and will not be set as an attribute.

__call__(f: Callable[[ndarray], ndarray]) ndarray[source]¶

Evaluate the callable f at the quadrature nodes and return its integral.

f is assumed to accept arrays of shape (dims, npts), or of shape (npts,) for 1D quadrature.

class modepy.TensorProductQuadrature(quads: Iterable[Quadrature])[source]¶

Bases: Quadrature

A tensor product quadrature of one-dimensional Quadratures.

__init__(quads: Iterable[Quadrature]) None[source]¶
Parameters:

quad – a iterable of Quadrature objects for one-dimensional intervals, one for each dimension of the tensor product.

class modepy.LegendreGaussTensorProductQuadrature(N: int, dims: int, backend: str | None = None)[source]¶

Bases: TensorProductQuadrature

A tensor product using only LegendreGaussQuadrature one-dimenisonal rules.