Skip to content
forked from sigma-py/quadpy

Numerical integration (quadrature, cubature) in Python

License

Notifications You must be signed in to change notification settings

basnijholt/quadpy

 
 

Repository files navigation

quadpy

Your one-stop shop for numerical integration in Python.

CircleCI codecov PyPi Version GitHub stars

Hundreds of numerical integration schemes for line segments, circles, disks, triangles, quadrilaterals, spheres, balls, tetrahedra, hexahedra, wedges, pyramids, n-spheres, n-balls, n-cubes, n-simplices, and the 2D/3D/nD spaces with weight functions exp(-r) and exp(-r2).

To numerically integrate any function over any given triangle, do

import numpy
import quadpy

def f(x):
    return numpy.sin(x[0]) * numpy.sin(x[1])

triangle = numpy.array([[0.0, 0.0], [1.0, 0.0], [0.7, 0.5]])

val = quadpy.triangle.integrate(f, triangle, quadpy.triangle.Strang(9))

This uses Strang's rule of degree 6.

quadpy is fully vectorized, so if you like to compute the integral of a function on many domains at once, you can provide them all in one integrate() call, e.g.,

triangles = numpy.stack([
    [[0.0, 0.0], [1.0, 0.0], [0.0, 1.0]],
    [[1.2, 0.6], [1.3, 0.7], [1.4, 0.8]],
    [[26.0, 31.0], [24.0, 27.0], [33.0, 28]],
    [[0.1, 0.3], [0.4, 0.4], [0.7, 0.1]],
    [[8.6, 6.0], [9.4, 5.6], [7.5, 7.4]]
    ], axis=-2)

The same goes for functions with vectorized output, e.g.,

def f(x):
    return [numpy.sin(x[0]), numpy.sin(x[1])]

More examples under test/examples_test.py.

Adaptive quadrature

quadpy can do adaptive quadrature for certain domains. Again, everything is fully vectorized, so you can provide multiple intervals and vector-valued functions.

Line segments

val, error_estimate = quadpy.line_segment.adaptive_integrate(
        lambda x: x * sin(5 * x),
        [0.0, pi],
        1.0e-10
        )

Triangles

val, error_estimate = quadpy.triangle.adaptive_integrate(
        lambda x: x[0] * sin(5 * x[1]),
        [[0.0, 0.0], [1.0, 0.0], [0.0, 1.0]],
        1.0e-10
        )

ProTip: You can provide many triangles that together form a domain to get an approximation of the integral over the domain.

Schemes

Line segment

  • Chebyshev-Gauss (both variants, arbitrary degree)
  • Clenshaw-Curtis (after Waldvogel, arbitrary degree)
  • Fejér-type-1 (after Waldvogel, arbitrary degree)
  • Fejér-type-2 (after Waldvogel, arbitrary degree)
  • Gauss-Hermite (via NumPy, arbitrary degree)
  • Gauss-Laguerre (via NumPy, arbitrary degree)
  • Gauss-Legendre (via NumPy, arbitrary degree)
  • Gauss-Lobatto (arbitrary degree)
  • Gauss-Kronrod (after Laurie, arbitrary degree)
  • Gauss-Patterson (7 schemes up to degree 191)
  • Gauss-Radau (arbitrary degree)
  • closed Newton-Cotes (arbitrary degree)
  • open Newton-Cotes (arbitrary degree)

You can use orthopy to generate Gauss formulas for your own weight functions.

Example:

val = quadpy.line_segment.integrate(
    lambda x: numpy.exp(x),
    [0.0, 1.0],
    quadpy.line_segment.GaussPatterson(5)
    )

Circle

  • Krylov (1959, arbitrary degree)

Example:

val = quadpy.circle.integrate(
    lambda x: numpy.exp(x[0]),
    [0.0, 0.0], 1.0,
    quadpy.circle.Krylov(7)
    )

Triangle

Apart from the classical centroid, vertex, and seven-point schemes we have

Example:

