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 dimension 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].

GaussLobatto(order, dim)

An arbitrary-order Gauss-Lobatto quadrature rule of dimension 1, 2 or 3 on the interval \([-1, 1]\).

GaussLobattoBoundary(order, dim)

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

Triangles and Tetrahedrons

quadrature.Triangle(order)

A quadrature scheme suitable for Triangles of order 1, 2 or 3 on the interval \([0, 1]\).

quadrature.Tetrahedron(order)

A quadrature scheme suitable for Tetrahedrons of order 1, 2 or 3 on the interval \([0, 1]\).

TriangleQuadrature

TetrahedronQuadrature

Sphere

BazantOh([n])

Quadrature scheme for a numeric integration of the surface of a unit sphere [1]_.

Detailed API Reference

class felupe.quadrature.Scheme(points, weights)[ソース]#

ベースクラス: 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.

メモ

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, point_size=20, weighted=False, **kwargs)[ソース]#

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

参考

felupe.Scene.plot

Plot method of a scene.

screenshot(*args, filename=None, transparent_background=None, scale=None, **kwargs)[ソース]#

Take a screenshot of the quadrature.

参考

pyvista.Plotter.screenshot

Take a screenshot of a PyVista plotter.

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

ベースクラス: Scheme

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

パラメータ:
  • 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). This is supported for two and three dimensions as well as first and second order schemes. Otherwise this flag is silently ignored.

メモ

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\).

サンプル

>>> import felupe as fem
>>>
>>> element = fem.Line()
>>> quadrature = fem.GaussLegendre(order=2, dim=1)
>>> quadrature.plot(
...     plotter=element.plot(add_point_labels=False, show_points=False),
...     weighted=True,
... ).show()
../_images/quadrature-2a44ff09216e9348_00_00.png
>>> import felupe as fem
>>>
>>> element = fem.QuadraticQuad()
>>> quadrature = fem.GaussLegendre(order=2, dim=2)
>>> quadrature.plot(
...     plotter=element.plot(add_point_labels=False, show_points=False),
...     weighted=True,
... ).show()
../_images/quadrature-1db5ecbba5ed37a1_00_00.png
>>> import felupe as fem
>>>
>>> element = fem.QuadraticHexahedron()
>>> quadrature = fem.GaussLegendre(order=2, dim=3)
>>> quadrature.plot(
...     plotter=element.plot(add_point_labels=False, show_points=False),
...     weighted=True,
... ).show()
../_images/quadrature-1acad856983d53b7_00_00.png
>>> import felupe as fem
>>>
>>> element = fem.ArbitraryOrderLagrangeElement(order=5, dim=2)
>>> quadrature = fem.GaussLegendre(order=5, dim=2)
>>> quadrature.plot(
...     plotter=element.plot(add_point_labels=False, show_points=False),
...     weighted=True,
... ).show()
../_images/quadrature-b98d30fe6cdf6b94_00_00.png
inv()[ソース]#

Return the inverse quadrature scheme.

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

ベースクラス: GaussLegendre

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

パラメータ:
  • 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).

メモ

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\).

サンプル

>>> import felupe as fem
>>>
>>> element = fem.QuadraticQuad()
>>> quadrature = fem.GaussLegendreBoundary(order=2, dim=2)
>>> quadrature.plot(
...     plotter=element.plot(add_point_labels=False, show_points=False),
...     weighted=True,
... ).show()
../_images/quadrature-0aa5a6f2fbd96b0c_00_00.png
>>> import felupe as fem
>>>
>>> element = fem.QuadraticHexahedron()
>>> quadrature = fem.GaussLegendreBoundary(order=2, dim=3)
>>> quadrature.plot(
...     plotter=element.plot(add_point_labels=False, show_points=False),
...     weighted=True,
... ).show()
../_images/quadrature-10a42b0f0c950bdc_00_00.png
>>> import felupe as fem
>>>
>>> element = fem.ArbitraryOrderLagrangeElement(order=5, dim=3)
>>> quadrature = fem.GaussLegendreBoundary(order=5, dim=3)
>>> quadrature.plot(
...     plotter=element.plot(add_point_labels=False, show_points=False),
...     weighted=True,
... ).show()
../_images/quadrature-9f84fb790ed91c2b_00_00.png
class felupe.GaussLobatto(order: int, dim: int)[ソース]#

