bundles / scipy latest / scipy / integrate / _cubature / cubature
function
scipy.integrate._cubature:cubature
Signature
def cubature ( f , a , b , * , rule = gk21 , rtol = 1e-08 , atol = 0 , max_subdivisions = 10000 , args = () , workers = 1 , points = None ) Summary
Adaptive cubature of multidimensional array-valued function.
Extended Summary
Given an arbitrary integration rule, this function returns an estimate of the integral to the requested tolerance over the region defined by the arrays a and b specifying the corners of a hypercube.
Convergence is not guaranteed for all integrals.
Parameters
f: callableFunction to integrate.
fmust have the signaturef(x : ndarray, *args) -> ndarrayfshould accept arraysxof shape(npoints, ndim)and output arrays of shape
(npoints, output_dim_1, ..., output_dim_n)In this case, cubature will return arrays of shape
(output_dim_1, ..., output_dim_n)a, b: array_likeLower and upper limits of integration as 1D arrays specifying the left and right endpoints of the intervals being integrated over. Limits can be infinite.
rule: str, optionalRule used to estimate the integral. If passing a string, the options are "gauss-kronrod" (21 node), or "genz-malik" (degree 7). If a rule like "gauss-kronrod" is specified for an
n-dim integrand, the corresponding Cartesian product rule is used. "gk21", "gk15" are also supported for compatibility with quad_vec. See Notes.rtol, atol: float, optionalRelative and absolute tolerances. Iterations are performed until the error is estimated to be less than
atol + rtol * abs(est). Herertolcontrols relative accuracy (number of correct digits), whileatolcontrols absolute accuracy (number of correct decimal places). To achieve the desiredrtol, setatolto be smaller than the smallest value that can be expected fromrtol * abs(y)so thatrtoldominates the allowable error. Ifatolis larger thanrtol * abs(y)the number of correct digits is not guaranteed. Conversely, to achieve the desiredatol, setrtolsuch thatrtol * abs(y)is always smaller thanatol. Default values are 1e-8 forrtoland 0 foratol.max_subdivisions: int, optionalUpper bound on the number of subdivisions to perform. Default is 10,000.
args: tuple, optionalAdditional positional args passed to
f, if any.workers: int or map-like callable, optionalIf
workersis an integer, part of the computation is done in parallel subdivided to this many tasks (usingpython:multiprocessing.pool.Pool). Supply-1to use all cores available to the Process. Alternatively, supply a map-like callable, such aspython:multiprocessing.pool.Pool.mapfor evaluating the population in parallel. This evaluation is carried out asworkers(func, iterable).points: list of array_like, optionalList of points to avoid evaluating
fat, under the condition that the rule being used does not evaluatefon the boundary of a region (which is the case for all Genz-Malik and Gauss-Kronrod rules). This can be useful iffhas a singularity at the specified point. This should be a list of array-likes where each element has lengthndim. Default is empty. See Examples.
Returns
res: objectObject containing the results of the estimation. It has the following attributes:
estimate
estimate
error
error
status
status
subdivisions
subdivisions
atol,rtol
atol, rtol
regions: list of object
List of objects containing the estimates of the integral over smaller regions of the domain.
Each object in
regionshas the following attributes:a,b
a, b
estimate
estimate
error
error
Notes
The algorithm uses a similar algorithm to quad_vec, which itself is based on the implementation of QUADPACK's DQAG* algorithms, implementing global error control and adaptive subdivision.
The source of the nodes and weights used for Gauss-Kronrod quadrature can be found in [1], and the algorithm for calculating the nodes and weights in Genz-Malik cubature can be found in [2].
The rules currently supported via the rule argument are:
"gauss-kronrod", 21-node Gauss-Kronrod"genz-malik", n-node Genz-Malik
If using Gauss-Kronrod for an n-dim integrand where n > 2, then the corresponding Cartesian product rule will be found by taking the Cartesian product of the nodes in the 1D case. This means that the number of nodes scales exponentially as 21^n in the Gauss-Kronrod case, which may be problematic in a moderate number of dimensions.
Genz-Malik is typically less accurate than Gauss-Kronrod but has much fewer nodes, so in this situation using "genz-malik" might be preferable.
Infinite limits are handled with an appropriate variable transformation. Assuming a = [a_1, ..., a_n] and b = [b_1, ..., b_n]:
If and , the i-th integration variable will use the transformation and .
If and , the i-th integration variable will use the transformation and .
If and , the i-th integration variable will use the transformation and .
Array API Standard Support
cubature has experimental support for Python Array API Standard compatible backends in addition to NumPy. Please consider testing these features by setting an environment variable SCIPY_ARRAY_API=1 and providing CuPy, PyTorch, JAX, or Dask arrays as array arguments. The following combinations of backend and device (or other capability) are supported.
==================== ==================== ==================== Library CPU GPU ==================== ==================== ==================== NumPy ✅ n/a CuPy n/a ✅ PyTorch ✅ ✅ JAX ⚠️ no JIT ⚠️ no JIT Dask ⚠️ computes graph n/a ==================== ==================== ====================
See
dev-arrayapifor more information.
Examples
**1D integral with vector output**: .. math:: \int^1_0 \mathbf f(x) \text dx Where ``f(x) = x^n`` and ``n = np.arange(10)`` is a vector. Since no rule is specified, the default "gk21" is used, which corresponds to Gauss-Kronrod integration with 21 nodes.import numpy as np from scipy.integrate import cubature def f(x, n): # Make sure x and n are broadcastable return x[:, np.newaxis]**n[np.newaxis, :] res = cubature( f, a=[0], b=[1], args=(np.arange(10),), )✓
res.estimate
✗import numpy as np from scipy.integrate import cubature def f(x, r, alphas): # f(x) = cos(2*pi*r + alphas @ x) # Need to allow r and alphas to be arbitrary shape npoints, ndim = x.shape[0], x.shape[-1] alphas = alphas[np.newaxis, ...] x = x.reshape(npoints, *([1]*(len(alphas.shape) - 1)), ndim) return np.cos(2*np.pi*r + np.sum(alphas * x, axis=-1)) rng = np.random.default_rng() r, alphas = rng.random((2, 3)), rng.random((2, 3, 7)) res = cubature( f=f, a=np.array([0, 0, 0, 0, 0, 0, 0]), b=np.array([1, 1, 1, 1, 1, 1, 1]), rtol=1e-5, rule="genz-malik", args=(r, alphas), )✓
res.estimate
✗from concurrent.futures import ThreadPoolExecutor with ThreadPoolExecutor() as executor: res = cubature( f=f, a=np.array([0, 0, 0, 0, 0, 0, 0]), b=np.array([1, 1, 1, 1, 1, 1, 1]), rtol=1e-5, rule="genz-malik", args=(r, alphas), workers=executor.map, )✓
res.estimate
✗def gaussian(x): return np.exp(-np.sum(x**2, axis=-1)) res = cubature(gaussian, [-np.inf, -np.inf], [np.inf, np.inf])✓
res.estimate
✗def sinc(x): return np.sin(x)/x res = cubature(sinc, [-1], [1], points=[[0]])✓
res.estimate
✗Aliases
-
scipy.integrate.cubature