Quadrature#

This module contains quadrature (numeric integration) schemes for different finite element formulations. The integration points of a boundary-quadrature are located on the first edge for 2d elements and on the first face for 3d elements.

Lines, Quads and Hexahedrons

GaussLegendre(order, dim[, permute])

An arbitrary-order Gauss-Legendre quadrature rule of dim 1, 2 or 3 on the interval [-1, 1].

GaussLegendreBoundary(order, dim[, permute])

An arbitrary-order Gauss-Legendre quadrature rule of dim 1, 2 or 3 on the interval [-1, 1].

Triangles and Tetrahedrons

TriangleQuadrature

alias of Triangle

TetrahedronQuadrature

alias of Tetrahedron

quadrature.Triangle(order)

An integration scheme for Triangles.

quadrature.Tetrahedron(order)

An integration scheme for Tetrahedrons.

Detailed API Reference

class felupe.quadrature.Scheme(points, weights)[source]#

Bases: object

A quadrature scheme with integration points \(x_q\) and weights \(w_q\). It approximates the integral of a function over a region \(V\) by a weighted sum of function values \(f_q = f(x_q)\), evaluated on the quadrature- points.

Notes

The approximation is given by

\[\int_V f(x) dV \approx \sum f(x_q) w_q\]

with quadrature points \(x_q\) and corresponding weights \(w_q\).

plot(plotter=None, add_axes=True, point_size=100, **kwargs)[source]#

Plot the quadrature points, scaled by their weights, into a (optionally provided) PyVista plotter.

See also

felupe.Scene.plot

Plot method of a scene.

screenshot(*args, filename=None, transparent_background=None, scale=None, **kwargs)[source]#

Take a screenshot of the quadrature.

See also

pyvista.Plotter.screenshot

Take a screenshot of a PyVista plotter.

class felupe.GaussLegendre(order: int, dim: int, permute: bool = True)[source]#

Bases: Scheme

An arbitrary-order Gauss-Legendre quadrature rule of dim 1, 2 or 3 on the interval [-1, 1].

Parameters:
  • order (int) – The number of sample points \(n\) minus one. The quadrature rule integrates degree \(2n-1\) polynomials exactly.

  • dim (int) – The dimension of the quadrature region.

  • permute (bool, optional) – Permute the quadrature points according to the cell point orderings (default is True).

Notes

The approximation is given by

\[\int_{-1}^1 f(x) dx \approx \sum f(x_q) w_q\]

with quadrature points \(x_q\) and corresponding weights \(w_q\).

Examples

>>> import felupe as fem
>>> fem.GaussLegendre(order=2, dim=3).screenshot()
../_images/quadrature.png
inv()[source]#

Return the inverse quadrature scheme.

class felupe.GaussLegendreBoundary(order: int, dim: int, permute: bool = True)[source]#

Bases: GaussLegendre

An arbitrary-order Gauss-Legendre quadrature rule of dim 1, 2 or 3 on the interval [-1, 1].

Parameters:
  • order (int) – The number of sample points \(n\) minus one. The quadrature rule integrates degree \(2n-1\) polynomials exactly.

  • dim (int) – The dimension of the quadrature region.

  • permute (bool, optional) – Permute the quadrature points according to the cell point orderings (default is True).

Notes

The approximation is given by

\[\int_{-1}^1 f(x) dx \approx \sum f(x_q) w_q\]

with quadrature points \(x_q\) and corresponding weights \(w_q\).

Examples

>>> import felupe as fem
>>> fem.GaussLegendreBoundary(order=2, dim=3).screenshot()
../_images/quadrature_boundary.png
class felupe.quadrature.Triangle(order: int)[source]#

Bases: Scheme

An integration scheme for Triangles.

class felupe.quadrature.Tetrahedron(order: int)[source]#

Bases: Scheme

An integration scheme for Tetrahedrons.