quadrature.py 17 KB

123456789101112131415161718192021222324252627282930313233343536373839404142434445464748495051525354555657585960616263646566676869707172737475767778798081828384858687888990919293949596979899100101102103104105106107108109110111112113114115116117118119120121122123124125126127128129130131132133134135136137138139140141142143144145146147148149150151152153154155156157158159160161162163164165166167168169170171172173174175176177178179180181182183184185186187188189190191192193194195196197198199200201202203204205206207208209210211212213214215216217218219220221222223224225226227228229230231232233234235236237238239240241242243244245246247248249250251252253254255256257258259260261262263264265266267268269270271272273274275276277278279280281282283284285286287288289290291292293294295296297298299300301302303304305306307308309310311312313314315316317318319320321322323324325326327328329330331332333334335336337338339340341342343344345346347348349350351352353354355356357358359360361362363364365366367368369370371372373374375376377378379380381382383384385386387388389390391392393394395396397398399400401402403404405406407408409410411412413414415416417418419420421422423424425426427428429430431432433434435436437438439440441442443444445446447448449450451452453454455456457458459460461462463464465466467468469470471472473474475476477478479480481482483484485486487488489490491492493494495496497498499500501502503504505506507508509510511512513514515516517518519520521522523524525526527528529530531532533534535536537538539540541542543544545546547548549550551552553554555556557558559560561562563564565566567568569570571572573574575576577578579580581582583584585586587588589590591592593594595596597598599600601602603604605606607608609610611612613614615616617
  1. from sympy.core import S, Dummy, pi
  2. from sympy.functions.combinatorial.factorials import factorial
  3. from sympy.functions.elementary.trigonometric import sin, cos
  4. from sympy.functions.elementary.miscellaneous import sqrt
  5. from sympy.functions.special.gamma_functions import gamma
  6. from sympy.polys.orthopolys import (legendre_poly, laguerre_poly,
  7. hermite_poly, jacobi_poly)
  8. from sympy.polys.rootoftools import RootOf
  9. def gauss_legendre(n, n_digits):
  10. r"""
  11. Computes the Gauss-Legendre quadrature [1]_ points and weights.
  12. Explanation
  13. ===========
  14. The Gauss-Legendre quadrature approximates the integral:
  15. .. math::
  16. \int_{-1}^1 f(x)\,dx \approx \sum_{i=1}^n w_i f(x_i)
  17. The nodes `x_i` of an order `n` quadrature rule are the roots of `P_n`
  18. and the weights `w_i` are given by:
  19. .. math::
  20. w_i = \frac{2}{\left(1-x_i^2\right) \left(P'_n(x_i)\right)^2}
  21. Parameters
  22. ==========
  23. n :
  24. The order of quadrature.
  25. n_digits :
  26. Number of significant digits of the points and weights to return.
  27. Returns
  28. =======
  29. (x, w) : the ``x`` and ``w`` are lists of points and weights as Floats.
  30. The points `x_i` and weights `w_i` are returned as ``(x, w)``
  31. tuple of lists.
  32. Examples
  33. ========
  34. >>> from sympy.integrals.quadrature import gauss_legendre
  35. >>> x, w = gauss_legendre(3, 5)
  36. >>> x
  37. [-0.7746, 0, 0.7746]
  38. >>> w
  39. [0.55556, 0.88889, 0.55556]
  40. >>> x, w = gauss_legendre(4, 5)
  41. >>> x
  42. [-0.86114, -0.33998, 0.33998, 0.86114]
  43. >>> w
  44. [0.34785, 0.65215, 0.65215, 0.34785]
  45. See Also
  46. ========
  47. gauss_laguerre, gauss_gen_laguerre, gauss_hermite, gauss_chebyshev_t, gauss_chebyshev_u, gauss_jacobi, gauss_lobatto
  48. References
  49. ==========
  50. .. [1] https://en.wikipedia.org/wiki/Gaussian_quadrature
  51. .. [2] http://people.sc.fsu.edu/~jburkardt/cpp_src/legendre_rule/legendre_rule.html
  52. """
  53. x = Dummy("x")
  54. p = legendre_poly(n, x, polys=True)
  55. pd = p.diff(x)
  56. xi = []
  57. w = []
  58. for r in p.real_roots():
  59. if isinstance(r, RootOf):
  60. r = r.eval_rational(S.One/10**(n_digits+2))
  61. xi.append(r.n(n_digits))
  62. w.append((2/((1-r**2) * pd.subs(x, r)**2)).n(n_digits))
  63. return xi, w
  64. def gauss_laguerre(n, n_digits):
  65. r"""
  66. Computes the Gauss-Laguerre quadrature [1]_ points and weights.
  67. Explanation
  68. ===========
  69. The Gauss-Laguerre quadrature approximates the integral:
  70. .. math::
  71. \int_0^{\infty} e^{-x} f(x)\,dx \approx \sum_{i=1}^n w_i f(x_i)
  72. The nodes `x_i` of an order `n` quadrature rule are the roots of `L_n`
  73. and the weights `w_i` are given by:
  74. .. math::
  75. w_i = \frac{x_i}{(n+1)^2 \left(L_{n+1}(x_i)\right)^2}
  76. Parameters
  77. ==========
  78. n :
  79. The order of quadrature.
  80. n_digits :
  81. Number of significant digits of the points and weights to return.
  82. Returns
  83. =======
  84. (x, w) : The ``x`` and ``w`` are lists of points and weights as Floats.
  85. The points `x_i` and weights `w_i` are returned as ``(x, w)``
  86. tuple of lists.
  87. Examples
  88. ========
  89. >>> from sympy.integrals.quadrature import gauss_laguerre
  90. >>> x, w = gauss_laguerre(3, 5)
  91. >>> x
  92. [0.41577, 2.2943, 6.2899]
  93. >>> w
  94. [0.71109, 0.27852, 0.010389]
  95. >>> x, w = gauss_laguerre(6, 5)
  96. >>> x
  97. [0.22285, 1.1889, 2.9927, 5.7751, 9.8375, 15.983]
  98. >>> w
  99. [0.45896, 0.417, 0.11337, 0.010399, 0.00026102, 8.9855e-7]
  100. See Also
  101. ========
  102. gauss_legendre, gauss_gen_laguerre, gauss_hermite, gauss_chebyshev_t, gauss_chebyshev_u, gauss_jacobi, gauss_lobatto
  103. References
  104. ==========
  105. .. [1] https://en.wikipedia.org/wiki/Gauss%E2%80%93Laguerre_quadrature
  106. .. [2] http://people.sc.fsu.edu/~jburkardt/cpp_src/laguerre_rule/laguerre_rule.html
  107. """
  108. x = Dummy("x")
  109. p = laguerre_poly(n, x, polys=True)
  110. p1 = laguerre_poly(n+1, x, polys=True)
  111. xi = []
  112. w = []
  113. for r in p.real_roots():
  114. if isinstance(r, RootOf):
  115. r = r.eval_rational(S.One/10**(n_digits+2))
  116. xi.append(r.n(n_digits))
  117. w.append((r/((n+1)**2 * p1.subs(x, r)**2)).n(n_digits))
  118. return xi, w
  119. def gauss_hermite(n, n_digits):
  120. r"""
  121. Computes the Gauss-Hermite quadrature [1]_ points and weights.
  122. Explanation
  123. ===========
  124. The Gauss-Hermite quadrature approximates the integral:
  125. .. math::
  126. \int_{-\infty}^{\infty} e^{-x^2} f(x)\,dx \approx
  127. \sum_{i=1}^n w_i f(x_i)
  128. The nodes `x_i` of an order `n` quadrature rule are the roots of `H_n`
  129. and the weights `w_i` are given by:
  130. .. math::
  131. w_i = \frac{2^{n-1} n! \sqrt{\pi}}{n^2 \left(H_{n-1}(x_i)\right)^2}
  132. Parameters
  133. ==========
  134. n :
  135. The order of quadrature.
  136. n_digits :
  137. Number of significant digits of the points and weights to return.
  138. Returns
  139. =======
  140. (x, w) : The ``x`` and ``w`` are lists of points and weights as Floats.
  141. The points `x_i` and weights `w_i` are returned as ``(x, w)``
  142. tuple of lists.
  143. Examples
  144. ========
  145. >>> from sympy.integrals.quadrature import gauss_hermite
  146. >>> x, w = gauss_hermite(3, 5)
  147. >>> x
  148. [-1.2247, 0, 1.2247]
  149. >>> w
  150. [0.29541, 1.1816, 0.29541]
  151. >>> x, w = gauss_hermite(6, 5)
  152. >>> x
  153. [-2.3506, -1.3358, -0.43608, 0.43608, 1.3358, 2.3506]
  154. >>> w
  155. [0.00453, 0.15707, 0.72463, 0.72463, 0.15707, 0.00453]
  156. See Also
  157. ========
  158. gauss_legendre, gauss_laguerre, gauss_gen_laguerre, gauss_chebyshev_t, gauss_chebyshev_u, gauss_jacobi, gauss_lobatto
  159. References
  160. ==========
  161. .. [1] https://en.wikipedia.org/wiki/Gauss-Hermite_Quadrature
  162. .. [2] http://people.sc.fsu.edu/~jburkardt/cpp_src/hermite_rule/hermite_rule.html
  163. .. [3] http://people.sc.fsu.edu/~jburkardt/cpp_src/gen_hermite_rule/gen_hermite_rule.html
  164. """
  165. x = Dummy("x")
  166. p = hermite_poly(n, x, polys=True)
  167. p1 = hermite_poly(n-1, x, polys=True)
  168. xi = []
  169. w = []
  170. for r in p.real_roots():
  171. if isinstance(r, RootOf):
  172. r = r.eval_rational(S.One/10**(n_digits+2))
  173. xi.append(r.n(n_digits))
  174. w.append(((2**(n-1) * factorial(n) * sqrt(pi)) /
  175. (n**2 * p1.subs(x, r)**2)).n(n_digits))
  176. return xi, w
  177. def gauss_gen_laguerre(n, alpha, n_digits):
  178. r"""
  179. Computes the generalized Gauss-Laguerre quadrature [1]_ points and weights.
  180. Explanation
  181. ===========
  182. The generalized Gauss-Laguerre quadrature approximates the integral:
  183. .. math::
  184. \int_{0}^\infty x^{\alpha} e^{-x} f(x)\,dx \approx
  185. \sum_{i=1}^n w_i f(x_i)
  186. The nodes `x_i` of an order `n` quadrature rule are the roots of
  187. `L^{\alpha}_n` and the weights `w_i` are given by:
  188. .. math::
  189. w_i = \frac{\Gamma(\alpha+n)}
  190. {n \Gamma(n) L^{\alpha}_{n-1}(x_i) L^{\alpha+1}_{n-1}(x_i)}
  191. Parameters
  192. ==========
  193. n :
  194. The order of quadrature.
  195. alpha :
  196. The exponent of the singularity, `\alpha > -1`.
  197. n_digits :
  198. Number of significant digits of the points and weights to return.
  199. Returns
  200. =======
  201. (x, w) : the ``x`` and ``w`` are lists of points and weights as Floats.
  202. The points `x_i` and weights `w_i` are returned as ``(x, w)``
  203. tuple of lists.
  204. Examples
  205. ========
  206. >>> from sympy import S
  207. >>> from sympy.integrals.quadrature import gauss_gen_laguerre
  208. >>> x, w = gauss_gen_laguerre(3, -S.Half, 5)
  209. >>> x
  210. [0.19016, 1.7845, 5.5253]
  211. >>> w
  212. [1.4493, 0.31413, 0.00906]
  213. >>> x, w = gauss_gen_laguerre(4, 3*S.Half, 5)
  214. >>> x
  215. [0.97851, 2.9904, 6.3193, 11.712]
  216. >>> w
  217. [0.53087, 0.67721, 0.11895, 0.0023152]
  218. See Also
  219. ========
  220. gauss_legendre, gauss_laguerre, gauss_hermite, gauss_chebyshev_t, gauss_chebyshev_u, gauss_jacobi, gauss_lobatto
  221. References
  222. ==========
  223. .. [1] https://en.wikipedia.org/wiki/Gauss%E2%80%93Laguerre_quadrature
  224. .. [2] http://people.sc.fsu.edu/~jburkardt/cpp_src/gen_laguerre_rule/gen_laguerre_rule.html
  225. """
  226. x = Dummy("x")
  227. p = laguerre_poly(n, x, alpha=alpha, polys=True)
  228. p1 = laguerre_poly(n-1, x, alpha=alpha, polys=True)
  229. p2 = laguerre_poly(n-1, x, alpha=alpha+1, polys=True)
  230. xi = []
  231. w = []
  232. for r in p.real_roots():
  233. if isinstance(r, RootOf):
  234. r = r.eval_rational(S.One/10**(n_digits+2))
  235. xi.append(r.n(n_digits))
  236. w.append((gamma(alpha+n) /
  237. (n*gamma(n)*p1.subs(x, r)*p2.subs(x, r))).n(n_digits))
  238. return xi, w
  239. def gauss_chebyshev_t(n, n_digits):
  240. r"""
  241. Computes the Gauss-Chebyshev quadrature [1]_ points and weights of
  242. the first kind.
  243. Explanation
  244. ===========
  245. The Gauss-Chebyshev quadrature of the first kind approximates the integral:
  246. .. math::
  247. \int_{-1}^{1} \frac{1}{\sqrt{1-x^2}} f(x)\,dx \approx
  248. \sum_{i=1}^n w_i f(x_i)
  249. The nodes `x_i` of an order `n` quadrature rule are the roots of `T_n`
  250. and the weights `w_i` are given by:
  251. .. math::
  252. w_i = \frac{\pi}{n}
  253. Parameters
  254. ==========
  255. n :
  256. The order of quadrature.
  257. n_digits :
  258. Number of significant digits of the points and weights to return.
  259. Returns
  260. =======
  261. (x, w) : the ``x`` and ``w`` are lists of points and weights as Floats.
  262. The points `x_i` and weights `w_i` are returned as ``(x, w)``
  263. tuple of lists.
  264. Examples
  265. ========
  266. >>> from sympy.integrals.quadrature import gauss_chebyshev_t
  267. >>> x, w = gauss_chebyshev_t(3, 5)
  268. >>> x
  269. [0.86602, 0, -0.86602]
  270. >>> w
  271. [1.0472, 1.0472, 1.0472]
  272. >>> x, w = gauss_chebyshev_t(6, 5)
  273. >>> x
  274. [0.96593, 0.70711, 0.25882, -0.25882, -0.70711, -0.96593]
  275. >>> w
  276. [0.5236, 0.5236, 0.5236, 0.5236, 0.5236, 0.5236]
  277. See Also
  278. ========
  279. gauss_legendre, gauss_laguerre, gauss_hermite, gauss_gen_laguerre, gauss_chebyshev_u, gauss_jacobi, gauss_lobatto
  280. References
  281. ==========
  282. .. [1] https://en.wikipedia.org/wiki/Chebyshev%E2%80%93Gauss_quadrature
  283. .. [2] http://people.sc.fsu.edu/~jburkardt/cpp_src/chebyshev1_rule/chebyshev1_rule.html
  284. """
  285. xi = []
  286. w = []
  287. for i in range(1, n+1):
  288. xi.append((cos((2*i-S.One)/(2*n)*S.Pi)).n(n_digits))
  289. w.append((S.Pi/n).n(n_digits))
  290. return xi, w
  291. def gauss_chebyshev_u(n, n_digits):
  292. r"""
  293. Computes the Gauss-Chebyshev quadrature [1]_ points and weights of
  294. the second kind.
  295. Explanation
  296. ===========
  297. The Gauss-Chebyshev quadrature of the second kind approximates the
  298. integral:
  299. .. math::
  300. \int_{-1}^{1} \sqrt{1-x^2} f(x)\,dx \approx \sum_{i=1}^n w_i f(x_i)
  301. The nodes `x_i` of an order `n` quadrature rule are the roots of `U_n`
  302. and the weights `w_i` are given by:
  303. .. math::
  304. w_i = \frac{\pi}{n+1} \sin^2 \left(\frac{i}{n+1}\pi\right)
  305. Parameters
  306. ==========
  307. n : the order of quadrature
  308. n_digits : number of significant digits of the points and weights to return
  309. Returns
  310. =======
  311. (x, w) : the ``x`` and ``w`` are lists of points and weights as Floats.
  312. The points `x_i` and weights `w_i` are returned as ``(x, w)``
  313. tuple of lists.
  314. Examples
  315. ========
  316. >>> from sympy.integrals.quadrature import gauss_chebyshev_u
  317. >>> x, w = gauss_chebyshev_u(3, 5)
  318. >>> x
  319. [0.70711, 0, -0.70711]
  320. >>> w
  321. [0.3927, 0.7854, 0.3927]
  322. >>> x, w = gauss_chebyshev_u(6, 5)
  323. >>> x
  324. [0.90097, 0.62349, 0.22252, -0.22252, -0.62349, -0.90097]
  325. >>> w
  326. [0.084489, 0.27433, 0.42658, 0.42658, 0.27433, 0.084489]
  327. See Also
  328. ========
  329. gauss_legendre, gauss_laguerre, gauss_hermite, gauss_gen_laguerre, gauss_chebyshev_t, gauss_jacobi, gauss_lobatto
  330. References
  331. ==========
  332. .. [1] https://en.wikipedia.org/wiki/Chebyshev%E2%80%93Gauss_quadrature
  333. .. [2] http://people.sc.fsu.edu/~jburkardt/cpp_src/chebyshev2_rule/chebyshev2_rule.html
  334. """
  335. xi = []
  336. w = []
  337. for i in range(1, n+1):
  338. xi.append((cos(i/(n+S.One)*S.Pi)).n(n_digits))
  339. w.append((S.Pi/(n+S.One)*sin(i*S.Pi/(n+S.One))**2).n(n_digits))
  340. return xi, w
  341. def gauss_jacobi(n, alpha, beta, n_digits):
  342. r"""
  343. Computes the Gauss-Jacobi quadrature [1]_ points and weights.
  344. Explanation
  345. ===========
  346. The Gauss-Jacobi quadrature of the first kind approximates the integral:
  347. .. math::
  348. \int_{-1}^1 (1-x)^\alpha (1+x)^\beta f(x)\,dx \approx
  349. \sum_{i=1}^n w_i f(x_i)
  350. The nodes `x_i` of an order `n` quadrature rule are the roots of
  351. `P^{(\alpha,\beta)}_n` and the weights `w_i` are given by:
  352. .. math::
  353. w_i = -\frac{2n+\alpha+\beta+2}{n+\alpha+\beta+1}
  354. \frac{\Gamma(n+\alpha+1)\Gamma(n+\beta+1)}
  355. {\Gamma(n+\alpha+\beta+1)(n+1)!}
  356. \frac{2^{\alpha+\beta}}{P'_n(x_i)
  357. P^{(\alpha,\beta)}_{n+1}(x_i)}
  358. Parameters
  359. ==========
  360. n : the order of quadrature
  361. alpha : the first parameter of the Jacobi Polynomial, `\alpha > -1`
  362. beta : the second parameter of the Jacobi Polynomial, `\beta > -1`
  363. n_digits : number of significant digits of the points and weights to return
  364. Returns
  365. =======
  366. (x, w) : the ``x`` and ``w`` are lists of points and weights as Floats.
  367. The points `x_i` and weights `w_i` are returned as ``(x, w)``
  368. tuple of lists.
  369. Examples
  370. ========
  371. >>> from sympy import S
  372. >>> from sympy.integrals.quadrature import gauss_jacobi
  373. >>> x, w = gauss_jacobi(3, S.Half, -S.Half, 5)
  374. >>> x
  375. [-0.90097, -0.22252, 0.62349]
  376. >>> w
  377. [1.7063, 1.0973, 0.33795]
  378. >>> x, w = gauss_jacobi(6, 1, 1, 5)
  379. >>> x
  380. [-0.87174, -0.5917, -0.2093, 0.2093, 0.5917, 0.87174]
  381. >>> w
  382. [0.050584, 0.22169, 0.39439, 0.39439, 0.22169, 0.050584]
  383. See Also
  384. ========
  385. gauss_legendre, gauss_laguerre, gauss_hermite, gauss_gen_laguerre,
  386. gauss_chebyshev_t, gauss_chebyshev_u, gauss_lobatto
  387. References
  388. ==========
  389. .. [1] https://en.wikipedia.org/wiki/Gauss%E2%80%93Jacobi_quadrature
  390. .. [2] http://people.sc.fsu.edu/~jburkardt/cpp_src/jacobi_rule/jacobi_rule.html
  391. .. [3] http://people.sc.fsu.edu/~jburkardt/cpp_src/gegenbauer_rule/gegenbauer_rule.html
  392. """
  393. x = Dummy("x")
  394. p = jacobi_poly(n, alpha, beta, x, polys=True)
  395. pd = p.diff(x)
  396. pn = jacobi_poly(n+1, alpha, beta, x, polys=True)
  397. xi = []
  398. w = []
  399. for r in p.real_roots():
  400. if isinstance(r, RootOf):
  401. r = r.eval_rational(S.One/10**(n_digits+2))
  402. xi.append(r.n(n_digits))
  403. w.append((
  404. - (2*n+alpha+beta+2) / (n+alpha+beta+S.One) *
  405. (gamma(n+alpha+1)*gamma(n+beta+1)) /
  406. (gamma(n+alpha+beta+S.One)*gamma(n+2)) *
  407. 2**(alpha+beta) / (pd.subs(x, r) * pn.subs(x, r))).n(n_digits))
  408. return xi, w
  409. def gauss_lobatto(n, n_digits):
  410. r"""
  411. Computes the Gauss-Lobatto quadrature [1]_ points and weights.
  412. Explanation
  413. ===========
  414. The Gauss-Lobatto quadrature approximates the integral:
  415. .. math::
  416. \int_{-1}^1 f(x)\,dx \approx \sum_{i=1}^n w_i f(x_i)
  417. The nodes `x_i` of an order `n` quadrature rule are the roots of `P'_(n-1)`
  418. and the weights `w_i` are given by:
  419. .. math::
  420. &w_i = \frac{2}{n(n-1) \left[P_{n-1}(x_i)\right]^2},\quad x\neq\pm 1\\
  421. &w_i = \frac{2}{n(n-1)},\quad x=\pm 1
  422. Parameters
  423. ==========
  424. n : the order of quadrature
  425. n_digits : number of significant digits of the points and weights to return
  426. Returns
  427. =======
  428. (x, w) : the ``x`` and ``w`` are lists of points and weights as Floats.
  429. The points `x_i` and weights `w_i` are returned as ``(x, w)``
  430. tuple of lists.
  431. Examples
  432. ========
  433. >>> from sympy.integrals.quadrature import gauss_lobatto
  434. >>> x, w = gauss_lobatto(3, 5)
  435. >>> x
  436. [-1, 0, 1]
  437. >>> w
  438. [0.33333, 1.3333, 0.33333]
  439. >>> x, w = gauss_lobatto(4, 5)
  440. >>> x
  441. [-1, -0.44721, 0.44721, 1]
  442. >>> w
  443. [0.16667, 0.83333, 0.83333, 0.16667]
  444. See Also
  445. ========
  446. gauss_legendre,gauss_laguerre, gauss_gen_laguerre, gauss_hermite, gauss_chebyshev_t, gauss_chebyshev_u, gauss_jacobi
  447. References
  448. ==========
  449. .. [1] https://en.wikipedia.org/wiki/Gaussian_quadrature#Gauss.E2.80.93Lobatto_rules
  450. .. [2] http://people.math.sfu.ca/~cbm/aands/page_888.htm
  451. """
  452. x = Dummy("x")
  453. p = legendre_poly(n-1, x, polys=True)
  454. pd = p.diff(x)
  455. xi = []
  456. w = []
  457. for r in pd.real_roots():
  458. if isinstance(r, RootOf):
  459. r = r.eval_rational(S.One/10**(n_digits+2))
  460. xi.append(r.n(n_digits))
  461. w.append((2/(n*(n-1) * p.subs(x, r)**2)).n(n_digits))
  462. xi.insert(0, -1)
  463. xi.append(1)
  464. w.insert(0, (S(2)/(n*(n-1))).n(n_digits))
  465. w.append((S(2)/(n*(n-1))).n(n_digits))
  466. return xi, w