val = quadpy.triangle.integrate(
    lambda x: numpy.exp(x[0]),
    [[0.0, 0.0], [1.0, 0.0], [0.5, 0.7]],
    quadpy.triangle.XiaoGimbutas(5)
    )

Disk

Example:

val = quadpy.disk.integrate(
    lambda x: numpy.exp(x[0]),
    [0.0, 0.0], 1.0,
    quadpy.disk.Lether(6)
    )

Quadrilateral

Example:

val = quadpy.quadrilateral.integrate(
    lambda x: numpy.exp(x[0]),
    [[[0.0, 0.0], [1.0, 0.0]], [[0.0, 1.0], [1.0, 1.0]]],
    quadpy.quadrilateral.Stroud('C2 7-2')
    )

The points are specified in an array of shape (2, 2, ...) such that arr[0][0] is the lower left corner, arr[1][1] the upper right. If your quadrilateral has its sides aligned with the coordinate axes, you can use the convenience function

quadpy.quadrilateral.rectangle_points([x0, x1], [y0, y1])

to generate the array.

2D space with weight function exp(-r)

Example:

val = quadpy.e2r.integrate(
    lambda x: numpy.exp(x[0]),
    quadpy.e2r.RabinowitzRichter(5)
    )

2D space with weight function exp(-r2)

Example:

val = quadpy.e2r2.integrate(
    lambda x: numpy.exp(x[0]),
    quadpy.e2r2.RabinowitzRichter(3)
    )

Sphere

  • via Stroud (1971):
  • Lebedev (1976, 32 schemes up to degree 131)
  • Heo-Xu (2001, 27 schemes up to degree 39, single-precision)

Example:

val = quadpy.sphere.integrate(
    lambda x: numpy.exp(x[0]),
    [0.0, 0.0, 0.0], 1.0,
    quadpy.sphere.Lebedev(19)
    )

Integration on the sphere can also be done for function defined in spherical coordinates:

val = quadpy.sphere.integrate_spherical(
    lambda phi_theta: numpy.sin(phi_theta[0])**2 * numpy.sin(phi_theta[1]),
    radius=1.0,
    rule=quadpy.sphere.Lebedev(19)
    )

Note that phi_theta[0] is the azimuthal, phi_theta[1] the polar angle here.

Ball

  • Hammer-Stroud (1958, 6 schemes up to degree 7)
  • via: Stroud (1971):
    • Ditkin (1948, 3 schemes up to degree 7)
    • Mysovskih (1964, degree 7)
    • 2 schemes up to degree 14

Example:

val = quadpy.ball.integrate(
    lambda x: numpy.exp(x[0]),
    [0.0, 0.0, 0.0], 1.0,
    quadpy.ball.HammerStroud('14-3a')
    )

Tetrahedron

Example:

val = quadpy.tetrahedron.integrate(
    lambda x: numpy.exp(x[0]),
    [[0.0, 0.0, 0.0], [1.0, 0.0, 0.0], [0.5, 0.7, 0.0], [0.3, 0.9, 1.0]],
    quadpy.tetrahedron.Keast(10)
    )

Hexahedron

Example:

val = quadpy.hexahedron.integrate(
    lambda x: numpy.exp(x[0]),
    quadpy.hexahedron.cube_points([0.0, 1.0], [-0.3, 0.4], [1.0, 2.1]),
    quadpy.hexahedron.Product(quadpy.line_segment.NewtonCotesClosed(3))
    )

Pyramid

  • Felippa (9 schemes up to degree 5)

Example:

val = quadpy.pyramid.integrate(
    lambda x: numpy.exp(x[0]),
    [
      [0.0, 0.0, 0.0], [1.0, 0.0, 0.0], [0.5, 0.7, 0.0], [0.3, 0.9, 0.0],
      [0.0, 0.1, 1.0],
    ],
    quadpy.pyramid.Felippa(5)
    )

Wedge

  • Felippa (6 schemes up to degree 6)

Example:

val = quadpy.wedge.integrate(
    lambda x: numpy.exp(x[0]),
    [
      [[0.0, 0.0, 0.0], [1.0, 0.0, 0.0], [0.5, 0.7, 0.0]],
      [[0.0, 0.0, 1.0], [1.0, 0.0, 1.0], [0.5, 0.7, 1.0]],
    ],
    quadpy.wedge.Felippa(3)
    )

