bundles / scipy 1.17.1 / scipy / integrate / _quadpack_py / quad
function
scipy.integrate._quadpack_py:quad
Signature
def quad ( func , a , b , args = () , full_output = 0 , epsabs = 1.49e-08 , epsrel = 1.49e-08 , limit = 50 , points = None , weight = None , wvar = None , wopts = None , maxp1 = 50 , limlst = 50 , complex_func = False ) Summary
Compute a definite integral.
Extended Summary
Integrate func from a to b (possibly infinite interval) using a technique from the Fortran library QUADPACK.
Parameters
func: {function, scipy.LowLevelCallable}A Python function or method to integrate. If
functakes many arguments, it is integrated along the axis corresponding to the first argument.If the user desires improved integration performance, then
fmay be a scipy.LowLevelCallable with one of the signaturesdouble func(double x) double func(double x, void *user_data) double func(int n, double *xx) double func(int n, double *xx, void *user_data)
The
user_datais the data contained in the scipy.LowLevelCallable. In the call forms withxx,nis the length of thexxarray which containsxx[0] == xand the rest of the items are numbers contained in theargsargument of quad.In addition, certain ctypes call signatures are supported for backward compatibility, but those should not be used in new code.
a: floatLower limit of integration (use -numpy.inf for -infinity).
b: floatUpper limit of integration (use numpy.inf for +infinity).
args: tuple, optionalExtra arguments to pass to
func.full_output: int, optionalNon-zero to return a dictionary of integration information. If non-zero, warning messages are also suppressed and the message is appended to the output tuple.
complex_func: bool, optionalIndicate if the function's (
func) return type is real (complex_func=False: default) or complex (complex_func=True). In both cases, the function's argument is real. If full_output is also non-zero, the infodict,message, andexplainfor the real and complex components are returned in a dictionary with keys "real output" and "imag output".
Returns
y: floatThe integral of func from
atob.abserr: floatAn estimate of the absolute error in the result.
infodict: dictA dictionary containing additional information.
: messageA convergence message.
: explainAppended only with 'cos' or 'sin' weighting and infinite integration limits, it contains an explanation of the codes in infodict['ierlst']
Other Parameters
epsabs: float or int, optionalAbsolute error tolerance. Default is 1.49e-8. quad tries to obtain an accuracy of
abs(i-result) <= max(epsabs, epsrel*abs(i))wherei= integral offuncfromatob, andresultis the numerical approximation. Seeepsrelbelow.epsrel: float or int, optionalRelative error tolerance. Default is 1.49e-8. If
epsabs <= 0,epsrelmust be greater than both 5e-29 and50 * (machine epsilon). Seeepsabsabove.limit: float or int, optionalAn upper bound on the number of subintervals used in the adaptive algorithm.
points: (sequence of floats,ints), optionalA sequence of break points in the bounded integration interval where local difficulties of the integrand may occur (e.g., singularities, discontinuities). The sequence does not have to be sorted. Note that this option cannot be used in conjunction with
weight.weight: float or int, optionalString indicating weighting function. Full explanation for this and the remaining arguments can be found below.
wvar: optionalVariables for use with weighting functions.
wopts: optionalOptional input for reusing Chebyshev moments.
maxp1: float or int, optionalAn upper bound on the number of Chebyshev moments.
limlst: int, optionalUpper bound on the number of cycles (>=3) for use with a sinusoidal weighting and an infinite end-point.
Notes
For valid results, the integral must converge; behavior for divergent integrals is not guaranteed.
Extra information for quad() inputs and outputs
If full_output is non-zero, then the third output argument (infodict) is a dictionary with entries as tabulated below. For infinite limits, the range is transformed to (0,1) and the optional outputs are given with respect to this transformed range. Let M be the input argument limit and let K be infodict['last']. The entries are:
'neval'
The number of function evaluations.
'last'
The number, K, of subintervals produced in the subdivision process.
'alist'
A rank-1 array of length M, the first K elements of which are the left end points of the subintervals in the partition of the integration range.
'blist'
A rank-1 array of length M, the first K elements of which are the right end points of the subintervals.
'rlist'
A rank-1 array of length M, the first K elements of which are the integral approximations on the subintervals.
'elist'
A rank-1 array of length M, the first K elements of which are the moduli of the absolute error estimates on the subintervals.
'iord'
A rank-1 integer array of length M, the first L elements of which are pointers to the error estimates over the subintervals with
L=KifK<=M/2+2orL=M+1-Kotherwise. Let I be the sequenceinfodict['iord']and let E be the sequenceinfodict['elist']. ThenE[I[1]], ..., E[I[L]]forms a decreasing sequence.
If the input argument points is provided (i.e., it is not None), the following additional outputs are placed in the output dictionary. Assume the points sequence is of length P.
'pts'
A rank-1 array of length P+2 containing the integration limits and the break points of the intervals in ascending order. This is an array giving the subintervals over which integration will occur.
'level'
A rank-1 integer array of length M (=limit), containing the subdivision levels of the subintervals, i.e., if (aa,bb) is a subinterval of
(pts[1], pts[2])wherepts[0]andpts[2]are adjacent elements ofinfodict['pts'], then (aa,bb) has level l if|bb-aa| = |pts[2]-pts[1]| * 2**(-l).'ndin'
A rank-1 integer array of length P+2. After the first integration over the intervals (pts[1], pts[2]), the error estimates over some of the intervals may have been increased artificially in order to put their subdivision forward. This array has ones in slots corresponding to the subintervals for which this happens.
Weighting the integrand
The input variables, weight and wvar, are used to weight the integrand by a select list of functions. Different integration methods are used to compute the integral with these weighting functions, and these do not support specifying break points. The possible values of weight and the corresponding weighting functions are.
========== =================================== ===================== ``weight`` Weight function used ``wvar`` ========== =================================== ===================== 'cos' cos(w*x) wvar = w 'sin' sin(w*x) wvar = w 'alg' g(x) = ((x-a)**alpha)*((b-x)**beta) wvar = (alpha, beta) 'alg-loga' g(x)*log(x-a) wvar = (alpha, beta) 'alg-logb' g(x)*log(b-x) wvar = (alpha, beta) 'alg-log' g(x)*log(x-a)*log(b-x) wvar = (alpha, beta) 'cauchy' 1/(x-c) wvar = c ========== =================================== =====================
wvar holds the parameter w, (alpha, beta), or c depending on the weight selected. In these expressions, a and b are the integration limits.
For the 'cos' and 'sin' weighting, additional inputs and outputs are available.
For weighted integrals with finite integration limits, the integration is performed using a Clenshaw-Curtis method, which uses Chebyshev moments. For repeated calculations, these moments are saved in the output dictionary:
'momcom'
The maximum level of Chebyshev moments that have been computed, i.e., if
M_cisinfodict['momcom']then the moments have been computed for intervals of length|b-a| * 2**(-l),l=0,1,...,M_c.'nnlog'
A rank-1 integer array of length M(=limit), containing the subdivision levels of the subintervals, i.e., an element of this array is equal to l if the corresponding subinterval is
|b-a|* 2**(-l).'chebmo'
A rank-2 array of shape (25, maxp1) containing the computed Chebyshev moments. These can be passed on to an integration over the same interval by passing this array as the second element of the sequence wopts and passing infodict['momcom'] as the first element.
If one of the integration limits is infinite, then a Fourier integral is computed (assuming w neq 0). If full_output is 1 and a numerical error is encountered, besides the error message attached to the output tuple, a dictionary is also appended to the output tuple which translates the error codes in the array info['ierlst'] to English messages. The output information dictionary contains the following entries instead of 'last', 'alist', 'blist', 'rlist', and 'elist':
'lst'
The number of subintervals needed for the integration (call it
K_f).'rslst'
A rank-1 array of length M_f=limlst, whose first
K_felements contain the integral contribution over the interval(a+(k-1)c, a+kc)wherec = (2*floor(|w|) + 1) * pi / |w|andk=1,2,...,K_f.'erlst'
A rank-1 array of length
M_fcontaining the error estimate corresponding to the interval in the same position ininfodict['rslist'].'ierlst'
A rank-1 integer array of length
M_fcontaining an error flag corresponding to the interval in the same position ininfodict['rslist']. See the explanation dictionary (last entry in the output tuple) for the meaning of the codes.
Details of QUADPACK level routines
quad 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 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 .
Integration of Complex Function of a Real Variable
A complex valued function, , of a real variable can be written as . Similarly, the integral of can be written as
assuming that the integrals of and exist over the interval [2]. Therefore, quad integrates complex-valued functions by integrating the real and imaginary components separately.
Array API Standard Support
quad 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^4_0 x^2 dx` and compare with an analytic resultfrom scipy import integrate import numpy as np x2 = lambda x: x**2✓
integrate.quad(x2, 0, 4) print(4**3 / 3.) # analytical result✗
invexp = lambda x: np.exp(-x)
✓integrate.quad(invexp, 0, np.inf)
✗f = lambda x, a: a*x y, err = integrate.quad(f, 0, 1, args=(1,)) y y, err = integrate.quad(f, 0, 1, args=(3,)) y✓
y = lambda x: 1 if x<=0 else 0 integrate.quad(y, -1, 1) integrate.quad(y, -1, 100) integrate.quad(y, -1, 10000)✓
See also
- dblquad
double integral
- fixed_quad
fixed-order Gaussian quadrature
- nquad
n-dimensional integrals (uses
quadrecursively)- romb
integrator for sampled data
- scipy.special
for coefficients and roots of orthogonal polynomials
- simpson
integrator for sampled data
- tplquad
triple integral
Aliases
-
scipy.integrate.quad