scipy.integrate.nquad¶
-
scipy.integrate.nquad(func, ranges, args=None, opts=None, full_output=False)[source]¶ Integration over multiple variables.
Wraps
quadto 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. Where`t0, ... tm`are 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 f may be a
scipy.LowLevelCallablewith one of the signatures:double 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 thescipy.LowLevelCallable.- rangesiterable object
Each 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).- argsiterable object, optional
Additional arguments
t0, ..., tn, required by func, ranges, andopts.- optsiterable object or dict, optional
Options 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 over x0, 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
quadandquad_explain.- full_outputbool, optional
Partial implementation of
full_outputfrom scipy.integrate.quad. The number of integrand function evaluationsnevalcan be obtained by settingfull_output=Truewhen calling nquad.
- Returns
- resultfloat
The result of the integration.
- abserrfloat
The maximum of the estimates of the absolute error in the various integration results.
- out_dictdict, optional
A dict containing additional information on the integration.
See also
quad1-D numerical integration
dblquad,tplquaddouble and triple integrals
fixed_quadfixed-order Gaussian quadrature
quadratureadaptive Gaussian quadrature
Examples
>>> 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) (1.5267454070738633, 2.9437360001402324e-14, {'neval': 388962})
>>> scale = .1 >>> def func2(x0, x1, x2, x3, t0, t1): ... return x0*x1*x3**2 + np.sin(x2) + 1 + (1 if x0+t1*x1-t0>0 else 0) >>> def lim0(x1, x2, x3, t0, t1): ... return [scale * (x1**2 + x2 + np.cos(x3)*t0*t1 + 1) - 1, ... scale * (x1**2 + x2 + np.cos(x3)*t0*t1 + 1) + 1] >>> def lim1(x2, x3, t0, t1): ... return [scale * (t0*x2 + t1*x3) - 1, ... scale * (t0*x2 + t1*x3) + 1] >>> def lim2(x3, t0, t1): ... return [scale * (x3 + t0**2*t1**3) - 1, ... scale * (x3 + t0**2*t1**3) + 1] >>> def lim3(t0, t1): ... return [scale * (t0+t1) - 1, scale * (t0+t1) + 1] >>> def opts0(x1, x2, x3, t0, t1): ... return {'points' : [t0 - t1*x1]} >>> def opts1(x2, x3, t0, t1): ... return {} >>> def opts2(x3, t0, t1): ... return {} >>> def opts3(t0, t1): ... return {} >>> integrate.nquad(func2, [lim0, lim1, lim2, lim3], args=(0,0), ... opts=[opts0, opts1, opts2, opts3]) (25.066666666666666, 2.7829590483937256e-13)