scipy.integrate.simpson#

scipy.integrate.simpson(y, *, x=None, dx=1.0, axis=-1, even=<object object>)[source]#

Integrate y(x) using samples along the given axis and the composite Simpson’s rule. If x is None, spacing of dx is assumed.

If there are an even number of samples, N, then there are an odd number of intervals (N-1), but Simpson’s rule requires an even number of intervals. The parameter ‘even’ controls how this is handled.

Parameters:
yarray_like

Array to be integrated.

xarray_like, optional

If given, the points at which y is sampled.

dxfloat, optional

Spacing of integration points along axis of x. Only used when x is None. Default is 1.

axisint, optional

Axis along which to integrate. Default is the last axis.

even{None, ‘simpson’, ‘avg’, ‘first’, ‘last’}, optional
‘avg’Average two results:
  1. use the first N-2 intervals with a trapezoidal rule on the last interval and

  2. use the last N-2 intervals with a trapezoidal rule on the first interval.

‘first’Use Simpson’s rule for the first N-2 intervals with

a trapezoidal rule on the last interval.

‘last’Use Simpson’s rule for the last N-2 intervals with a

trapezoidal rule on the first interval.

None : equivalent to ‘simpson’ (default)

‘simpson’Use Simpson’s rule for the first N-2 intervals with the

addition of a 3-point parabolic segment for the last interval using equations outlined by Cartwright [1]. If the axis to be integrated over only has two points then the integration falls back to a trapezoidal integration.

Added in version 1.11.0.

Changed in version 1.11.0: The newly added ‘simpson’ option is now the default as it is more accurate in most situations.

Deprecated since version 1.11.0: Parameter even is deprecated and will be removed in SciPy 1.14.0. After this time the behaviour for an even number of points will follow that of even=’simpson’.

Returns:
float

The estimated integral computed with the composite Simpson’s rule.

See also

quad

adaptive quadrature using QUADPACK

romberg

adaptive Romberg quadrature

quadrature

adaptive Gaussian quadrature

fixed_quad

fixed-order Gaussian quadrature

dblquad

double integrals

tplquad

triple integrals

romb

integrators for sampled data

cumulative_trapezoid

cumulative integration for sampled data

cumulative_simpson

cumulative integration using Simpson’s 1/3 rule

ode

ODE integrators

odeint

ODE integrators

Notes

For an odd number of samples that are equally spaced the result is exact if the function is a polynomial of order 3 or less. If the samples are not equally spaced, then the result is exact only if the function is a polynomial of order 2 or less.

References

[1]

Cartwright, Kenneth V. Simpson’s Rule Cumulative Integration with MS Excel and Irregularly-spaced Data. Journal of Mathematical Sciences and Mathematics Education. 12 (2): 1-9

Examples

>>> from scipy import integrate
>>> import numpy as np
>>> x = np.arange(0, 10)
>>> y = np.arange(0, 10)
>>> integrate.simpson(y, x=x)
40.5
>>> y = np.power(x, 3)
>>> integrate.simpson(y, x=x)
1640.5
>>> integrate.quad(lambda x: x**3, 0, 9)[0]
1640.25
>>> integrate.simpson(y, x=x, even='first')
1644.5