quadrature.py 30 KB

123456789101112131415161718192021222324252627282930313233343536373839404142434445464748495051525354555657585960616263646566676869707172737475767778798081828384858687888990919293949596979899100101102103104105106107108109110111112113114115116117118119120121122123124125126127128129130131132133134135136137138139140141142143144145146147148149150151152153154155156157158159160161162163164165166167168169170171172173174175176177178179180181182183184185186187188189190191192193194195196197198199200201202203204205206207208209210211212213214215216217218219220221222223224225226227228229230231232233234235236237238239240241242243244245246247248249250251252253254255256257258259260261262263264265266267268269270271272273274275276277278279280281282283284285286287288289290291292293294295296297298299300301302303304305306307308309310311312313314315316317318319320321322323324325326327328329330331332333334335336337338339340341342343344345346347348349350351352353354355356357358359360361362363364365366367368369370371372373374375376377378379380381382383384385386387388389390391392393394395396397398399400401402403404405406407408409410411412413414415416417418419420421422423424425426427428429430431432433434435436437438439440441442443444445446447448449450451452453454455456457458459460461462463464465466467468469470471472473474475476477478479480481482483484485486487488489490491492493494495496497498499500501502503504505506507508509510511512513514515516517518519520521522523524525526527528529530531532533534535536537538539540541542543544545546547548549550551552553554555556557558559560561562563564565566567568569570571572573574575576577578579580581582583584585586587588589590591592593594595596597598599600601602603604605606607608609610611612613614615616617618619620621622623624625626627628629630631632633634635636637638639640641642643644645646647648649650651652653654655656657658659660661662663664665666667668669670671672673674675676677678679680681682683684685686687688689690691692693694695696697698699700701702703704705706707708709710711712713714715716717718719720721722723724725726727728729730731732733734735736737738739740741742743744745746747748749750751752753754755756757758759760761762763764765766767768769770771772773774775776777778779780781782783784785786787788789790791792793794795796797798799800801802803804805806807808809810811812813814815816817818819820821822823824825826827828829830831832833834835836837838839840841842843844845846847848849850851852853854855856857858859860861862863864865866867868869870871872873874875876877878879880881882883884885886887888889890891892893894895896897898899900901902903904905906907908909910911912913914915916917918919920921922923924925926927928929930931932933934935936937938939940941942943944945946947948949950951952953954955956957
  1. from __future__ import division, print_function, absolute_import
  2. import numpy as np
  3. import math
  4. import warnings
  5. # trapz is a public function for scipy.integrate,
  6. # even though it's actually a numpy function.
  7. from numpy import trapz
  8. from scipy.special import roots_legendre
  9. from scipy.special import gammaln
  10. from scipy._lib.six import xrange
  11. __all__ = ['fixed_quad', 'quadrature', 'romberg', 'trapz', 'simps', 'romb',
  12. 'cumtrapz', 'newton_cotes']
  13. class AccuracyWarning(Warning):
  14. pass
  15. def _cached_roots_legendre(n):
  16. """
  17. Cache roots_legendre results to speed up calls of the fixed_quad
  18. function.
  19. """
  20. if n in _cached_roots_legendre.cache:
  21. return _cached_roots_legendre.cache[n]
  22. _cached_roots_legendre.cache[n] = roots_legendre(n)
  23. return _cached_roots_legendre.cache[n]
  24. _cached_roots_legendre.cache = dict()
  25. def fixed_quad(func, a, b, args=(), n=5):
  26. """
  27. Compute a definite integral using fixed-order Gaussian quadrature.
  28. Integrate `func` from `a` to `b` using Gaussian quadrature of
  29. order `n`.
  30. Parameters
  31. ----------
  32. func : callable
  33. A Python function or method to integrate (must accept vector inputs).
  34. If integrating a vector-valued function, the returned array must have
  35. shape ``(..., len(x))``.
  36. a : float
  37. Lower limit of integration.
  38. b : float
  39. Upper limit of integration.
  40. args : tuple, optional
  41. Extra arguments to pass to function, if any.
  42. n : int, optional
  43. Order of quadrature integration. Default is 5.
  44. Returns
  45. -------
  46. val : float
  47. Gaussian quadrature approximation to the integral
  48. none : None
  49. Statically returned value of None
  50. See Also
  51. --------
  52. quad : adaptive quadrature using QUADPACK
  53. dblquad : double integrals
  54. tplquad : triple integrals
  55. romberg : adaptive Romberg quadrature
  56. quadrature : adaptive Gaussian quadrature
  57. romb : integrators for sampled data
  58. simps : integrators for sampled data
  59. cumtrapz : cumulative integration for sampled data
  60. ode : ODE integrator
  61. odeint : ODE integrator
  62. Examples
  63. --------
  64. >>> from scipy import integrate
  65. >>> f = lambda x: x**8
  66. >>> integrate.fixed_quad(f, 0.0, 1.0, n=4)
  67. (0.1110884353741496, None)
  68. >>> integrate.fixed_quad(f, 0.0, 1.0, n=5)
  69. (0.11111111111111102, None)
  70. >>> print(1/9.0) # analytical result
  71. 0.1111111111111111
  72. >>> integrate.fixed_quad(np.cos, 0.0, np.pi/2, n=4)
  73. (0.9999999771971152, None)
  74. >>> integrate.fixed_quad(np.cos, 0.0, np.pi/2, n=5)
  75. (1.000000000039565, None)
  76. >>> np.sin(np.pi/2)-np.sin(0) # analytical result
  77. 1.0
  78. """
  79. x, w = _cached_roots_legendre(n)
  80. x = np.real(x)
  81. if np.isinf(a) or np.isinf(b):
  82. raise ValueError("Gaussian quadrature is only available for "
  83. "finite limits.")
  84. y = (b-a)*(x+1)/2.0 + a
  85. return (b-a)/2.0 * np.sum(w*func(y, *args), axis=-1), None
  86. def vectorize1(func, args=(), vec_func=False):
  87. """Vectorize the call to a function.
  88. This is an internal utility function used by `romberg` and
  89. `quadrature` to create a vectorized version of a function.
  90. If `vec_func` is True, the function `func` is assumed to take vector
  91. arguments.
  92. Parameters
  93. ----------
  94. func : callable
  95. User defined function.
  96. args : tuple, optional
  97. Extra arguments for the function.
  98. vec_func : bool, optional
  99. True if the function func takes vector arguments.
  100. Returns
  101. -------
  102. vfunc : callable
  103. A function that will take a vector argument and return the
  104. result.
  105. """
  106. if vec_func:
  107. def vfunc(x):
  108. return func(x, *args)
  109. else:
  110. def vfunc(x):
  111. if np.isscalar(x):
  112. return func(x, *args)
  113. x = np.asarray(x)
  114. # call with first point to get output type
  115. y0 = func(x[0], *args)
  116. n = len(x)
  117. dtype = getattr(y0, 'dtype', type(y0))
  118. output = np.empty((n,), dtype=dtype)
  119. output[0] = y0
  120. for i in xrange(1, n):
  121. output[i] = func(x[i], *args)
  122. return output
  123. return vfunc
  124. def quadrature(func, a, b, args=(), tol=1.49e-8, rtol=1.49e-8, maxiter=50,
  125. vec_func=True, miniter=1):
  126. """
  127. Compute a definite integral using fixed-tolerance Gaussian quadrature.
  128. Integrate `func` from `a` to `b` using Gaussian quadrature
  129. with absolute tolerance `tol`.
  130. Parameters
  131. ----------
  132. func : function
  133. A Python function or method to integrate.
  134. a : float
  135. Lower limit of integration.
  136. b : float
  137. Upper limit of integration.
  138. args : tuple, optional
  139. Extra arguments to pass to function.
  140. tol, rtol : float, optional
  141. Iteration stops when error between last two iterates is less than
  142. `tol` OR the relative change is less than `rtol`.
  143. maxiter : int, optional
  144. Maximum order of Gaussian quadrature.
  145. vec_func : bool, optional
  146. True or False if func handles arrays as arguments (is
  147. a "vector" function). Default is True.
  148. miniter : int, optional
  149. Minimum order of Gaussian quadrature.
  150. Returns
  151. -------
  152. val : float
  153. Gaussian quadrature approximation (within tolerance) to integral.
  154. err : float
  155. Difference between last two estimates of the integral.
  156. See also
  157. --------
  158. romberg: adaptive Romberg quadrature
  159. fixed_quad: fixed-order Gaussian quadrature
  160. quad: adaptive quadrature using QUADPACK
  161. dblquad: double integrals
  162. tplquad: triple integrals
  163. romb: integrator for sampled data
  164. simps: integrator for sampled data
  165. cumtrapz: cumulative integration for sampled data
  166. ode: ODE integrator
  167. odeint: ODE integrator
  168. Examples
  169. --------
  170. >>> from scipy import integrate
  171. >>> f = lambda x: x**8
  172. >>> integrate.quadrature(f, 0.0, 1.0)
  173. (0.11111111111111106, 4.163336342344337e-17)
  174. >>> print(1/9.0) # analytical result
  175. 0.1111111111111111
  176. >>> integrate.quadrature(np.cos, 0.0, np.pi/2)
  177. (0.9999999999999536, 3.9611425250996035e-11)
  178. >>> np.sin(np.pi/2)-np.sin(0) # analytical result
  179. 1.0
  180. """
  181. if not isinstance(args, tuple):
  182. args = (args,)
  183. vfunc = vectorize1(func, args, vec_func=vec_func)
  184. val = np.inf
  185. err = np.inf
  186. maxiter = max(miniter+1, maxiter)
  187. for n in xrange(miniter, maxiter+1):
  188. newval = fixed_quad(vfunc, a, b, (), n)[0]
  189. err = abs(newval-val)
  190. val = newval
  191. if err < tol or err < rtol*abs(val):
  192. break
  193. else:
  194. warnings.warn(
  195. "maxiter (%d) exceeded. Latest difference = %e" % (maxiter, err),
  196. AccuracyWarning)
  197. return val, err
  198. def tupleset(t, i, value):
  199. l = list(t)
  200. l[i] = value
  201. return tuple(l)
  202. def cumtrapz(y, x=None, dx=1.0, axis=-1, initial=None):
  203. """
  204. Cumulatively integrate y(x) using the composite trapezoidal rule.
  205. Parameters
  206. ----------
  207. y : array_like
  208. Values to integrate.
  209. x : array_like, optional
  210. The coordinate to integrate along. If None (default), use spacing `dx`
  211. between consecutive elements in `y`.
  212. dx : float, optional
  213. Spacing between elements of `y`. Only used if `x` is None.
  214. axis : int, optional
  215. Specifies the axis to cumulate. Default is -1 (last axis).
  216. initial : scalar, optional
  217. If given, insert this value at the beginning of the returned result.
  218. Typically this value should be 0. Default is None, which means no
  219. value at ``x[0]`` is returned and `res` has one element less than `y`
  220. along the axis of integration.
  221. Returns
  222. -------
  223. res : ndarray
  224. The result of cumulative integration of `y` along `axis`.
  225. If `initial` is None, the shape is such that the axis of integration
  226. has one less value than `y`. If `initial` is given, the shape is equal
  227. to that of `y`.
  228. See Also
  229. --------
  230. numpy.cumsum, numpy.cumprod
  231. quad: adaptive quadrature using QUADPACK
  232. romberg: adaptive Romberg quadrature
  233. quadrature: adaptive Gaussian quadrature
  234. fixed_quad: fixed-order Gaussian quadrature
  235. dblquad: double integrals
  236. tplquad: triple integrals
  237. romb: integrators for sampled data
  238. ode: ODE integrators
  239. odeint: ODE integrators
  240. Examples
  241. --------
  242. >>> from scipy import integrate
  243. >>> import matplotlib.pyplot as plt
  244. >>> x = np.linspace(-2, 2, num=20)
  245. >>> y = x
  246. >>> y_int = integrate.cumtrapz(y, x, initial=0)
  247. >>> plt.plot(x, y_int, 'ro', x, y[0] + 0.5 * x**2, 'b-')
  248. >>> plt.show()
  249. """
  250. y = np.asarray(y)
  251. if x is None:
  252. d = dx
  253. else:
  254. x = np.asarray(x)
  255. if x.ndim == 1:
  256. d = np.diff(x)
  257. # reshape to correct shape
  258. shape = [1] * y.ndim
  259. shape[axis] = -1
  260. d = d.reshape(shape)
  261. elif len(x.shape) != len(y.shape):
  262. raise ValueError("If given, shape of x must be 1-d or the "
  263. "same as y.")
  264. else:
  265. d = np.diff(x, axis=axis)
  266. if d.shape[axis] != y.shape[axis] - 1:
  267. raise ValueError("If given, length of x along axis must be the "
  268. "same as y.")
  269. nd = len(y.shape)
  270. slice1 = tupleset((slice(None),)*nd, axis, slice(1, None))
  271. slice2 = tupleset((slice(None),)*nd, axis, slice(None, -1))
  272. res = np.cumsum(d * (y[slice1] + y[slice2]) / 2.0, axis=axis)
  273. if initial is not None:
  274. if not np.isscalar(initial):
  275. raise ValueError("`initial` parameter should be a scalar.")
  276. shape = list(res.shape)
  277. shape[axis] = 1
  278. res = np.concatenate([np.full(shape, initial, dtype=res.dtype), res],
  279. axis=axis)
  280. return res
  281. def _basic_simps(y, start, stop, x, dx, axis):
  282. nd = len(y.shape)
  283. if start is None:
  284. start = 0
  285. step = 2
  286. slice_all = (slice(None),)*nd
  287. slice0 = tupleset(slice_all, axis, slice(start, stop, step))
  288. slice1 = tupleset(slice_all, axis, slice(start+1, stop+1, step))
  289. slice2 = tupleset(slice_all, axis, slice(start+2, stop+2, step))
  290. if x is None: # Even spaced Simpson's rule.
  291. result = np.sum(dx/3.0 * (y[slice0]+4*y[slice1]+y[slice2]),
  292. axis=axis)
  293. else:
  294. # Account for possibly different spacings.
  295. # Simpson's rule changes a bit.
  296. h = np.diff(x, axis=axis)
  297. sl0 = tupleset(slice_all, axis, slice(start, stop, step))
  298. sl1 = tupleset(slice_all, axis, slice(start+1, stop+1, step))
  299. h0 = h[sl0]
  300. h1 = h[sl1]
  301. hsum = h0 + h1
  302. hprod = h0 * h1
  303. h0divh1 = h0 / h1
  304. tmp = hsum/6.0 * (y[slice0]*(2-1.0/h0divh1) +
  305. y[slice1]*hsum*hsum/hprod +
  306. y[slice2]*(2-h0divh1))
  307. result = np.sum(tmp, axis=axis)
  308. return result
  309. def simps(y, x=None, dx=1, axis=-1, even='avg'):
  310. """
  311. Integrate y(x) using samples along the given axis and the composite
  312. Simpson's rule. If x is None, spacing of dx is assumed.
  313. If there are an even number of samples, N, then there are an odd
  314. number of intervals (N-1), but Simpson's rule requires an even number
  315. of intervals. The parameter 'even' controls how this is handled.
  316. Parameters
  317. ----------
  318. y : array_like
  319. Array to be integrated.
  320. x : array_like, optional
  321. If given, the points at which `y` is sampled.
  322. dx : int, optional
  323. Spacing of integration points along axis of `y`. Only used when
  324. `x` is None. Default is 1.
  325. axis : int, optional
  326. Axis along which to integrate. Default is the last axis.
  327. even : str {'avg', 'first', 'last'}, optional
  328. 'avg' : Average two results:1) use the first N-2 intervals with
  329. a trapezoidal rule on the last interval and 2) use the last
  330. N-2 intervals with a trapezoidal rule on the first interval.
  331. 'first' : Use Simpson's rule for the first N-2 intervals with
  332. a trapezoidal rule on the last interval.
  333. 'last' : Use Simpson's rule for the last N-2 intervals with a
  334. trapezoidal rule on the first interval.
  335. See Also
  336. --------
  337. quad: adaptive quadrature using QUADPACK
  338. romberg: adaptive Romberg quadrature
  339. quadrature: adaptive Gaussian quadrature
  340. fixed_quad: fixed-order Gaussian quadrature
  341. dblquad: double integrals
  342. tplquad: triple integrals
  343. romb: integrators for sampled data
  344. cumtrapz: cumulative integration for sampled data
  345. ode: ODE integrators
  346. odeint: ODE integrators
  347. Notes
  348. -----
  349. For an odd number of samples that are equally spaced the result is
  350. exact if the function is a polynomial of order 3 or less. If
  351. the samples are not equally spaced, then the result is exact only
  352. if the function is a polynomial of order 2 or less.
  353. Examples
  354. --------
  355. >>> from scipy import integrate
  356. >>> x = np.arange(0, 10)
  357. >>> y = np.arange(0, 10)
  358. >>> integrate.simps(y, x)
  359. 40.5
  360. >>> y = np.power(x, 3)
  361. >>> integrate.simps(y, x)
  362. 1642.5
  363. >>> integrate.quad(lambda x: x**3, 0, 9)[0]
  364. 1640.25
  365. >>> integrate.simps(y, x, even='first')
  366. 1644.5
  367. """
  368. y = np.asarray(y)
  369. nd = len(y.shape)
  370. N = y.shape[axis]
  371. last_dx = dx
  372. first_dx = dx
  373. returnshape = 0
  374. if x is not None:
  375. x = np.asarray(x)
  376. if len(x.shape) == 1:
  377. shapex = [1] * nd
  378. shapex[axis] = x.shape[0]
  379. saveshape = x.shape
  380. returnshape = 1
  381. x = x.reshape(tuple(shapex))
  382. elif len(x.shape) != len(y.shape):
  383. raise ValueError("If given, shape of x must be 1-d or the "
  384. "same as y.")
  385. if x.shape[axis] != N:
  386. raise ValueError("If given, length of x along axis must be the "
  387. "same as y.")
  388. if N % 2 == 0:
  389. val = 0.0
  390. result = 0.0
  391. slice1 = (slice(None),)*nd
  392. slice2 = (slice(None),)*nd
  393. if even not in ['avg', 'last', 'first']:
  394. raise ValueError("Parameter 'even' must be "
  395. "'avg', 'last', or 'first'.")
  396. # Compute using Simpson's rule on first intervals
  397. if even in ['avg', 'first']:
  398. slice1 = tupleset(slice1, axis, -1)
  399. slice2 = tupleset(slice2, axis, -2)
  400. if x is not None:
  401. last_dx = x[slice1] - x[slice2]
  402. val += 0.5*last_dx*(y[slice1]+y[slice2])
  403. result = _basic_simps(y, 0, N-3, x, dx, axis)
  404. # Compute using Simpson's rule on last set of intervals
  405. if even in ['avg', 'last']:
  406. slice1 = tupleset(slice1, axis, 0)
  407. slice2 = tupleset(slice2, axis, 1)
  408. if x is not None:
  409. first_dx = x[tuple(slice2)] - x[tuple(slice1)]
  410. val += 0.5*first_dx*(y[slice2]+y[slice1])
  411. result += _basic_simps(y, 1, N-2, x, dx, axis)
  412. if even == 'avg':
  413. val /= 2.0
  414. result /= 2.0
  415. result = result + val
  416. else:
  417. result = _basic_simps(y, 0, N-2, x, dx, axis)
  418. if returnshape:
  419. x = x.reshape(saveshape)
  420. return result
  421. def romb(y, dx=1.0, axis=-1, show=False):
  422. """
  423. Romberg integration using samples of a function.
  424. Parameters
  425. ----------
  426. y : array_like
  427. A vector of ``2**k + 1`` equally-spaced samples of a function.
  428. dx : float, optional
  429. The sample spacing. Default is 1.
  430. axis : int, optional
  431. The axis along which to integrate. Default is -1 (last axis).
  432. show : bool, optional
  433. When `y` is a single 1-D array, then if this argument is True
  434. print the table showing Richardson extrapolation from the
  435. samples. Default is False.
  436. Returns
  437. -------
  438. romb : ndarray
  439. The integrated result for `axis`.
  440. See also
  441. --------
  442. quad : adaptive quadrature using QUADPACK
  443. romberg : adaptive Romberg quadrature
  444. quadrature : adaptive Gaussian quadrature
  445. fixed_quad : fixed-order Gaussian quadrature
  446. dblquad : double integrals
  447. tplquad : triple integrals
  448. simps : integrators for sampled data
  449. cumtrapz : cumulative integration for sampled data
  450. ode : ODE integrators
  451. odeint : ODE integrators
  452. Examples
  453. --------
  454. >>> from scipy import integrate
  455. >>> x = np.arange(10, 14.25, 0.25)
  456. >>> y = np.arange(3, 12)
  457. >>> integrate.romb(y)
  458. 56.0
  459. >>> y = np.sin(np.power(x, 2.5))
  460. >>> integrate.romb(y)
  461. -0.742561336672229
  462. >>> integrate.romb(y, show=True)
  463. Richardson Extrapolation Table for Romberg Integration
  464. ====================================================================
  465. -0.81576
  466. 4.63862 6.45674
  467. -1.10581 -3.02062 -3.65245
  468. -2.57379 -3.06311 -3.06595 -3.05664
  469. -1.34093 -0.92997 -0.78776 -0.75160 -0.74256
  470. ====================================================================
  471. -0.742561336672229
  472. """
  473. y = np.asarray(y)
  474. nd = len(y.shape)
  475. Nsamps = y.shape[axis]
  476. Ninterv = Nsamps-1
  477. n = 1
  478. k = 0
  479. while n < Ninterv:
  480. n <<= 1
  481. k += 1
  482. if n != Ninterv:
  483. raise ValueError("Number of samples must be one plus a "
  484. "non-negative power of 2.")
  485. R = {}
  486. slice_all = (slice(None),) * nd
  487. slice0 = tupleset(slice_all, axis, 0)
  488. slicem1 = tupleset(slice_all, axis, -1)
  489. h = Ninterv * np.asarray(dx, dtype=float)
  490. R[(0, 0)] = (y[slice0] + y[slicem1])/2.0*h
  491. slice_R = slice_all
  492. start = stop = step = Ninterv
  493. for i in xrange(1, k+1):
  494. start >>= 1
  495. slice_R = tupleset(slice_R, axis, slice(start, stop, step))
  496. step >>= 1
  497. R[(i, 0)] = 0.5*(R[(i-1, 0)] + h*y[slice_R].sum(axis=axis))
  498. for j in xrange(1, i+1):
  499. prev = R[(i, j-1)]
  500. R[(i, j)] = prev + (prev-R[(i-1, j-1)]) / ((1 << (2*j))-1)
  501. h /= 2.0
  502. if show:
  503. if not np.isscalar(R[(0, 0)]):
  504. print("*** Printing table only supported for integrals" +
  505. " of a single data set.")
  506. else:
  507. try:
  508. precis = show[0]
  509. except (TypeError, IndexError):
  510. precis = 5
  511. try:
  512. width = show[1]
  513. except (TypeError, IndexError):
  514. width = 8
  515. formstr = "%%%d.%df" % (width, precis)
  516. title = "Richardson Extrapolation Table for Romberg Integration"
  517. print("", title.center(68), "=" * 68, sep="\n", end="\n")
  518. for i in xrange(k+1):
  519. for j in xrange(i+1):
  520. print(formstr % R[(i, j)], end=" ")
  521. print()
  522. print("=" * 68)
  523. print()
  524. return R[(k, k)]
  525. # Romberg quadratures for numeric integration.
  526. #
  527. # Written by Scott M. Ransom <ransom@cfa.harvard.edu>
  528. # last revision: 14 Nov 98
  529. #
  530. # Cosmetic changes by Konrad Hinsen <hinsen@cnrs-orleans.fr>
  531. # last revision: 1999-7-21
  532. #
  533. # Adapted to scipy by Travis Oliphant <oliphant.travis@ieee.org>
  534. # last revision: Dec 2001
  535. def _difftrap(function, interval, numtraps):
  536. """
  537. Perform part of the trapezoidal rule to integrate a function.
  538. Assume that we had called difftrap with all lower powers-of-2
  539. starting with 1. Calling difftrap only returns the summation
  540. of the new ordinates. It does _not_ multiply by the width
  541. of the trapezoids. This must be performed by the caller.
  542. 'function' is the function to evaluate (must accept vector arguments).
  543. 'interval' is a sequence with lower and upper limits
  544. of integration.
  545. 'numtraps' is the number of trapezoids to use (must be a
  546. power-of-2).
  547. """
  548. if numtraps <= 0:
  549. raise ValueError("numtraps must be > 0 in difftrap().")
  550. elif numtraps == 1:
  551. return 0.5*(function(interval[0])+function(interval[1]))
  552. else:
  553. numtosum = numtraps/2
  554. h = float(interval[1]-interval[0])/numtosum
  555. lox = interval[0] + 0.5 * h
  556. points = lox + h * np.arange(numtosum)
  557. s = np.sum(function(points), axis=0)
  558. return s
  559. def _romberg_diff(b, c, k):
  560. """
  561. Compute the differences for the Romberg quadrature corrections.
  562. See Forman Acton's "Real Computing Made Real," p 143.
  563. """
  564. tmp = 4.0**k
  565. return (tmp * c - b)/(tmp - 1.0)
  566. def _printresmat(function, interval, resmat):
  567. # Print the Romberg result matrix.
  568. i = j = 0
  569. print('Romberg integration of', repr(function), end=' ')
  570. print('from', interval)
  571. print('')
  572. print('%6s %9s %9s' % ('Steps', 'StepSize', 'Results'))
  573. for i in xrange(len(resmat)):
  574. print('%6d %9f' % (2**i, (interval[1]-interval[0])/(2.**i)), end=' ')
  575. for j in xrange(i+1):
  576. print('%9f' % (resmat[i][j]), end=' ')
  577. print('')
  578. print('')
  579. print('The final result is', resmat[i][j], end=' ')
  580. print('after', 2**(len(resmat)-1)+1, 'function evaluations.')
  581. def romberg(function, a, b, args=(), tol=1.48e-8, rtol=1.48e-8, show=False,
  582. divmax=10, vec_func=False):
  583. """
  584. Romberg integration of a callable function or method.
  585. Returns the integral of `function` (a function of one variable)
  586. over the interval (`a`, `b`).
  587. If `show` is 1, the triangular array of the intermediate results
  588. will be printed. If `vec_func` is True (default is False), then
  589. `function` is assumed to support vector arguments.
  590. Parameters
  591. ----------
  592. function : callable
  593. Function to be integrated.
  594. a : float
  595. Lower limit of integration.
  596. b : float
  597. Upper limit of integration.
  598. Returns
  599. -------
  600. results : float
  601. Result of the integration.
  602. Other Parameters
  603. ----------------
  604. args : tuple, optional
  605. Extra arguments to pass to function. Each element of `args` will
  606. be passed as a single argument to `func`. Default is to pass no
  607. extra arguments.
  608. tol, rtol : float, optional
  609. The desired absolute and relative tolerances. Defaults are 1.48e-8.
  610. show : bool, optional
  611. Whether to print the results. Default is False.
  612. divmax : int, optional
  613. Maximum order of extrapolation. Default is 10.
  614. vec_func : bool, optional
  615. Whether `func` handles arrays as arguments (i.e whether it is a
  616. "vector" function). Default is False.
  617. See Also
  618. --------
  619. fixed_quad : Fixed-order Gaussian quadrature.
  620. quad : Adaptive quadrature using QUADPACK.
  621. dblquad : Double integrals.
  622. tplquad : Triple integrals.
  623. romb : Integrators for sampled data.
  624. simps : Integrators for sampled data.
  625. cumtrapz : Cumulative integration for sampled data.
  626. ode : ODE integrator.
  627. odeint : ODE integrator.
  628. References
  629. ----------
  630. .. [1] 'Romberg's method' https://en.wikipedia.org/wiki/Romberg%27s_method
  631. Examples
  632. --------
  633. Integrate a gaussian from 0 to 1 and compare to the error function.
  634. >>> from scipy import integrate
  635. >>> from scipy.special import erf
  636. >>> gaussian = lambda x: 1/np.sqrt(np.pi) * np.exp(-x**2)
  637. >>> result = integrate.romberg(gaussian, 0, 1, show=True)
  638. Romberg integration of <function vfunc at ...> from [0, 1]
  639. ::
  640. Steps StepSize Results
  641. 1 1.000000 0.385872
  642. 2 0.500000 0.412631 0.421551
  643. 4 0.250000 0.419184 0.421368 0.421356
  644. 8 0.125000 0.420810 0.421352 0.421350 0.421350
  645. 16 0.062500 0.421215 0.421350 0.421350 0.421350 0.421350
  646. 32 0.031250 0.421317 0.421350 0.421350 0.421350 0.421350 0.421350
  647. The final result is 0.421350396475 after 33 function evaluations.
  648. >>> print("%g %g" % (2*result, erf(1)))
  649. 0.842701 0.842701
  650. """
  651. if np.isinf(a) or np.isinf(b):
  652. raise ValueError("Romberg integration only available "
  653. "for finite limits.")
  654. vfunc = vectorize1(function, args, vec_func=vec_func)
  655. n = 1
  656. interval = [a, b]
  657. intrange = b - a
  658. ordsum = _difftrap(vfunc, interval, n)
  659. result = intrange * ordsum
  660. resmat = [[result]]
  661. err = np.inf
  662. last_row = resmat[0]
  663. for i in xrange(1, divmax+1):
  664. n *= 2
  665. ordsum += _difftrap(vfunc, interval, n)
  666. row = [intrange * ordsum / n]
  667. for k in xrange(i):
  668. row.append(_romberg_diff(last_row[k], row[k], k+1))
  669. result = row[i]
  670. lastresult = last_row[i-1]
  671. if show:
  672. resmat.append(row)
  673. err = abs(result - lastresult)
  674. if err < tol or err < rtol * abs(result):
  675. break
  676. last_row = row
  677. else:
  678. warnings.warn(
  679. "divmax (%d) exceeded. Latest difference = %e" % (divmax, err),
  680. AccuracyWarning)
  681. if show:
  682. _printresmat(vfunc, interval, resmat)
  683. return result
  684. # Coefficients for Netwon-Cotes quadrature
  685. #
  686. # These are the points being used
  687. # to construct the local interpolating polynomial
  688. # a are the weights for Newton-Cotes integration
  689. # B is the error coefficient.
  690. # error in these coefficients grows as N gets larger.
  691. # or as samples are closer and closer together
  692. # You can use maxima to find these rational coefficients
  693. # for equally spaced data using the commands
  694. # a(i,N) := integrate(product(r-j,j,0,i-1) * product(r-j,j,i+1,N),r,0,N) / ((N-i)! * i!) * (-1)^(N-i);
  695. # Be(N) := N^(N+2)/(N+2)! * (N/(N+3) - sum((i/N)^(N+2)*a(i,N),i,0,N));
  696. # Bo(N) := N^(N+1)/(N+1)! * (N/(N+2) - sum((i/N)^(N+1)*a(i,N),i,0,N));
  697. # B(N) := (if (mod(N,2)=0) then Be(N) else Bo(N));
  698. #
  699. # pre-computed for equally-spaced weights
  700. #
  701. # num_a, den_a, int_a, num_B, den_B = _builtincoeffs[N]
  702. #
  703. # a = num_a*array(int_a)/den_a
  704. # B = num_B*1.0 / den_B
  705. #
  706. # integrate(f(x),x,x_0,x_N) = dx*sum(a*f(x_i)) + B*(dx)^(2k+3) f^(2k+2)(x*)
  707. # where k = N // 2
  708. #
  709. _builtincoeffs = {
  710. 1: (1,2,[1,1],-1,12),
  711. 2: (1,3,[1,4,1],-1,90),
  712. 3: (3,8,[1,3,3,1],-3,80),
  713. 4: (2,45,[7,32,12,32,7],-8,945),
  714. 5: (5,288,[19,75,50,50,75,19],-275,12096),
  715. 6: (1,140,[41,216,27,272,27,216,41],-9,1400),
  716. 7: (7,17280,[751,3577,1323,2989,2989,1323,3577,751],-8183,518400),
  717. 8: (4,14175,[989,5888,-928,10496,-4540,10496,-928,5888,989],
  718. -2368,467775),
  719. 9: (9,89600,[2857,15741,1080,19344,5778,5778,19344,1080,
  720. 15741,2857], -4671, 394240),
  721. 10: (5,299376,[16067,106300,-48525,272400,-260550,427368,
  722. -260550,272400,-48525,106300,16067],
  723. -673175, 163459296),
  724. 11: (11,87091200,[2171465,13486539,-3237113, 25226685,-9595542,
  725. 15493566,15493566,-9595542,25226685,-3237113,
  726. 13486539,2171465], -2224234463, 237758976000),
  727. 12: (1, 5255250, [1364651,9903168,-7587864,35725120,-51491295,
  728. 87516288,-87797136,87516288,-51491295,35725120,
  729. -7587864,9903168,1364651], -3012, 875875),
  730. 13: (13, 402361344000,[8181904909, 56280729661, -31268252574,
  731. 156074417954,-151659573325,206683437987,
  732. -43111992612,-43111992612,206683437987,
  733. -151659573325,156074417954,-31268252574,
  734. 56280729661,8181904909], -2639651053,
  735. 344881152000),
  736. 14: (7, 2501928000, [90241897,710986864,-770720657,3501442784,
  737. -6625093363,12630121616,-16802270373,19534438464,
  738. -16802270373,12630121616,-6625093363,3501442784,
  739. -770720657,710986864,90241897], -3740727473,
  740. 1275983280000)
  741. }
  742. def newton_cotes(rn, equal=0):
  743. r"""
  744. Return weights and error coefficient for Newton-Cotes integration.
  745. Suppose we have (N+1) samples of f at the positions
  746. x_0, x_1, ..., x_N. Then an N-point Newton-Cotes formula for the
  747. integral between x_0 and x_N is:
  748. :math:`\int_{x_0}^{x_N} f(x)dx = \Delta x \sum_{i=0}^{N} a_i f(x_i)
  749. + B_N (\Delta x)^{N+2} f^{N+1} (\xi)`
  750. where :math:`\xi \in [x_0,x_N]`
  751. and :math:`\Delta x = \frac{x_N-x_0}{N}` is the average samples spacing.
  752. If the samples are equally-spaced and N is even, then the error
  753. term is :math:`B_N (\Delta x)^{N+3} f^{N+2}(\xi)`.
  754. Parameters
  755. ----------
  756. rn : int
  757. The integer order for equally-spaced data or the relative positions of
  758. the samples with the first sample at 0 and the last at N, where N+1 is
  759. the length of `rn`. N is the order of the Newton-Cotes integration.
  760. equal : int, optional
  761. Set to 1 to enforce equally spaced data.
  762. Returns
  763. -------
  764. an : ndarray
  765. 1-D array of weights to apply to the function at the provided sample
  766. positions.
  767. B : float
  768. Error coefficient.
  769. Examples
  770. --------
  771. Compute the integral of sin(x) in [0, :math:`\pi`]:
  772. >>> from scipy.integrate import newton_cotes
  773. >>> def f(x):
  774. ... return np.sin(x)
  775. >>> a = 0
  776. >>> b = np.pi
  777. >>> exact = 2
  778. >>> for N in [2, 4, 6, 8, 10]:
  779. ... x = np.linspace(a, b, N + 1)
  780. ... an, B = newton_cotes(N, 1)
  781. ... dx = (b - a) / N
  782. ... quad = dx * np.sum(an * f(x))
  783. ... error = abs(quad - exact)
  784. ... print('{:2d} {:10.9f} {:.5e}'.format(N, quad, error))
  785. ...
  786. 2 2.094395102 9.43951e-02
  787. 4 1.998570732 1.42927e-03
  788. 6 2.000017814 1.78136e-05
  789. 8 1.999999835 1.64725e-07
  790. 10 2.000000001 1.14677e-09
  791. Notes
  792. -----
  793. Normally, the Newton-Cotes rules are used on smaller integration
  794. regions and a composite rule is used to return the total integral.
  795. """
  796. try:
  797. N = len(rn)-1
  798. if equal:
  799. rn = np.arange(N+1)
  800. elif np.all(np.diff(rn) == 1):
  801. equal = 1
  802. except Exception:
  803. N = rn
  804. rn = np.arange(N+1)
  805. equal = 1
  806. if equal and N in _builtincoeffs:
  807. na, da, vi, nb, db = _builtincoeffs[N]
  808. an = na * np.array(vi, dtype=float) / da
  809. return an, float(nb)/db
  810. if (rn[0] != 0) or (rn[-1] != N):
  811. raise ValueError("The sample positions must start at 0"
  812. " and end at N")
  813. yi = rn / float(N)
  814. ti = 2 * yi - 1
  815. nvec = np.arange(N+1)
  816. C = ti ** nvec[:, np.newaxis]
  817. Cinv = np.linalg.inv(C)
  818. # improve precision of result
  819. for i in range(2):
  820. Cinv = 2*Cinv - Cinv.dot(C).dot(Cinv)
  821. vec = 2.0 / (nvec[::2]+1)
  822. ai = Cinv[:, ::2].dot(vec) * (N / 2.)
  823. if (N % 2 == 0) and equal:
  824. BN = N/(N+3.)
  825. power = N+2
  826. else:
  827. BN = N/(N+2.)
  828. power = N+1
  829. BN = BN - np.dot(yi**power, ai)
  830. p1 = power+1
  831. fac = power*math.log(N) - gammaln(p1)
  832. fac = math.exp(fac)
  833. return ai, BN*fac