ベースクラス: Scheme

An arbitrary-order Gauss-Lobatto quadrature rule of dimension 1, 2 or 3 on the interval \([-1, 1]\).

パラメータ:
  • order (int) -- The number of sample points \(n\) minus two. The quadrature rule integrates degree \(2n-3\) 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). This is supported for two and three dimensions as well as first and second order schemes. Otherwise this flag is silently ignored.

メモ

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\).

サンプル

注釈

The quadrature points are not weighted in the plots because the end-points would be almost invisible for the Gauss-Lobatto quadrature rule.

>>> import felupe as fem
>>>
>>> element = fem.Line()
>>> quadrature = fem.GaussLobatto(order=2, dim=1)
>>> quadrature.plot(
...     plotter=element.plot(add_point_labels=False, show_points=False),
...     weighted=False,
... ).show()
../_images/quadrature-151cbbc28d36d292_00_00.png
>>> import felupe as fem
>>>
>>> element = fem.QuadraticQuad()
>>> quadrature = fem.GaussLobatto(order=2, dim=2)
>>> quadrature.plot(
...     plotter=element.plot(add_point_labels=False, show_points=False),
...     weighted=False,
... ).show()
../_images/quadrature-851f588ebef23299_00_00.png
>>> import felupe as fem
>>>
>>> element = fem.QuadraticHexahedron()
>>> quadrature = fem.GaussLobatto(order=2, dim=3)
>>> quadrature.plot(
...     plotter=element.plot(add_point_labels=False, show_points=False),
...     weighted=False,
... ).show()
../_images/quadrature-22cbf00c7427ef4c_00_00.png
>>> import felupe as fem
>>>
>>> element = fem.ArbitraryOrderLagrangeElement(order=5, dim=2)
>>> quadrature = fem.GaussLobatto(order=5, dim=2)
>>> quadrature.plot(
...     plotter=element.plot(add_point_labels=False, show_points=False),
...     weighted=False,
... ).show()
../_images/quadrature-e014e05812b7232f_00_00.png
class felupe.GaussLobattoBoundary(order: int, dim: int)[ソース]#

ベースクラス: GaussLobatto

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

パラメータ:
  • order (int) -- The number of sample points \(n\) minus two. The quadrature rule integrates degree \(2n-3\) 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).

メモ

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\).

サンプル

注釈

The quadrature points are not weighted in the plots because the end-points would be almost invisible for the Gauss-Lobatto quadrature rule.

>>> import felupe as fem
>>>
>>> element = fem.QuadraticQuad()
>>> quadrature = fem.GaussLobattoBoundary(order=2, dim=2)
>>> quadrature.plot(
...     plotter=element.plot(add_point_labels=False, show_points=False),
...     weighted=False,
... ).show()
../_images/quadrature-79070e9ec3c56ed2_00_00.png
>>> import felupe as fem
>>>
>>> element = fem.QuadraticHexahedron()
>>> quadrature = fem.GaussLobattoBoundary(order=2, dim=3)
>>> quadrature.plot(
...     plotter=element.plot(add_point_labels=False, show_points=False),
...     weighted=False,
... ).show()
../_images/quadrature-510cc264fbde0c9b_00_00.png
>>> import felupe as fem
>>>
>>> element = fem.ArbitraryOrderLagrangeElement(order=5, dim=3)
>>> quadrature = fem.GaussLobattoBoundary(order=5, dim=3)
>>> quadrature.plot(
...     plotter=element.plot(add_point_labels=False, show_points=False),
...     weighted=False,
... ).show()
../_images/quadrature-833f4be4a2632c4c_00_00.png
class felupe.quadrature.Triangle(order: int)[ソース]#

