polynomials_on_simplices.calculus.quadrature module

Numerical evaluation of integrals over triangles (or simplices).

For triangle quadrature rules, see D.A. Dunavant, High degree efficient symmetrical Gaussian quadrature rules for the triangle. For tetrahedron quadrature rules, see Yu Jinyun Symmetric Gaussian Quadrature Formulae for Tetrahedronal Regions, and the PHG FEM software: http://lsec.cc.ac.cn/phg/index_en.htm.

gauss_legendre_unit_interval(n, n_digits=10)[source]

Compute quadrature points and weights for Gauss-Legendre quadrature modified from the interval [-1, 1] to the unit interval.

\[\int_0^1 f(x)\,dx \approx \sum_{i=0}^{n-1} w_i f(x_i).\]
Parameters:
  • n – The order of quadrature (number of points and weights).
  • n_digits – Number of significant digits of the points and weights to return.
Returns:

Tuple (x, w) of lists containing the points x and the weights w.

quadrature_interval(f, a, b)[source]

Numerically compute the integral \(\int_a^b f(x) \, dx\).

Parameters:
  • f (Callable f(x)) – Function to integrate.
  • a – Start point of the interval over which we should integrate f.
  • b – End point of the interval over which we should integrate f.
Returns:

Approximate value of the integral.

quadrature_interval_fixed(f, a, b, n)[source]

Numerically compute the integral \(\int_a^b f(x) \, dx\).

Parameters:
  • f (Callable f(x)) – Function to integrate.
  • a – Start point of the interval over which we should integrate f.
  • b – End point of the interval over which we should integrate f.
  • n – Quadrature degree. The result is guaranteed to be exact for all polynomials of degree <= n.
Returns:

Approximate value of the integral.

quadrature_tetrahedron(f, vertices)[source]

Numerically compute the integral \(\int_{T} f(x) \, dx\), for a given tetrahedron T.

Parameters:
  • f (Callable f(x, y, z)) – Function to integrate.
  • vertices – Vertices of the tetrahedron (4 by 3 array).
Returns:

Approximate value of the integral.

quadrature_tetrahedron_fixed(f, vertices, n)[source]

Numerically compute the integral \(\int_{T} f(x) \, dx\), for a given tetrahedron T.

Parameters:
  • f (Callable f(x, y, z)) – Function to integrate.
  • vertices – Vertices of the tetrahedron (4 by 3 array).
  • n – Quadrature degree. The result is guaranteed to be exact for all polynomials of degree <= n.
Returns:

Approximate value of the integral.

quadrature_tetrahedron_midpoint_rule(f, vertices, n)[source]

Numerically compute the integral \(\int_T f \, dx\) using the midpoint rule.

Note

Only useful for testing, since quadrature_tetrahedron() and quadrature_tetrahedron_fixed() are much more efficient.

Parameters:
  • f (Callable f(x, y, z)) – Function to integrate.
  • vertices – Vertices of the tetrahedron to integrate over (4 by 3 array).
  • n – Number of smaller tetrahedrons along each edge of the full tetrahedron, over which the function is approximated with a constant value.
Returns:

Approximate value of the integral.

quadrature_triangle(f, vertices)[source]

Numerically compute the integral \(\int_{T} f(x) \, dx\), for a given triangle T.

Parameters:
  • f (Callable f(x, y)) – Function to integrate.
  • vertices – Vertices of the triangle (3 by 2 array).
Returns:

Approximate value of the integral.

quadrature_triangle_fixed(f, vertices, n)[source]

Numerically compute the integral \(\int_{T} f(x) \, dx\), for a given triangle T.

Parameters:
  • f (Callable f(x, y)) – Function to integrate.
  • vertices – Vertices of the triangle (3 by 2 array).
  • n – Quadrature degree. The result is guaranteed to be exact for all polynomials of degree <= n.
Returns:

Approximate value of the integral.

quadrature_triangle_midpoint_rule(f, vertices, n)[source]

Numerically compute the integral \(\int_T f \, dx\) using the midpoint rule.

Note

Only useful for testing, since quadrature_triangle() and quadrature_triangle_fixed() are much more efficient.

Parameters:
  • f (Callable f(x, y)) – Function to integrate.
  • vertices – Vertices of the triangle to integrate over (3 by 2 array).
  • n – Number of smaller triangles along each edge of the full triangle, over which the function is approximated with a constant value.
Returns:

Approximate value of the integral.

quadrature_unit_interval(f)[source]

Numerically compute the integral \(\int_0^1 f(x) \, dx\).

Parameters:f (Callable f(x)) – Function to integrate.
Returns:Approximate value of the integral.
quadrature_unit_interval_fixed(f, r)[source]

Numerically compute the integral \(\int_0^1 f(x) \, dx\).

Parameters:
  • f (Callable f(x)) – Function to integrate.
  • r – Quadrature degree. The result is guaranteed to be exact for all polynomials of degree <= r.
Returns:

Approximate value of the integral.

quadrature_unit_tetrahedron(f)[source]

Numerically compute the integral \(\int_{\Delta^3_c} f(x) \, dx\).

Parameters:f (Callable f(x, y, z)) – Function to integrate.
Returns:Approximate value of the integral.
quadrature_unit_tetrahedron_fixed(f, n)[source]

Numerically compute the integral \(\int_{\Delta^3_c} f(x) \, dx\).

Parameters:
  • f (Callable f(x, y, z)) – Function to integrate.
  • n – Quadrature degree. The result is guaranteed to be exact for all polynomials of degree <= n.
Returns:

Approximate value of the integral.

quadrature_unit_triangle(f)[source]

Numerically compute the integral \(\int_{\Delta^2_c} f(x) \, dx\).

Parameters:f (Callable f(x, y)) – Function to integrate.
Returns:Approximate value of the integral.
quadrature_unit_triangle_fixed(f, n)[source]

Numerically compute the integral \(\int_{\Delta^2_c} f(x) \, dx\).

Parameters:
  • f (Callable f(x, y)) – Function to integrate.
  • n – Quadrature degree. The result is guaranteed to be exact for all polynomials of degree <= n.
Returns:

Approximate value of the integral.