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
|
An arbitrary-order Gauss-Legendre quadrature rule of dim 1, 2 or 3 on the interval |
|
An arbitrary-order Gauss-Legendre quadrature rule of dim 1, 2 or 3 on the interval |
Triangles and Tetrahedrons
|
alias of |
|
alias of |
|
An integration scheme for Triangles. |
|
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\).
- 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()
- 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()