ベースクラス: Scheme

A quadrature scheme suitable for Triangles of order 1, 2 or 3 on the interval \([0, 1]\).

パラメータ:

order (int) -- The order of the quadrature scheme.

メモ

The approximation is given by

\[\int_A f(x) dA \approx \sum f(x_q) w_q\]

with quadrature points \(x_q\) and corresponding weights \(w_q\) [1]_ [2].

サンプル

>>> import felupe as fem
>>>
>>> element = fem.Triangle()
>>> quadrature = fem.TriangleQuadrature(order=1)
>>> quadrature.plot(
...     plotter=element.plot(add_point_labels=False, show_points=False),
...     weighted=True,
... ).show()
../_images/quadrature-ccc70f5d403152f1_00_00.png
>>> import felupe as fem
>>>
>>> element = fem.QuadraticTriangle()
>>> quadrature = fem.TriangleQuadrature(order=2)
>>> quadrature.plot(
...     plotter=element.plot(add_point_labels=False, show_points=False),
...     weighted=True,
... ).show()
../_images/quadrature-ae42ca563d597536_00_00.png
>>> import felupe as fem
>>>
>>> element = fem.QuadraticTriangle()
>>> quadrature = fem.TriangleQuadrature(order=5)
>>> quadrature.plot(
...     plotter=element.plot(add_point_labels=False, show_points=False),
...     weighted=True,
... ).show()
../_images/quadrature-2b87b96fe61e3997_00_00.png

参照

class felupe.quadrature.Tetrahedron(order: int)[ソース]#

ベースクラス: Scheme

A quadrature scheme suitable for Tetrahedrons of order 1, 2 or 3 on the interval \([0, 1]\).

パラメータ:

order (int) -- The order of the quadrature scheme.

メモ

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\) [1]_.

サンプル

>>> import felupe as fem
>>>
>>> element = fem.Tetra()
>>> quadrature = fem.TetrahedronQuadrature(order=1)
>>> quadrature.plot(
...     plotter=element.plot(add_point_labels=False, show_points=False),
...     weighted=True,
... ).show()
../_images/quadrature-57170068689e4120_00_00.png
>>> import felupe as fem
>>>
>>> element = fem.QuadraticTetra()
>>> quadrature = fem.TetrahedronQuadrature(order=2)
>>> quadrature.plot(
...     plotter=element.plot(add_point_labels=False, show_points=False),
...     weighted=True,
... ).show()
../_images/quadrature-60935262927cc49c_00_00.png
>>> import felupe as fem
>>>
>>> element = fem.QuadraticTetra()
>>> quadrature = fem.TetrahedronQuadrature(order=5)
>>> quadrature.plot(
...     plotter=element.plot(add_point_labels=False, show_points=False),
...     weighted=True,
... ).show()
../_images/quadrature-958779c285dd11cb_00_00.png

参照

class felupe.BazantOh(n: int = 21)[ソース]#

ベースクラス: Scheme

Quadrature scheme for a numeric integration of the surface of a unit sphere [1]_.

パラメータ:

n (int, optional) -- The number of quadrature points (default is 21).

メモ

The approximation is given by

\[\int_{\partial V} f(x) dA \approx \sum f(x_q) w_q\]

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

サンプル

>>> import felupe as fem
>>> import pyvista as pv
>>>
>>> plotter = pv.Plotter()
>>> sphere = pv.Sphere(radius=1).clip(normal="z", invert=False)
>>> actor = plotter.add_mesh(sphere, opacity=0.3, color="white")
>>>
>>> quadrature = fem.BazantOh(n=21)
>>> quadrature.plot(plotter=plotter, weighted=True).show()
../_images/quadrature-572ac8783a70dfb9_00_00.png

参照