orthopolys.py 9.9 KB

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