orthopolys.py 9.9 KB

123456789101112131415161718192021222324252627282930313233343536373839404142434445464748495051525354555657585960616263646566676869707172737475767778798081828384858687888990919293949596979899100101102103104105106107108109110111112113114115116117118119120121122123124125126127128129130131132133134135136137138139140141142143144145146147148149150151152153154155156157158159160161162163164165166167168169170171172173174175176177178179180181182183184185186187188189190191192193194195196197198199200201202203204205206207208209210211212213214215216217218219220221222223224225226227228229230231232233234235236237238239240241242243244245246247248249250251252253254255256257258259260261262263264265266267268269270271272273274275276277278279280281282283284285286287288289290291292293294295296297298299300301302303304305306307308309310311312313314315316317318319320321322323324325326327328329330331332333334335336337338339340341342343
  1. """Efficient functions for generating orthogonal polynomials."""
  2. from sympy.core.symbol import Dummy
  3. from sympy.polys.densearith import (dup_mul, dup_mul_ground,
  4. dup_lshift, dup_sub, dup_add, dup_sub_term, dup_sub_ground, dup_sqr)
  5. from sympy.polys.domains import ZZ, QQ
  6. from sympy.polys.polytools import named_poly
  7. from sympy.utilities import public
  8. def dup_jacobi(n, a, b, K):
  9. """Low-level implementation of Jacobi polynomials."""
  10. if n < 1:
  11. return [K.one]
  12. m2, m1 = [K.one], [(a+b)/K(2) + K.one, (a-b)/K(2)]
  13. for i in range(2, n+1):
  14. den = K(i)*(a + b + i)*(a + b + K(2)*i - K(2))
  15. f0 = (a + b + K(2)*i - K.one) * (a*a - b*b) / (K(2)*den)
  16. f1 = (a + b + K(2)*i - K.one) * (a + b + K(2)*i - K(2)) * (a + b + K(2)*i) / (K(2)*den)
  17. f2 = (a + i - K.one)*(b + i - K.one)*(a + b + K(2)*i) / den
  18. p0 = dup_mul_ground(m1, f0, K)
  19. p1 = dup_mul_ground(dup_lshift(m1, 1, K), f1, K)
  20. p2 = dup_mul_ground(m2, f2, K)
  21. m2, m1 = m1, dup_sub(dup_add(p0, p1, K), p2, K)
  22. return m1
  23. @public
  24. def jacobi_poly(n, a, b, x=None, polys=False):
  25. r"""Generates the Jacobi polynomial `P_n^{(a,b)}(x)`.
  26. Parameters
  27. ==========
  28. n : int
  29. Degree of the polynomial.
  30. a
  31. Lower limit of minimal domain for the list of coefficients.
  32. b
  33. Upper limit of minimal domain for the list of coefficients.
  34. x : optional
  35. polys : bool, optional
  36. If True, return a Poly, otherwise (default) return an expression.
  37. """
  38. return named_poly(n, dup_jacobi, None, "Jacobi polynomial", (x, a, b), polys)
  39. def dup_gegenbauer(n, a, K):
  40. """Low-level implementation of Gegenbauer polynomials."""
  41. if n < 1:
  42. return [K.one]
  43. m2, m1 = [K.one], [K(2)*a, K.zero]
  44. for i in range(2, n+1):
  45. p1 = dup_mul_ground(dup_lshift(m1, 1, K), K(2)*(a-K.one)/K(i) + K(2), K)
  46. p2 = dup_mul_ground(m2, K(2)*(a-K.one)/K(i) + K.one, K)
  47. m2, m1 = m1, dup_sub(p1, p2, K)
  48. return m1
  49. def gegenbauer_poly(n, a, x=None, polys=False):
  50. r"""Generates the Gegenbauer polynomial `C_n^{(a)}(x)`.
  51. Parameters
  52. ==========
  53. n : int
  54. Degree of the polynomial.
  55. x : optional
  56. a
  57. Decides minimal domain for the list of coefficients.
  58. polys : bool, optional
  59. If True, return a Poly, otherwise (default) return an expression.
  60. """
  61. return named_poly(n, dup_gegenbauer, None, "Gegenbauer polynomial", (x, a), polys)
  62. def dup_chebyshevt(n, K):
  63. """Low-level implementation of Chebyshev polynomials of the first kind."""
  64. if n < 1:
  65. return [K.one]
  66. # When n is small, it is faster to directly calculate the recurrence relation.
  67. if n < 64: # The threshold serves as a heuristic
  68. return _dup_chebyshevt_rec(n, K)
  69. return _dup_chebyshevt_prod(n, K)
  70. def _dup_chebyshevt_rec(n, K):
  71. r""" Chebyshev polynomials of the first kind using recurrence.
  72. Explanation
  73. ===========
  74. Chebyshev polynomials of the first kind are defined by the recurrence
  75. relation:
  76. .. math::
  77. T_0(x) &= 1\\
  78. T_1(x) &= x\\
  79. T_n(x) &= 2xT_{n-1}(x) - T_{n-2}(x)
  80. This function calculates the Chebyshev polynomial of the first kind using
  81. the above recurrence relation.
  82. Parameters
  83. ==========
  84. n : int
  85. n is a nonnegative integer.
  86. K : domain
  87. """
  88. m2, m1 = [K.one], [K.one, K.zero]
  89. for _ in range(n - 1):
  90. m2, m1 = m1, dup_sub(dup_mul_ground(dup_lshift(m1, 1, K), K(2), K), m2, K)
  91. return m1
  92. def _dup_chebyshevt_prod(n, K):
  93. r""" Chebyshev polynomials of the first kind using recursive products.
  94. Explanation
  95. ===========
  96. Computes Chebyshev polynomials of the first kind using
  97. .. math::
  98. T_{2n}(x) &= 2T_n^2(x) - 1\\
  99. T_{2n+1}(x) &= 2T_{n+1}(x)T_n(x) - x
  100. This is faster than ``_dup_chebyshevt_rec`` for large ``n``.
  101. Parameters
  102. ==========
  103. n : int
  104. n is a nonnegative integer.
  105. K : domain
  106. """
  107. m2, m1 = [K.one, K.zero], [K(2), K.zero, -K.one]
  108. for i in bin(n)[3:]:
  109. c = dup_sub_term(dup_mul_ground(dup_mul(m1, m2, K), K(2), K), K.one, 1, K)
  110. if i == '1':
  111. m2, m1 = c, dup_sub_ground(dup_mul_ground(dup_sqr(m1, K), K(2), K), K.one, K)
  112. else:
  113. m2, m1 = dup_sub_ground(dup_mul_ground(dup_sqr(m2, K), K(2), K), K.one, K), c
  114. return m2
  115. def dup_chebyshevu(n, K):
  116. """Low-level implementation of Chebyshev polynomials of the second kind."""
  117. if n < 1:
  118. return [K.one]
  119. m2, m1 = [K.one], [K(2), K.zero]
  120. for i in range(2, n+1):
  121. m2, m1 = m1, dup_sub(dup_mul_ground(dup_lshift(m1, 1, K), K(2), K), m2, K)
  122. return m1
  123. @public
  124. def chebyshevt_poly(n, x=None, polys=False):
  125. r"""Generates the Chebyshev polynomial of the first kind `T_n(x)`.
  126. Parameters
  127. ==========
  128. n : int
  129. Degree of the polynomial.
  130. x : optional
  131. polys : bool, optional
  132. If True, return a Poly, otherwise (default) return an expression.
  133. """
  134. return named_poly(n, dup_chebyshevt, ZZ,
  135. "Chebyshev polynomial of the first kind", (x,), polys)
  136. @public
  137. def chebyshevu_poly(n, x=None, polys=False):
  138. r"""Generates the Chebyshev polynomial of the second kind `U_n(x)`.
  139. Parameters
  140. ==========
  141. n : int
  142. Degree of the polynomial.
  143. x : optional
  144. polys : bool, optional
  145. If True, return a Poly, otherwise (default) return an expression.
  146. """
  147. return named_poly(n, dup_chebyshevu, ZZ,
  148. "Chebyshev polynomial of the second kind", (x,), polys)
  149. def dup_hermite(n, K):
  150. """Low-level implementation of Hermite polynomials."""
  151. if n < 1:
  152. return [K.one]
  153. m2, m1 = [K.one], [K(2), K.zero]
  154. for i in range(2, n+1):
  155. a = dup_lshift(m1, 1, K)
  156. b = dup_mul_ground(m2, K(i-1), K)
  157. m2, m1 = m1, dup_mul_ground(dup_sub(a, b, K), K(2), K)
  158. return m1
  159. def dup_hermite_prob(n, K):
  160. """Low-level implementation of probabilist's Hermite polynomials."""
  161. if n < 1:
  162. return [K.one]
  163. m2, m1 = [K.one], [K.one, K.zero]
  164. for i in range(2, n+1):
  165. a = dup_lshift(m1, 1, K)
  166. b = dup_mul_ground(m2, K(i-1), K)
  167. m2, m1 = m1, dup_sub(a, b, K)
  168. return m1
  169. @public
  170. def hermite_poly(n, x=None, polys=False):
  171. r"""Generates the Hermite polynomial `H_n(x)`.
  172. Parameters
  173. ==========
  174. n : int
  175. Degree of the polynomial.
  176. x : optional
  177. polys : bool, optional
  178. If True, return a Poly, otherwise (default) return an expression.
  179. """
  180. return named_poly(n, dup_hermite, ZZ, "Hermite polynomial", (x,), polys)
  181. @public
  182. def hermite_prob_poly(n, x=None, polys=False):
  183. r"""Generates the probabilist's Hermite polynomial `He_n(x)`.
  184. Parameters
  185. ==========
  186. n : int
  187. Degree of the polynomial.
  188. x : optional
  189. polys : bool, optional
  190. If True, return a Poly, otherwise (default) return an expression.
  191. """
  192. return named_poly(n, dup_hermite_prob, ZZ,
  193. "probabilist's Hermite polynomial", (x,), polys)
  194. def dup_legendre(n, K):
  195. """Low-level implementation of Legendre polynomials."""
  196. if n < 1:
  197. return [K.one]
  198. m2, m1 = [K.one], [K.one, K.zero]
  199. for i in range(2, n+1):
  200. a = dup_mul_ground(dup_lshift(m1, 1, K), K(2*i-1, i), K)
  201. b = dup_mul_ground(m2, K(i-1, i), K)
  202. m2, m1 = m1, dup_sub(a, b, K)
  203. return m1
  204. @public
  205. def legendre_poly(n, x=None, polys=False):
  206. r"""Generates the Legendre polynomial `P_n(x)`.
  207. Parameters
  208. ==========
  209. n : int
  210. Degree of the polynomial.
  211. x : optional
  212. polys : bool, optional
  213. If True, return a Poly, otherwise (default) return an expression.
  214. """
  215. return named_poly(n, dup_legendre, QQ, "Legendre polynomial", (x,), polys)
  216. def dup_laguerre(n, alpha, K):
  217. """Low-level implementation of Laguerre polynomials."""
  218. m2, m1 = [K.zero], [K.one]
  219. for i in range(1, n+1):
  220. a = dup_mul(m1, [-K.one/K(i), (alpha-K.one)/K(i) + K(2)], K)
  221. b = dup_mul_ground(m2, (alpha-K.one)/K(i) + K.one, K)
  222. m2, m1 = m1, dup_sub(a, b, K)
  223. return m1
  224. @public
  225. def laguerre_poly(n, x=None, alpha=0, polys=False):
  226. r"""Generates the Laguerre polynomial `L_n^{(\alpha)}(x)`.
  227. Parameters
  228. ==========
  229. n : int
  230. Degree of the polynomial.
  231. x : optional
  232. alpha : optional
  233. Decides minimal domain for the list of coefficients.
  234. polys : bool, optional
  235. If True, return a Poly, otherwise (default) return an expression.
  236. """
  237. return named_poly(n, dup_laguerre, None, "Laguerre polynomial", (x, alpha), polys)
  238. def dup_spherical_bessel_fn(n, K):
  239. """Low-level implementation of fn(n, x)."""
  240. if n < 1:
  241. return [K.one, K.zero]
  242. m2, m1 = [K.one], [K.one, K.zero]
  243. for i in range(2, n+1):
  244. m2, m1 = m1, dup_sub(dup_mul_ground(dup_lshift(m1, 1, K), K(2*i-1), K), m2, K)
  245. return dup_lshift(m1, 1, K)
  246. def dup_spherical_bessel_fn_minus(n, K):
  247. """Low-level implementation of fn(-n, x)."""
  248. m2, m1 = [K.one, K.zero], [K.zero]
  249. for i in range(2, n+1):
  250. m2, m1 = m1, dup_sub(dup_mul_ground(dup_lshift(m1, 1, K), K(3-2*i), K), m2, K)
  251. return m1
  252. def spherical_bessel_fn(n, x=None, polys=False):
  253. """
  254. Coefficients for the spherical Bessel functions.
  255. These are only needed in the jn() function.
  256. The coefficients are calculated from:
  257. fn(0, z) = 1/z
  258. fn(1, z) = 1/z**2
  259. fn(n-1, z) + fn(n+1, z) == (2*n+1)/z * fn(n, z)
  260. Parameters
  261. ==========
  262. n : int
  263. Degree of the polynomial.
  264. x : optional
  265. polys : bool, optional
  266. If True, return a Poly, otherwise (default) return an expression.
  267. Examples
  268. ========
  269. >>> from sympy.polys.orthopolys import spherical_bessel_fn as fn
  270. >>> from sympy import Symbol
  271. >>> z = Symbol("z")
  272. >>> fn(1, z)
  273. z**(-2)
  274. >>> fn(2, z)
  275. -1/z + 3/z**3
  276. >>> fn(3, z)
  277. -6/z**2 + 15/z**4
  278. >>> fn(4, z)
  279. 1/z - 45/z**3 + 105/z**5
  280. """
  281. if x is None:
  282. x = Dummy("x")
  283. f = dup_spherical_bessel_fn_minus if n < 0 else dup_spherical_bessel_fn
  284. return named_poly(abs(n), f, ZZ, "", (QQ(1)/x,), polys)