bundles / scipy 1.17.1 / scipy / integrate / _quadpack_py / nquad
function
scipy.integrate._quadpack_py:nquad
Signature
def nquad ( func , ranges , args = None , opts = None , full_output = False ) Summary
Integration over multiple variables.
Extended Summary
Wraps quad to enable integration over multiple variables. Various options allow improved integration of discontinuous functions, as well as the use of weighted integration, and generally finer control of the integration process.
Parameters
func: {callable, scipy.LowLevelCallable}The function to be integrated. Has arguments of
x0, ..., xn,t0, ..., tm, where integration is carried out overx0, ..., xn, which must be floats. Wheret0, ..., tmare extra arguments passed in args. Function signature should befunc(x0, x1, ..., xn, t0, t1, ..., tm). Integration is carried out in order. That is, integration overx0is the innermost integral, andxnis the outermost.If the user desires improved integration performance, then
fmay be a scipy.LowLevelCallable with one of the signaturesdouble func(int n, double *xx) double func(int n, double *xx, void *user_data)
where
nis the number of variables and args. Thexxarray contains the coordinates and extra arguments.user_datais the data contained in the scipy.LowLevelCallable.ranges: iterable objectEach element of ranges may be either a sequence of 2 numbers, or else a callable that returns such a sequence.
ranges[0]corresponds to integration over x0, and so on. If an element of ranges is a callable, then it will be called with all of the integration arguments available, as well as any parametric arguments. e.g., iffunc = f(x0, x1, x2, t0, t1), thenranges[0]may be defined as either(a, b)or else as(a, b) = range0(x1, x2, t0, t1).args: iterable object, optionalAdditional arguments
t0, ..., tn, required byfunc,ranges, andopts.opts: iterable object or dict, optionalOptions to be passed to quad. May be empty, a dict, or a sequence of dicts or functions that return a dict. If empty, the default options from scipy.integrate.quad are used. If a dict, the same options are used for all levels of integraion. If a sequence, then each element of the sequence corresponds to a particular integration. e.g.,
opts[0]corresponds to integration overx0, and so on. If a callable, the signature must be the same as forranges. The available options together with their default values are:epsabs = 1.49e-08
epsrel = 1.49e-08
limit = 50
points = None
weight = None
wvar = None
wopts = None
For more information on these options, see quad.
full_output: bool, optionalPartial implementation of
full_outputfrom scipy.integrate.quad. The number of integrand function evaluationsnevalcan be obtained by settingfull_output=Truewhen calling nquad.
Returns
result: floatThe result of the integration.
abserr: floatThe maximum of the estimates of the absolute error in the various integration results.
out_dict: dict, optionalA dict containing additional information on the integration.
Notes
For valid results, the integral must converge; behavior for divergent integrals is not guaranteed.
Details of QUADPACK level routines
nquad calls routines from the FORTRAN library QUADPACK. This section provides details on the conditions for each routine to be called and a short description of each routine. The routine called depends on weight, points and the integration limits a and b.
================ ============== ========== ===================== QUADPACK routine `weight` `points` infinite bounds ================ ============== ========== ===================== qagse None No No qagie None No Yes qagpe None Yes No qawoe 'sin', 'cos' No No qawfe 'sin', 'cos' No either `a` or `b` qawse 'alg*' No No qawce 'cauchy' No No ================ ============== ========== =====================
The following provides a short description from [1] for each routine.
qagse
is an integrator based on globally adaptive interval subdivision in connection with extrapolation, which will eliminate the effects of integrand singularities of several types. The integration is is performed using a 21-point Gauss-Kronrod quadrature within each subinterval.
qagie
handles integration over infinite intervals. The infinite range is mapped onto a finite interval and subsequently the same strategy as in
QAGSis applied.qagpe
serves the same purposes as QAGS, but also allows the user to provide explicit information about the location and type of trouble-spots i.e. the abscissae of internal singularities, discontinuities and other difficulties of the integrand function.
qawoe
is an integrator for the evaluation of or over a finite interval [a,b], where and are specified by the user. The rule evaluation component is based on the modified Clenshaw-Curtis technique
An adaptive subdivision scheme is used in connection with an extrapolation procedure, which is a modification of that in
QAGSand allows the algorithm to deal with singularities in .qawfe
calculates the Fourier transform or for user-provided and . The procedure of
QAWOis applied on successive finite intervals, and convergence acceleration by means of the -algorithm is applied to the series of integral approximations.qawse
approximate , with where with , where may be one of the following functions: , , , .
The user specifies , and the type of the function . A globally adaptive subdivision strategy is applied, with modified Clenshaw-Curtis integration on those subintervals which contain
aorb.qawce
compute where the integral must be interpreted as a Cauchy principal value integral, for user specified and . The strategy is globally adaptive. Modified Clenshaw-Curtis integration is used on those intervals containing the point .
Array API Standard Support
nquad 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 ⛔ ⛔ Dask ⛔ n/a ==================== ==================== ====================
See
dev-arrayapifor more information.
Examples
Calculate .. math:: \int^{1}_{-0.15} \int^{0.8}_{0.13} \int^{1}_{-1} \int^{1}_{0} f(x_0, x_1, x_2, x_3) \,dx_0 \,dx_1 \,dx_2 \,dx_3 , where .. math:: f(x_0, x_1, x_2, x_3) = \begin{cases} x_0^2+x_1 x_2-x_3^3+ \sin{x_0}+1 & (x_0-0.2 x_3-0.5-0.25 x_1 > 0) \\ x_0^2+x_1 x_2-x_3^3+ \sin{x_0}+0 & (x_0-0.2 x_3-0.5-0.25 x_1 \leq 0) \end{cases} .import numpy as np from scipy import integrate func = lambda x0,x1,x2,x3 : x0**2 + x1*x2 - x3**3 + np.sin(x0) + ( 1 if (x0-.2*x3-.5-.25*x1>0) else 0) def opts0(*args, **kwargs): return {'points':[0.2*args[2] + 0.5 + 0.25*args[0]]}✓
integrate.nquad(func, [[0,1], [-1,1], [.13,.8], [-.15,1]], opts=[opts0,{},{},{}], full_output=True)✗
def func2(x0, x1, x2, t0, t1): return x0*x2**2 + np.sin(x1) + 1 + (1 if x0+t1*x1-t0>0 else 0) def lim0(x1, x2, t0, t1): return [t0*x1 + t1*x2 - 1, t0*x1 + t1*x2 + 1] def lim1(x2, t0, t1): return [x2 + t0**2*t1**3 - 1, x2 + t0**2*t1**3 + 1] def lim2(t0, t1): return [t0 + t1 - 1, t0 + t1 + 1] def opts0(x1, x2, t0, t1): return {'points' : [t0 - t1*x1]} def opts1(x2, t0, t1): return {} def opts2(t0, t1): return {}✓
integrate.nquad(func2, [lim0, lim1, lim2], args=(0,1), opts=[opts0, opts1, opts2])✗
See also
- dblquad
double and triple integrals
- fixed_quad
fixed-order Gaussian quadrature
- quad
1-D numerical integration
- tplquad
double and triple integrals
Aliases
-
scipy.integrate.nquad