3D space with weight function exp(-r)

Example:

val = quadpy.e2r.integrate(
    lambda x: numpy.exp(x[0]),
    quadpy.e2r.StroudSecrest('IX')
    )

3D space with weight function exp(-r2)

Example:

val = quadpy.e2r2.integrate(
    lambda x: numpy.exp(x[0]),
    quadpy.e2r2.RabinowitzRichter(3)
    )

n-Simplex

Example:

dim = 4
val = quadpy.simplex.integrate(
    lambda x: numpy.exp(x[0]),
    numpy.array([
        [0.0, 0.0, 0.0, 0.0],
        [1.0, 2.0, 0.0, 0.0],
        [0.0, 1.0, 0.0, 0.0],
        [0.0, 3.0, 1.0, 0.0],
        [0.0, 0.0, 4.0, 1.0],
        ]),
    quadpy.simplex.GrundmannMoeller(dim, 3)
    )

n-Sphere

  • via Stroud (1971):
    • Stroud (1967, degree 7)
    • Stroud (1969, 3 <= n <= 16, degree 11)
    • 6 schemes up to degree 5
  • Dobrodeev (1978, n >= 2, degree 5)

Example:

dim = 4
quadpy.nsphere.integrate(
    lambda x: numpy.exp(x[0]),
    numpy.zeros(dim), 1.0,
    quadpy.nsphere.Dobrodeev1978(dim)
    )

n-Ball

  • Dobrodeev (1970, n >= 3, degree 7)
  • via Stroud (1971):
    • Stroud (1957, degree 2)
    • Hammer-Stroud (1958, 2 schemes up to degree 5)
    • Stroud (1966, 4 schemes of degree 5)
    • Stroud (1967, 4 <= n <= 7, 2 schemes of degree 5)
    • Stroud (1967, n >= 3, 3 schemes of degree 7)
    • Stenger (1967, 6 schemes up to degree 11)
  • Dobrodeev (1978, 2 <= n <= 20, degree 5)

Example:

dim = 4
quadpy.nball.integrate(
    lambda x: numpy.exp(x[0]),
    numpy.zeros(dim), 1.0,
    quadpy.nball.Dobrodeev1970(dim)
    )

n-Cube

Example:

dim = 4
quadpy.ncube.integrate(
    lambda x: numpy.exp(x[0]),
    quadpy.ncube.ncube_points(
        [0.0, 1.0], [0.1, 0.9], [-1.0, 1.0], [-1.0, -0.5]
        ),
    quadpy.ncube.Stroud(dim, 'Cn 3-3')
    )

nD space with weight function exp(-r)

Example:

dim = 4
val = quadpy.enr.integrate(
    lambda x: numpy.exp(x[0]),
    quadpy.enr.Stroud(dim, '5-4')
    )

nD space with weight function exp(-r2)

  • via Stroud (1971):
    • Stroud-Secrest (1963, 4 schemes up to degree 5)
    • Stroud (1967, 2 schemes of degree 5)
    • Stroud (1967, 3 schemes of degree 7)
    • Stenger (1971, 6 schemes up to degree 11, varying dimensionality restrictions)
    • 5 schemes up to degree 5

Example:

dim = 4
val = quadpy.enr2.integrate(
    lambda x: numpy.exp(x[0]),
    quadpy.enr2.Stroud(dim, '5-2')
    )

Installation

quadpy is available from the Python Package Index, so with

pip install -U quadpy

you can install/upgrade.

Testing

To run the tests, just check out this repository and type

MPLBACKEND=Agg pytest

Distribution

To create a new release

  1. bump the __version__ number,

  2. publish to PyPi and GitHub:

    $ make publish
    

License

quadpy is published under the MIT license.

About

Numerical integration (quadrature, cubature) in Python

Resources

License

Stars

Watchers

Forks

Packages

No packages published

Languages

  • Python 99.9%
  • Makefile 0.1%