quadrature.py 38 KB

123456789101112131415161718192021222324252627282930313233343536373839404142434445464748495051525354555657585960616263646566676869707172737475767778798081828384858687888990919293949596979899100101102103104105106107108109110111112113114115116117118119120121122123124125126127128129130131132133134135136137138139140141142143144145146147148149150151152153154155156157158159160161162163164165166167168169170171172173174175176177178179180181182183184185186187188189190191192193194195196197198199200201202203204205206207208209210211212213214215216217218219220221222223224225226227228229230231232233234235236237238239240241242243244245246247248249250251252253254255256257258259260261262263264265266267268269270271272273274275276277278279280281282283284285286287288289290291292293294295296297298299300301302303304305306307308309310311312313314315316317318319320321322323324325326327328329330331332333334335336337338339340341342343344345346347348349350351352353354355356357358359360361362363364365366367368369370371372373374375376377378379380381382383384385386387388389390391392393394395396397398399400401402403404405406407408409410411412413414415416417418419420421422423424425426427428429430431432433434435436437438439440441442443444445446447448449450451452453454455456457458459460461462463464465466467468469470471472473474475476477478479480481482483484485486487488489490491492493494495496497498499500501502503504505506507508509510511512513514515516517518519520521522523524525526527528529530531532533534535536537538539540541542543544545546547548549550551552553554555556557558559560561562563564565566567568569570571572573574575576577578579580581582583584585586587588589590591592593594595596597598599600601602603604605606607608609610611612613614615616617618619620621622623624625626627628629630631632633634635636637638639640641642643644645646647648649650651652653654655656657658659660661662663664665666667668669670671672673674675676677678679680681682683684685686687688689690691692693694695696697698699700701702703704705706707708709710711712713714715716717718719720721722723724725726727728729730731732733734735736737738739740741742743744745746747748749750751752753754755756757758759760761762763764765766767768769770771772773774775776777778779780781782783784785786787788789790791792793794795796797798799800801802803804805806807808809810811812813814815816817818819820821822823824825826827828829830831832833834835836837838839840841842843844845846847848849850851852853854855856857858859860861862863864865866867868869870871872873874875876877878879880881882883884885886887888889890891892893894895896897898899900901902903904905906907908909910911912913914915916917918919920921922923924925926927928929930931932933934935936937938939940941942943944945946947948949950951952953954955956957958959960961962963964965966967968969970971972973974975976977978979980981982983984985986987988989990991992993994995996997998999100010011002100310041005100610071008100910101011
  1. import math
  2. from ..libmp.backend import xrange
  3. class QuadratureRule(object):
  4. """
  5. Quadrature rules are implemented using this class, in order to
  6. simplify the code and provide a common infrastructure
  7. for tasks such as error estimation and node caching.
  8. You can implement a custom quadrature rule by subclassing
  9. :class:`QuadratureRule` and implementing the appropriate
  10. methods. The subclass can then be used by :func:`~mpmath.quad` by
  11. passing it as the *method* argument.
  12. :class:`QuadratureRule` instances are supposed to be singletons.
  13. :class:`QuadratureRule` therefore implements instance caching
  14. in :func:`~mpmath.__new__`.
  15. """
  16. def __init__(self, ctx):
  17. self.ctx = ctx
  18. self.standard_cache = {}
  19. self.transformed_cache = {}
  20. self.interval_count = {}
  21. def clear(self):
  22. """
  23. Delete cached node data.
  24. """
  25. self.standard_cache = {}
  26. self.transformed_cache = {}
  27. self.interval_count = {}
  28. def calc_nodes(self, degree, prec, verbose=False):
  29. r"""
  30. Compute nodes for the standard interval `[-1, 1]`. Subclasses
  31. should probably implement only this method, and use
  32. :func:`~mpmath.get_nodes` method to retrieve the nodes.
  33. """
  34. raise NotImplementedError
  35. def get_nodes(self, a, b, degree, prec, verbose=False):
  36. """
  37. Return nodes for given interval, degree and precision. The
  38. nodes are retrieved from a cache if already computed;
  39. otherwise they are computed by calling :func:`~mpmath.calc_nodes`
  40. and are then cached.
  41. Subclasses should probably not implement this method,
  42. but just implement :func:`~mpmath.calc_nodes` for the actual
  43. node computation.
  44. """
  45. key = (a, b, degree, prec)
  46. if key in self.transformed_cache:
  47. return self.transformed_cache[key]
  48. orig = self.ctx.prec
  49. try:
  50. self.ctx.prec = prec+20
  51. # Get nodes on standard interval
  52. if (degree, prec) in self.standard_cache:
  53. nodes = self.standard_cache[degree, prec]
  54. else:
  55. nodes = self.calc_nodes(degree, prec, verbose)
  56. self.standard_cache[degree, prec] = nodes
  57. # Transform to general interval
  58. nodes = self.transform_nodes(nodes, a, b, verbose)
  59. if key in self.interval_count:
  60. self.transformed_cache[key] = nodes
  61. else:
  62. self.interval_count[key] = True
  63. finally:
  64. self.ctx.prec = orig
  65. return nodes
  66. def transform_nodes(self, nodes, a, b, verbose=False):
  67. r"""
  68. Rescale standardized nodes (for `[-1, 1]`) to a general
  69. interval `[a, b]`. For a finite interval, a simple linear
  70. change of variables is used. Otherwise, the following
  71. transformations are used:
  72. .. math ::
  73. \lbrack a, \infty \rbrack : t = \frac{1}{x} + (a-1)
  74. \lbrack -\infty, b \rbrack : t = (b+1) - \frac{1}{x}
  75. \lbrack -\infty, \infty \rbrack : t = \frac{x}{\sqrt{1-x^2}}
  76. """
  77. ctx = self.ctx
  78. a = ctx.convert(a)
  79. b = ctx.convert(b)
  80. one = ctx.one
  81. if (a, b) == (-one, one):
  82. return nodes
  83. half = ctx.mpf(0.5)
  84. new_nodes = []
  85. if ctx.isinf(a) or ctx.isinf(b):
  86. if (a, b) == (ctx.ninf, ctx.inf):
  87. p05 = -half
  88. for x, w in nodes:
  89. x2 = x*x
  90. px1 = one-x2
  91. spx1 = px1**p05
  92. x = x*spx1
  93. w *= spx1/px1
  94. new_nodes.append((x, w))
  95. elif a == ctx.ninf:
  96. b1 = b+1
  97. for x, w in nodes:
  98. u = 2/(x+one)
  99. x = b1-u
  100. w *= half*u**2
  101. new_nodes.append((x, w))
  102. elif b == ctx.inf:
  103. a1 = a-1
  104. for x, w in nodes:
  105. u = 2/(x+one)
  106. x = a1+u
  107. w *= half*u**2
  108. new_nodes.append((x, w))
  109. elif a == ctx.inf or b == ctx.ninf:
  110. return [(x,-w) for (x,w) in self.transform_nodes(nodes, b, a, verbose)]
  111. else:
  112. raise NotImplementedError
  113. else:
  114. # Simple linear change of variables
  115. C = (b-a)/2
  116. D = (b+a)/2
  117. for x, w in nodes:
  118. new_nodes.append((D+C*x, C*w))
  119. return new_nodes
  120. def guess_degree(self, prec):
  121. """
  122. Given a desired precision `p` in bits, estimate the degree `m`
  123. of the quadrature required to accomplish full accuracy for
  124. typical integrals. By default, :func:`~mpmath.quad` will perform up
  125. to `m` iterations. The value of `m` should be a slight
  126. overestimate, so that "slightly bad" integrals can be dealt
  127. with automatically using a few extra iterations. On the
  128. other hand, it should not be too big, so :func:`~mpmath.quad` can
  129. quit within a reasonable amount of time when it is given
  130. an "unsolvable" integral.
  131. The default formula used by :func:`~mpmath.guess_degree` is tuned
  132. for both :class:`TanhSinh` and :class:`GaussLegendre`.
  133. The output is roughly as follows:
  134. +---------+---------+
  135. | `p` | `m` |
  136. +=========+=========+
  137. | 50 | 6 |
  138. +---------+---------+
  139. | 100 | 7 |
  140. +---------+---------+
  141. | 500 | 10 |
  142. +---------+---------+
  143. | 3000 | 12 |
  144. +---------+---------+
  145. This formula is based purely on a limited amount of
  146. experimentation and will sometimes be wrong.
  147. """
  148. # Expected degree
  149. # XXX: use mag
  150. g = int(4 + max(0, self.ctx.log(prec/30.0, 2)))
  151. # Reasonable "worst case"
  152. g += 2
  153. return g
  154. def estimate_error(self, results, prec, epsilon):
  155. r"""
  156. Given results from integrations `[I_1, I_2, \ldots, I_k]` done
  157. with a quadrature of rule of degree `1, 2, \ldots, k`, estimate
  158. the error of `I_k`.
  159. For `k = 2`, we estimate `|I_{\infty}-I_2|` as `|I_2-I_1|`.
  160. For `k > 2`, we extrapolate `|I_{\infty}-I_k| \approx |I_{k+1}-I_k|`
  161. from `|I_k-I_{k-1}|` and `|I_k-I_{k-2}|` under the assumption
  162. that each degree increment roughly doubles the accuracy of
  163. the quadrature rule (this is true for both :class:`TanhSinh`
  164. and :class:`GaussLegendre`). The extrapolation formula is given
  165. by Borwein, Bailey & Girgensohn. Although not very conservative,
  166. this method seems to be very robust in practice.
  167. """
  168. if len(results) == 2:
  169. return abs(results[0]-results[1])
  170. try:
  171. if results[-1] == results[-2] == results[-3]:
  172. return self.ctx.zero
  173. D1 = self.ctx.log(abs(results[-1]-results[-2]), 10)
  174. D2 = self.ctx.log(abs(results[-1]-results[-3]), 10)
  175. except ValueError:
  176. return epsilon
  177. D3 = -prec
  178. D4 = min(0, max(D1**2/D2, 2*D1, D3))
  179. return self.ctx.mpf(10) ** int(D4)
  180. def summation(self, f, points, prec, epsilon, max_degree, verbose=False):
  181. """
  182. Main integration function. Computes the 1D integral over
  183. the interval specified by *points*. For each subinterval,
  184. performs quadrature of degree from 1 up to *max_degree*
  185. until :func:`~mpmath.estimate_error` signals convergence.
  186. :func:`~mpmath.summation` transforms each subintegration to
  187. the standard interval and then calls :func:`~mpmath.sum_next`.
  188. """
  189. ctx = self.ctx
  190. I = total_err = ctx.zero
  191. for i in xrange(len(points)-1):
  192. a, b = points[i], points[i+1]
  193. if a == b:
  194. continue
  195. # XXX: we could use a single variable transformation,
  196. # but this is not good in practice. We get better accuracy
  197. # by having 0 as an endpoint.
  198. if (a, b) == (ctx.ninf, ctx.inf):
  199. _f = f
  200. f = lambda x: _f(-x) + _f(x)
  201. a, b = (ctx.zero, ctx.inf)
  202. results = []
  203. err = ctx.zero
  204. for degree in xrange(1, max_degree+1):
  205. nodes = self.get_nodes(a, b, degree, prec, verbose)
  206. if verbose:
  207. print("Integrating from %s to %s (degree %s of %s)" % \
  208. (ctx.nstr(a), ctx.nstr(b), degree, max_degree))
  209. result = self.sum_next(f, nodes, degree, prec, results, verbose)
  210. results.append(result)
  211. if degree > 1:
  212. err = self.estimate_error(results, prec, epsilon)
  213. if verbose:
  214. print("Estimated error:", ctx.nstr(err), " epsilon:", ctx.nstr(epsilon), " result: ", ctx.nstr(result))
  215. if err <= epsilon:
  216. break
  217. I += results[-1]
  218. total_err += err
  219. if total_err > epsilon:
  220. if verbose:
  221. print("Failed to reach full accuracy. Estimated error:", ctx.nstr(total_err))
  222. return I, total_err
  223. def sum_next(self, f, nodes, degree, prec, previous, verbose=False):
  224. r"""
  225. Evaluates the step sum `\sum w_k f(x_k)` where the *nodes* list
  226. contains the `(w_k, x_k)` pairs.
  227. :func:`~mpmath.summation` will supply the list *results* of
  228. values computed by :func:`~mpmath.sum_next` at previous degrees, in
  229. case the quadrature rule is able to reuse them.
  230. """
  231. return self.ctx.fdot((w, f(x)) for (x,w) in nodes)
  232. class TanhSinh(QuadratureRule):
  233. r"""
  234. This class implements "tanh-sinh" or "doubly exponential"
  235. quadrature. This quadrature rule is based on the Euler-Maclaurin
  236. integral formula. By performing a change of variables involving
  237. nested exponentials / hyperbolic functions (hence the name), the
  238. derivatives at the endpoints vanish rapidly. Since the error term
  239. in the Euler-Maclaurin formula depends on the derivatives at the
  240. endpoints, a simple step sum becomes extremely accurate. In
  241. practice, this means that doubling the number of evaluation
  242. points roughly doubles the number of accurate digits.
  243. Comparison to Gauss-Legendre:
  244. * Initial computation of nodes is usually faster
  245. * Handles endpoint singularities better
  246. * Handles infinite integration intervals better
  247. * Is slower for smooth integrands once nodes have been computed
  248. The implementation of the tanh-sinh algorithm is based on the
  249. description given in Borwein, Bailey & Girgensohn, "Experimentation
  250. in Mathematics - Computational Paths to Discovery", A K Peters,
  251. 2003, pages 312-313. In the present implementation, a few
  252. improvements have been made:
  253. * A more efficient scheme is used to compute nodes (exploiting
  254. recurrence for the exponential function)
  255. * The nodes are computed successively instead of all at once
  256. Various documents describing the algorithm are available online, e.g.:
  257. * http://crd.lbl.gov/~dhbailey/dhbpapers/dhb-tanh-sinh.pdf
  258. * http://users.cs.dal.ca/~jborwein/tanh-sinh.pdf
  259. """
  260. def sum_next(self, f, nodes, degree, prec, previous, verbose=False):
  261. """
  262. Step sum for tanh-sinh quadrature of degree `m`. We exploit the
  263. fact that half of the abscissas at degree `m` are precisely the
  264. abscissas from degree `m-1`. Thus reusing the result from
  265. the previous level allows a 2x speedup.
  266. """
  267. h = self.ctx.mpf(2)**(-degree)
  268. # Abscissas overlap, so reusing saves half of the time
  269. if previous:
  270. S = previous[-1]/(h*2)
  271. else:
  272. S = self.ctx.zero
  273. S += self.ctx.fdot((w,f(x)) for (x,w) in nodes)
  274. return h*S
  275. def calc_nodes(self, degree, prec, verbose=False):
  276. r"""
  277. The abscissas and weights for tanh-sinh quadrature of degree
  278. `m` are given by
  279. .. math::
  280. x_k = \tanh(\pi/2 \sinh(t_k))
  281. w_k = \pi/2 \cosh(t_k) / \cosh(\pi/2 \sinh(t_k))^2
  282. where `t_k = t_0 + hk` for a step length `h \sim 2^{-m}`. The
  283. list of nodes is actually infinite, but the weights die off so
  284. rapidly that only a few are needed.
  285. """
  286. ctx = self.ctx
  287. nodes = []
  288. extra = 20
  289. ctx.prec += extra
  290. tol = ctx.ldexp(1, -prec-10)
  291. pi4 = ctx.pi/4
  292. # For simplicity, we work in steps h = 1/2^n, with the first point
  293. # offset so that we can reuse the sum from the previous degree
  294. # We define degree 1 to include the "degree 0" steps, including
  295. # the point x = 0. (It doesn't work well otherwise; not sure why.)
  296. t0 = ctx.ldexp(1, -degree)
  297. if degree == 1:
  298. #nodes.append((mpf(0), pi4))
  299. #nodes.append((-mpf(0), pi4))
  300. nodes.append((ctx.zero, ctx.pi/2))
  301. h = t0
  302. else:
  303. h = t0*2
  304. # Since h is fixed, we can compute the next exponential
  305. # by simply multiplying by exp(h)
  306. expt0 = ctx.exp(t0)
  307. a = pi4 * expt0
  308. b = pi4 / expt0
  309. udelta = ctx.exp(h)
  310. urdelta = 1/udelta
  311. for k in xrange(0, 20*2**degree+1):
  312. # Reference implementation:
  313. # t = t0 + k*h
  314. # x = tanh(pi/2 * sinh(t))
  315. # w = pi/2 * cosh(t) / cosh(pi/2 * sinh(t))**2
  316. # Fast implementation. Note that c = exp(pi/2 * sinh(t))
  317. c = ctx.exp(a-b)
  318. d = 1/c
  319. co = (c+d)/2
  320. si = (c-d)/2
  321. x = si / co
  322. w = (a+b) / co**2
  323. diff = abs(x-1)
  324. if diff <= tol:
  325. break
  326. nodes.append((x, w))
  327. nodes.append((-x, w))
  328. a *= udelta
  329. b *= urdelta
  330. if verbose and k % 300 == 150:
  331. # Note: the number displayed is rather arbitrary. Should
  332. # figure out how to print something that looks more like a
  333. # percentage
  334. print("Calculating nodes:", ctx.nstr(-ctx.log(diff, 10) / prec))
  335. ctx.prec -= extra
  336. return nodes
  337. class GaussLegendre(QuadratureRule):
  338. r"""
  339. This class implements Gauss-Legendre quadrature, which is
  340. exceptionally efficient for polynomials and polynomial-like (i.e.
  341. very smooth) integrands.
  342. The abscissas and weights are given by roots and values of
  343. Legendre polynomials, which are the orthogonal polynomials
  344. on `[-1, 1]` with respect to the unit weight
  345. (see :func:`~mpmath.legendre`).
  346. In this implementation, we take the "degree" `m` of the quadrature
  347. to denote a Gauss-Legendre rule of degree `3 \cdot 2^m` (following
  348. Borwein, Bailey & Girgensohn). This way we get quadratic, rather
  349. than linear, convergence as the degree is incremented.
  350. Comparison to tanh-sinh quadrature:
  351. * Is faster for smooth integrands once nodes have been computed
  352. * Initial computation of nodes is usually slower
  353. * Handles endpoint singularities worse
  354. * Handles infinite integration intervals worse
  355. """
  356. def calc_nodes(self, degree, prec, verbose=False):
  357. r"""
  358. Calculates the abscissas and weights for Gauss-Legendre
  359. quadrature of degree of given degree (actually `3 \cdot 2^m`).
  360. """
  361. ctx = self.ctx
  362. # It is important that the epsilon is set lower than the
  363. # "real" epsilon
  364. epsilon = ctx.ldexp(1, -prec-8)
  365. # Fairly high precision might be required for accurate
  366. # evaluation of the roots
  367. orig = ctx.prec
  368. ctx.prec = int(prec*1.5)
  369. if degree == 1:
  370. x = ctx.sqrt(ctx.mpf(3)/5)
  371. w = ctx.mpf(5)/9
  372. nodes = [(-x,w),(ctx.zero,ctx.mpf(8)/9),(x,w)]
  373. ctx.prec = orig
  374. return nodes
  375. nodes = []
  376. n = 3*2**(degree-1)
  377. upto = n//2 + 1
  378. for j in xrange(1, upto):
  379. # Asymptotic formula for the roots
  380. r = ctx.mpf(math.cos(math.pi*(j-0.25)/(n+0.5)))
  381. # Newton iteration
  382. while 1:
  383. t1, t2 = 1, 0
  384. # Evaluates the Legendre polynomial using its defining
  385. # recurrence relation
  386. for j1 in xrange(1,n+1):
  387. t3, t2, t1 = t2, t1, ((2*j1-1)*r*t1 - (j1-1)*t2)/j1
  388. t4 = n*(r*t1-t2)/(r**2-1)
  389. a = t1/t4
  390. r = r - a
  391. if abs(a) < epsilon:
  392. break
  393. x = r
  394. w = 2/((1-r**2)*t4**2)
  395. if verbose and j % 30 == 15:
  396. print("Computing nodes (%i of %i)" % (j, upto))
  397. nodes.append((x, w))
  398. nodes.append((-x, w))
  399. ctx.prec = orig
  400. return nodes
  401. class QuadratureMethods(object):
  402. def __init__(ctx, *args, **kwargs):
  403. ctx._gauss_legendre = GaussLegendre(ctx)
  404. ctx._tanh_sinh = TanhSinh(ctx)
  405. def quad(ctx, f, *points, **kwargs):
  406. r"""
  407. Computes a single, double or triple integral over a given
  408. 1D interval, 2D rectangle, or 3D cuboid. A basic example::
  409. >>> from mpmath import *
  410. >>> mp.dps = 15; mp.pretty = True
  411. >>> quad(sin, [0, pi])
  412. 2.0
  413. A basic 2D integral::
  414. >>> f = lambda x, y: cos(x+y/2)
  415. >>> quad(f, [-pi/2, pi/2], [0, pi])
  416. 4.0
  417. **Interval format**
  418. The integration range for each dimension may be specified
  419. using a list or tuple. Arguments are interpreted as follows:
  420. ``quad(f, [x1, x2])`` -- calculates
  421. `\int_{x_1}^{x_2} f(x) \, dx`
  422. ``quad(f, [x1, x2], [y1, y2])`` -- calculates
  423. `\int_{x_1}^{x_2} \int_{y_1}^{y_2} f(x,y) \, dy \, dx`
  424. ``quad(f, [x1, x2], [y1, y2], [z1, z2])`` -- calculates
  425. `\int_{x_1}^{x_2} \int_{y_1}^{y_2} \int_{z_1}^{z_2} f(x,y,z)
  426. \, dz \, dy \, dx`
  427. Endpoints may be finite or infinite. An interval descriptor
  428. may also contain more than two points. In this
  429. case, the integration is split into subintervals, between
  430. each pair of consecutive points. This is useful for
  431. dealing with mid-interval discontinuities, or integrating
  432. over large intervals where the function is irregular or
  433. oscillates.
  434. **Options**
  435. :func:`~mpmath.quad` recognizes the following keyword arguments:
  436. *method*
  437. Chooses integration algorithm (described below).
  438. *error*
  439. If set to true, :func:`~mpmath.quad` returns `(v, e)` where `v` is the
  440. integral and `e` is the estimated error.
  441. *maxdegree*
  442. Maximum degree of the quadrature rule to try before
  443. quitting.
  444. *verbose*
  445. Print details about progress.
  446. **Algorithms**
  447. Mpmath presently implements two integration algorithms: tanh-sinh
  448. quadrature and Gauss-Legendre quadrature. These can be selected
  449. using *method='tanh-sinh'* or *method='gauss-legendre'* or by
  450. passing the classes *method=TanhSinh*, *method=GaussLegendre*.
  451. The functions :func:`~mpmath.quadts` and :func:`~mpmath.quadgl` are also available
  452. as shortcuts.
  453. Both algorithms have the property that doubling the number of
  454. evaluation points roughly doubles the accuracy, so both are ideal
  455. for high precision quadrature (hundreds or thousands of digits).
  456. At high precision, computing the nodes and weights for the
  457. integration can be expensive (more expensive than computing the
  458. function values). To make repeated integrations fast, nodes
  459. are automatically cached.
  460. The advantages of the tanh-sinh algorithm are that it tends to
  461. handle endpoint singularities well, and that the nodes are cheap
  462. to compute on the first run. For these reasons, it is used by
  463. :func:`~mpmath.quad` as the default algorithm.
  464. Gauss-Legendre quadrature often requires fewer function
  465. evaluations, and is therefore often faster for repeated use, but
  466. the algorithm does not handle endpoint singularities as well and
  467. the nodes are more expensive to compute. Gauss-Legendre quadrature
  468. can be a better choice if the integrand is smooth and repeated
  469. integrations are required (e.g. for multiple integrals).
  470. See the documentation for :class:`TanhSinh` and
  471. :class:`GaussLegendre` for additional details.
  472. **Examples of 1D integrals**
  473. Intervals may be infinite or half-infinite. The following two
  474. examples evaluate the limits of the inverse tangent function
  475. (`\int 1/(1+x^2) = \tan^{-1} x`), and the Gaussian integral
  476. `\int_{\infty}^{\infty} \exp(-x^2)\,dx = \sqrt{\pi}`::
  477. >>> mp.dps = 15
  478. >>> quad(lambda x: 2/(x**2+1), [0, inf])
  479. 3.14159265358979
  480. >>> quad(lambda x: exp(-x**2), [-inf, inf])**2
  481. 3.14159265358979
  482. Integrals can typically be resolved to high precision.
  483. The following computes 50 digits of `\pi` by integrating the
  484. area of the half-circle defined by `x^2 + y^2 \le 1`,
  485. `-1 \le x \le 1`, `y \ge 0`::
  486. >>> mp.dps = 50
  487. >>> 2*quad(lambda x: sqrt(1-x**2), [-1, 1])
  488. 3.1415926535897932384626433832795028841971693993751
  489. One can just as well compute 1000 digits (output truncated)::
  490. >>> mp.dps = 1000
  491. >>> 2*quad(lambda x: sqrt(1-x**2), [-1, 1]) #doctest:+ELLIPSIS
  492. 3.141592653589793238462643383279502884...216420199
  493. Complex integrals are supported. The following computes
  494. a residue at `z = 0` by integrating counterclockwise along the
  495. diamond-shaped path from `1` to `+i` to `-1` to `-i` to `1`::
  496. >>> mp.dps = 15
  497. >>> chop(quad(lambda z: 1/z, [1,j,-1,-j,1]))
  498. (0.0 + 6.28318530717959j)
  499. **Examples of 2D and 3D integrals**
  500. Here are several nice examples of analytically solvable
  501. 2D integrals (taken from MathWorld [1]) that can be evaluated
  502. to high precision fairly rapidly by :func:`~mpmath.quad`::
  503. >>> mp.dps = 30
  504. >>> f = lambda x, y: (x-1)/((1-x*y)*log(x*y))
  505. >>> quad(f, [0, 1], [0, 1])
  506. 0.577215664901532860606512090082
  507. >>> +euler
  508. 0.577215664901532860606512090082
  509. >>> f = lambda x, y: 1/sqrt(1+x**2+y**2)
  510. >>> quad(f, [-1, 1], [-1, 1])
  511. 3.17343648530607134219175646705
  512. >>> 4*log(2+sqrt(3))-2*pi/3
  513. 3.17343648530607134219175646705
  514. >>> f = lambda x, y: 1/(1-x**2 * y**2)
  515. >>> quad(f, [0, 1], [0, 1])
  516. 1.23370055013616982735431137498
  517. >>> pi**2 / 8
  518. 1.23370055013616982735431137498
  519. >>> quad(lambda x, y: 1/(1-x*y), [0, 1], [0, 1])
  520. 1.64493406684822643647241516665
  521. >>> pi**2 / 6
  522. 1.64493406684822643647241516665
  523. Multiple integrals may be done over infinite ranges::
  524. >>> mp.dps = 15
  525. >>> print(quad(lambda x,y: exp(-x-y), [0, inf], [1, inf]))
  526. 0.367879441171442
  527. >>> print(1/e)
  528. 0.367879441171442
  529. For nonrectangular areas, one can call :func:`~mpmath.quad` recursively.
  530. For example, we can replicate the earlier example of calculating
  531. `\pi` by integrating over the unit-circle, and actually use double
  532. quadrature to actually measure the area circle::
  533. >>> f = lambda x: quad(lambda y: 1, [-sqrt(1-x**2), sqrt(1-x**2)])
  534. >>> quad(f, [-1, 1])
  535. 3.14159265358979
  536. Here is a simple triple integral::
  537. >>> mp.dps = 15
  538. >>> f = lambda x,y,z: x*y/(1+z)
  539. >>> quad(f, [0,1], [0,1], [1,2], method='gauss-legendre')
  540. 0.101366277027041
  541. >>> (log(3)-log(2))/4
  542. 0.101366277027041
  543. **Singularities**
  544. Both tanh-sinh and Gauss-Legendre quadrature are designed to
  545. integrate smooth (infinitely differentiable) functions. Neither
  546. algorithm copes well with mid-interval singularities (such as
  547. mid-interval discontinuities in `f(x)` or `f'(x)`).
  548. The best solution is to split the integral into parts::
  549. >>> mp.dps = 15
  550. >>> quad(lambda x: abs(sin(x)), [0, 2*pi]) # Bad
  551. 3.99900894176779
  552. >>> quad(lambda x: abs(sin(x)), [0, pi, 2*pi]) # Good
  553. 4.0
  554. The tanh-sinh rule often works well for integrands having a
  555. singularity at one or both endpoints::
  556. >>> mp.dps = 15
  557. >>> quad(log, [0, 1], method='tanh-sinh') # Good
  558. -1.0
  559. >>> quad(log, [0, 1], method='gauss-legendre') # Bad
  560. -0.999932197413801
  561. However, the result may still be inaccurate for some functions::
  562. >>> quad(lambda x: 1/sqrt(x), [0, 1], method='tanh-sinh')
  563. 1.99999999946942
  564. This problem is not due to the quadrature rule per se, but to
  565. numerical amplification of errors in the nodes. The problem can be
  566. circumvented by temporarily increasing the precision::
  567. >>> mp.dps = 30
  568. >>> a = quad(lambda x: 1/sqrt(x), [0, 1], method='tanh-sinh')
  569. >>> mp.dps = 15
  570. >>> +a
  571. 2.0
  572. **Highly variable functions**
  573. For functions that are smooth (in the sense of being infinitely
  574. differentiable) but contain sharp mid-interval peaks or many
  575. "bumps", :func:`~mpmath.quad` may fail to provide full accuracy. For
  576. example, with default settings, :func:`~mpmath.quad` is able to integrate
  577. `\sin(x)` accurately over an interval of length 100 but not over
  578. length 1000::
  579. >>> quad(sin, [0, 100]); 1-cos(100) # Good
  580. 0.137681127712316
  581. 0.137681127712316
  582. >>> quad(sin, [0, 1000]); 1-cos(1000) # Bad
  583. -37.8587612408485
  584. 0.437620923709297
  585. One solution is to break the integration into 10 intervals of
  586. length 100::
  587. >>> quad(sin, linspace(0, 1000, 10)) # Good
  588. 0.437620923709297
  589. Another is to increase the degree of the quadrature::
  590. >>> quad(sin, [0, 1000], maxdegree=10) # Also good
  591. 0.437620923709297
  592. Whether splitting the interval or increasing the degree is
  593. more efficient differs from case to case. Another example is the
  594. function `1/(1+x^2)`, which has a sharp peak centered around
  595. `x = 0`::
  596. >>> f = lambda x: 1/(1+x**2)
  597. >>> quad(f, [-100, 100]) # Bad
  598. 3.64804647105268
  599. >>> quad(f, [-100, 100], maxdegree=10) # Good
  600. 3.12159332021646
  601. >>> quad(f, [-100, 0, 100]) # Also good
  602. 3.12159332021646
  603. **References**
  604. 1. http://mathworld.wolfram.com/DoubleIntegral.html
  605. """
  606. rule = kwargs.get('method', 'tanh-sinh')
  607. if type(rule) is str:
  608. if rule == 'tanh-sinh':
  609. rule = ctx._tanh_sinh
  610. elif rule == 'gauss-legendre':
  611. rule = ctx._gauss_legendre
  612. else:
  613. raise ValueError("unknown quadrature rule: %s" % rule)
  614. else:
  615. rule = rule(ctx)
  616. verbose = kwargs.get('verbose')
  617. dim = len(points)
  618. orig = prec = ctx.prec
  619. epsilon = ctx.eps/8
  620. m = kwargs.get('maxdegree') or rule.guess_degree(prec)
  621. points = [ctx._as_points(p) for p in points]
  622. try:
  623. ctx.prec += 20
  624. if dim == 1:
  625. v, err = rule.summation(f, points[0], prec, epsilon, m, verbose)
  626. elif dim == 2:
  627. v, err = rule.summation(lambda x: \
  628. rule.summation(lambda y: f(x,y), \
  629. points[1], prec, epsilon, m)[0],
  630. points[0], prec, epsilon, m, verbose)
  631. elif dim == 3:
  632. v, err = rule.summation(lambda x: \
  633. rule.summation(lambda y: \
  634. rule.summation(lambda z: f(x,y,z), \
  635. points[2], prec, epsilon, m)[0],
  636. points[1], prec, epsilon, m)[0],
  637. points[0], prec, epsilon, m, verbose)
  638. else:
  639. raise NotImplementedError("quadrature must have dim 1, 2 or 3")
  640. finally:
  641. ctx.prec = orig
  642. if kwargs.get("error"):
  643. return +v, err
  644. return +v
  645. def quadts(ctx, *args, **kwargs):
  646. """
  647. Performs tanh-sinh quadrature. The call
  648. quadts(func, *points, ...)
  649. is simply a shortcut for:
  650. quad(func, *points, ..., method=TanhSinh)
  651. For example, a single integral and a double integral:
  652. quadts(lambda x: exp(cos(x)), [0, 1])
  653. quadts(lambda x, y: exp(cos(x+y)), [0, 1], [0, 1])
  654. See the documentation for quad for information about how points
  655. arguments and keyword arguments are parsed.
  656. See documentation for TanhSinh for algorithmic information about
  657. tanh-sinh quadrature.
  658. """
  659. kwargs['method'] = 'tanh-sinh'
  660. return ctx.quad(*args, **kwargs)
  661. def quadgl(ctx, *args, **kwargs):
  662. """
  663. Performs Gauss-Legendre quadrature. The call
  664. quadgl(func, *points, ...)
  665. is simply a shortcut for:
  666. quad(func, *points, ..., method=GaussLegendre)
  667. For example, a single integral and a double integral:
  668. quadgl(lambda x: exp(cos(x)), [0, 1])
  669. quadgl(lambda x, y: exp(cos(x+y)), [0, 1], [0, 1])
  670. See the documentation for quad for information about how points
  671. arguments and keyword arguments are parsed.
  672. See documentation for TanhSinh for algorithmic information about
  673. tanh-sinh quadrature.
  674. """
  675. kwargs['method'] = 'gauss-legendre'
  676. return ctx.quad(*args, **kwargs)
  677. def quadosc(ctx, f, interval, omega=None, period=None, zeros=None):
  678. r"""
  679. Calculates
  680. .. math ::
  681. I = \int_a^b f(x) dx
  682. where at least one of `a` and `b` is infinite and where
  683. `f(x) = g(x) \cos(\omega x + \phi)` for some slowly
  684. decreasing function `g(x)`. With proper input, :func:`~mpmath.quadosc`
  685. can also handle oscillatory integrals where the oscillation
  686. rate is different from a pure sine or cosine wave.
  687. In the standard case when `|a| < \infty, b = \infty`,
  688. :func:`~mpmath.quadosc` works by evaluating the infinite series
  689. .. math ::
  690. I = \int_a^{x_1} f(x) dx +
  691. \sum_{k=1}^{\infty} \int_{x_k}^{x_{k+1}} f(x) dx
  692. where `x_k` are consecutive zeros (alternatively
  693. some other periodic reference point) of `f(x)`.
  694. Accordingly, :func:`~mpmath.quadosc` requires information about the
  695. zeros of `f(x)`. For a periodic function, you can specify
  696. the zeros by either providing the angular frequency `\omega`
  697. (*omega*) or the *period* `2 \pi/\omega`. In general, you can
  698. specify the `n`-th zero by providing the *zeros* arguments.
  699. Below is an example of each::
  700. >>> from mpmath import *
  701. >>> mp.dps = 15; mp.pretty = True
  702. >>> f = lambda x: sin(3*x)/(x**2+1)
  703. >>> quadosc(f, [0,inf], omega=3)
  704. 0.37833007080198
  705. >>> quadosc(f, [0,inf], period=2*pi/3)
  706. 0.37833007080198
  707. >>> quadosc(f, [0,inf], zeros=lambda n: pi*n/3)
  708. 0.37833007080198
  709. >>> (ei(3)*exp(-3)-exp(3)*ei(-3))/2 # Computed by Mathematica
  710. 0.37833007080198
  711. Note that *zeros* was specified to multiply `n` by the
  712. *half-period*, not the full period. In theory, it does not matter
  713. whether each partial integral is done over a half period or a full
  714. period. However, if done over half-periods, the infinite series
  715. passed to :func:`~mpmath.nsum` becomes an *alternating series* and this
  716. typically makes the extrapolation much more efficient.
  717. Here is an example of an integration over the entire real line,
  718. and a half-infinite integration starting at `-\infty`::
  719. >>> quadosc(lambda x: cos(x)/(1+x**2), [-inf, inf], omega=1)
  720. 1.15572734979092
  721. >>> pi/e
  722. 1.15572734979092
  723. >>> quadosc(lambda x: cos(x)/x**2, [-inf, -1], period=2*pi)
  724. -0.0844109505595739
  725. >>> cos(1)+si(1)-pi/2
  726. -0.0844109505595738
  727. Of course, the integrand may contain a complex exponential just as
  728. well as a real sine or cosine::
  729. >>> quadosc(lambda x: exp(3*j*x)/(1+x**2), [-inf,inf], omega=3)
  730. (0.156410688228254 + 0.0j)
  731. >>> pi/e**3
  732. 0.156410688228254
  733. >>> quadosc(lambda x: exp(3*j*x)/(2+x+x**2), [-inf,inf], omega=3)
  734. (0.00317486988463794 - 0.0447701735209082j)
  735. >>> 2*pi/sqrt(7)/exp(3*(j+sqrt(7))/2)
  736. (0.00317486988463794 - 0.0447701735209082j)
  737. **Non-periodic functions**
  738. If `f(x) = g(x) h(x)` for some function `h(x)` that is not
  739. strictly periodic, *omega* or *period* might not work, and it might
  740. be necessary to use *zeros*.
  741. A notable exception can be made for Bessel functions which, though not
  742. periodic, are "asymptotically periodic" in a sufficiently strong sense
  743. that the sum extrapolation will work out::
  744. >>> quadosc(j0, [0, inf], period=2*pi)
  745. 1.0
  746. >>> quadosc(j1, [0, inf], period=2*pi)
  747. 1.0
  748. More properly, one should provide the exact Bessel function zeros::
  749. >>> j0zero = lambda n: findroot(j0, pi*(n-0.25))
  750. >>> quadosc(j0, [0, inf], zeros=j0zero)
  751. 1.0
  752. For an example where *zeros* becomes necessary, consider the
  753. complete Fresnel integrals
  754. .. math ::
  755. \int_0^{\infty} \cos x^2\,dx = \int_0^{\infty} \sin x^2\,dx
  756. = \sqrt{\frac{\pi}{8}}.
  757. Although the integrands do not decrease in magnitude as
  758. `x \to \infty`, the integrals are convergent since the oscillation
  759. rate increases (causing consecutive periods to asymptotically
  760. cancel out). These integrals are virtually impossible to calculate
  761. to any kind of accuracy using standard quadrature rules. However,
  762. if one provides the correct asymptotic distribution of zeros
  763. (`x_n \sim \sqrt{n}`), :func:`~mpmath.quadosc` works::
  764. >>> mp.dps = 30
  765. >>> f = lambda x: cos(x**2)
  766. >>> quadosc(f, [0,inf], zeros=lambda n:sqrt(pi*n))
  767. 0.626657068657750125603941321203
  768. >>> f = lambda x: sin(x**2)
  769. >>> quadosc(f, [0,inf], zeros=lambda n:sqrt(pi*n))
  770. 0.626657068657750125603941321203
  771. >>> sqrt(pi/8)
  772. 0.626657068657750125603941321203
  773. (Interestingly, these integrals can still be evaluated if one
  774. places some other constant than `\pi` in the square root sign.)
  775. In general, if `f(x) \sim g(x) \cos(h(x))`, the zeros follow
  776. the inverse-function distribution `h^{-1}(x)`::
  777. >>> mp.dps = 15
  778. >>> f = lambda x: sin(exp(x))
  779. >>> quadosc(f, [1,inf], zeros=lambda n: log(n))
  780. -0.25024394235267
  781. >>> pi/2-si(e)
  782. -0.250243942352671
  783. **Non-alternating functions**
  784. If the integrand oscillates around a positive value, without
  785. alternating signs, the extrapolation might fail. A simple trick
  786. that sometimes works is to multiply or divide the frequency by 2::
  787. >>> f = lambda x: 1/x**2+sin(x)/x**4
  788. >>> quadosc(f, [1,inf], omega=1) # Bad
  789. 1.28642190869861
  790. >>> quadosc(f, [1,inf], omega=0.5) # Perfect
  791. 1.28652953559617
  792. >>> 1+(cos(1)+ci(1)+sin(1))/6
  793. 1.28652953559617
  794. **Fast decay**
  795. :func:`~mpmath.quadosc` is primarily useful for slowly decaying
  796. integrands. If the integrand decreases exponentially or faster,
  797. :func:`~mpmath.quad` will likely handle it without trouble (and generally be
  798. much faster than :func:`~mpmath.quadosc`)::
  799. >>> quadosc(lambda x: cos(x)/exp(x), [0, inf], omega=1)
  800. 0.5
  801. >>> quad(lambda x: cos(x)/exp(x), [0, inf])
  802. 0.5
  803. """
  804. a, b = ctx._as_points(interval)
  805. a = ctx.convert(a)
  806. b = ctx.convert(b)
  807. if [omega, period, zeros].count(None) != 2:
  808. raise ValueError( \
  809. "must specify exactly one of omega, period, zeros")
  810. if a == ctx.ninf and b == ctx.inf:
  811. s1 = ctx.quadosc(f, [a, 0], omega=omega, zeros=zeros, period=period)
  812. s2 = ctx.quadosc(f, [0, b], omega=omega, zeros=zeros, period=period)
  813. return s1 + s2
  814. if a == ctx.ninf:
  815. if zeros:
  816. return ctx.quadosc(lambda x:f(-x), [-b,-a], lambda n: zeros(-n))
  817. else:
  818. return ctx.quadosc(lambda x:f(-x), [-b,-a], omega=omega, period=period)
  819. if b != ctx.inf:
  820. raise ValueError("quadosc requires an infinite integration interval")
  821. if not zeros:
  822. if omega:
  823. period = 2*ctx.pi/omega
  824. zeros = lambda n: n*period/2
  825. #for n in range(1,10):
  826. # p = zeros(n)
  827. # if p > a:
  828. # break
  829. #if n >= 9:
  830. # raise ValueError("zeros do not appear to be correctly indexed")
  831. n = 1
  832. s = ctx.quadgl(f, [a, zeros(n)])
  833. def term(k):
  834. return ctx.quadgl(f, [zeros(k), zeros(k+1)])
  835. s += ctx.nsum(term, [n, ctx.inf])
  836. return s
  837. if __name__ == '__main__':
  838. import doctest
  839. doctest.testmod()