bessel.py 62 KB

1234567891011121314151617181920212223242526272829303132333435363738394041424344454647484950515253545556575859606162636465666768697071727374757677787980818283848586878889909192939495969798991001011021031041051061071081091101111121131141151161171181191201211221231241251261271281291301311321331341351361371381391401411421431441451461471481491501511521531541551561571581591601611621631641651661671681691701711721731741751761771781791801811821831841851861871881891901911921931941951961971981992002012022032042052062072082092102112122132142152162172182192202212222232242252262272282292302312322332342352362372382392402412422432442452462472482492502512522532542552562572582592602612622632642652662672682692702712722732742752762772782792802812822832842852862872882892902912922932942952962972982993003013023033043053063073083093103113123133143153163173183193203213223233243253263273283293303313323333343353363373383393403413423433443453463473483493503513523533543553563573583593603613623633643653663673683693703713723733743753763773783793803813823833843853863873883893903913923933943953963973983994004014024034044054064074084094104114124134144154164174184194204214224234244254264274284294304314324334344354364374384394404414424434444454464474484494504514524534544554564574584594604614624634644654664674684694704714724734744754764774784794804814824834844854864874884894904914924934944954964974984995005015025035045055065075085095105115125135145155165175185195205215225235245255265275285295305315325335345355365375385395405415425435445455465475485495505515525535545555565575585595605615625635645655665675685695705715725735745755765775785795805815825835845855865875885895905915925935945955965975985996006016026036046056066076086096106116126136146156166176186196206216226236246256266276286296306316326336346356366376386396406416426436446456466476486496506516526536546556566576586596606616626636646656666676686696706716726736746756766776786796806816826836846856866876886896906916926936946956966976986997007017027037047057067077087097107117127137147157167177187197207217227237247257267277287297307317327337347357367377387397407417427437447457467477487497507517527537547557567577587597607617627637647657667677687697707717727737747757767777787797807817827837847857867877887897907917927937947957967977987998008018028038048058068078088098108118128138148158168178188198208218228238248258268278288298308318328338348358368378388398408418428438448458468478488498508518528538548558568578588598608618628638648658668678688698708718728738748758768778788798808818828838848858868878888898908918928938948958968978988999009019029039049059069079089099109119129139149159169179189199209219229239249259269279289299309319329339349359369379389399409419429439449459469479489499509519529539549559569579589599609619629639649659669679689699709719729739749759769779789799809819829839849859869879889899909919929939949959969979989991000100110021003100410051006100710081009101010111012101310141015101610171018101910201021102210231024102510261027102810291030103110321033103410351036103710381039104010411042104310441045104610471048104910501051105210531054105510561057105810591060106110621063106410651066106710681069107010711072107310741075107610771078107910801081108210831084108510861087108810891090109110921093109410951096109710981099110011011102110311041105110611071108110911101111111211131114111511161117111811191120112111221123112411251126112711281129113011311132113311341135113611371138113911401141114211431144114511461147114811491150115111521153115411551156115711581159116011611162116311641165116611671168116911701171117211731174117511761177117811791180118111821183118411851186118711881189119011911192119311941195119611971198119912001201120212031204120512061207120812091210121112121213121412151216121712181219122012211222122312241225122612271228122912301231123212331234123512361237123812391240124112421243124412451246124712481249125012511252125312541255125612571258125912601261126212631264126512661267126812691270127112721273127412751276127712781279128012811282128312841285128612871288128912901291129212931294129512961297129812991300130113021303130413051306130713081309131013111312131313141315131613171318131913201321132213231324132513261327132813291330133113321333133413351336133713381339134013411342134313441345134613471348134913501351135213531354135513561357135813591360136113621363136413651366136713681369137013711372137313741375137613771378137913801381138213831384138513861387138813891390139113921393139413951396139713981399140014011402140314041405140614071408140914101411141214131414141514161417141814191420142114221423142414251426142714281429143014311432143314341435143614371438143914401441144214431444144514461447144814491450145114521453145414551456145714581459146014611462146314641465146614671468146914701471147214731474147514761477147814791480148114821483148414851486148714881489149014911492149314941495149614971498149915001501150215031504150515061507150815091510151115121513151415151516151715181519152015211522152315241525152615271528152915301531153215331534153515361537153815391540154115421543154415451546154715481549155015511552155315541555155615571558155915601561156215631564156515661567156815691570157115721573157415751576157715781579158015811582158315841585158615871588158915901591159215931594159515961597159815991600160116021603160416051606160716081609161016111612161316141615161616171618161916201621162216231624162516261627162816291630163116321633163416351636163716381639164016411642164316441645164616471648164916501651165216531654165516561657165816591660166116621663166416651666166716681669167016711672167316741675167616771678167916801681168216831684168516861687168816891690169116921693169416951696169716981699170017011702170317041705170617071708170917101711171217131714171517161717171817191720172117221723172417251726172717281729173017311732173317341735173617371738173917401741174217431744174517461747174817491750175117521753175417551756175717581759176017611762176317641765176617671768176917701771177217731774177517761777177817791780178117821783178417851786178717881789179017911792179317941795179617971798179918001801180218031804180518061807180818091810181118121813181418151816181718181819182018211822182318241825182618271828182918301831183218331834183518361837183818391840184118421843184418451846184718481849185018511852185318541855185618571858185918601861186218631864186518661867186818691870187118721873187418751876187718781879188018811882188318841885188618871888188918901891189218931894189518961897189818991900190119021903190419051906190719081909191019111912191319141915191619171918191919201921192219231924192519261927192819291930193119321933193419351936193719381939194019411942194319441945194619471948194919501951195219531954195519561957195819591960196119621963196419651966196719681969197019711972197319741975197619771978197919801981198219831984198519861987198819891990199119921993199419951996199719981999200020012002200320042005200620072008200920102011201220132014201520162017201820192020202120222023202420252026202720282029203020312032203320342035203620372038203920402041204220432044204520462047204820492050205120522053205420552056205720582059206020612062206320642065206620672068206920702071207220732074207520762077207820792080208120822083208420852086208720882089
  1. from functools import wraps
  2. from sympy.core import S
  3. from sympy.core.add import Add
  4. from sympy.core.cache import cacheit
  5. from sympy.core.expr import Expr
  6. from sympy.core.function import Function, ArgumentIndexError, _mexpand
  7. from sympy.core.logic import fuzzy_or, fuzzy_not
  8. from sympy.core.numbers import Rational, pi, I
  9. from sympy.core.power import Pow
  10. from sympy.core.symbol import Dummy, Wild
  11. from sympy.core.sympify import sympify
  12. from sympy.functions.combinatorial.factorials import factorial
  13. from sympy.functions.elementary.trigonometric import sin, cos, csc, cot
  14. from sympy.functions.elementary.integers import ceiling
  15. from sympy.functions.elementary.exponential import exp, log
  16. from sympy.functions.elementary.miscellaneous import cbrt, sqrt, root
  17. from sympy.functions.elementary.complexes import (Abs, re, im, polar_lift, unpolarify)
  18. from sympy.functions.special.gamma_functions import gamma, digamma, uppergamma
  19. from sympy.functions.special.hyper import hyper
  20. from sympy.polys.orthopolys import spherical_bessel_fn
  21. from mpmath import mp, workprec
  22. # TODO
  23. # o Scorer functions G1 and G2
  24. # o Asymptotic expansions
  25. # These are possible, e.g. for fixed order, but since the bessel type
  26. # functions are oscillatory they are not actually tractable at
  27. # infinity, so this is not particularly useful right now.
  28. # o Nicer series expansions.
  29. # o More rewriting.
  30. # o Add solvers to ode.py (or rather add solvers for the hypergeometric equation).
  31. class BesselBase(Function):
  32. """
  33. Abstract base class for Bessel-type functions.
  34. This class is meant to reduce code duplication.
  35. All Bessel-type functions can 1) be differentiated, with the derivatives
  36. expressed in terms of similar functions, and 2) be rewritten in terms
  37. of other Bessel-type functions.
  38. Here, Bessel-type functions are assumed to have one complex parameter.
  39. To use this base class, define class attributes ``_a`` and ``_b`` such that
  40. ``2*F_n' = -_a*F_{n+1} + b*F_{n-1}``.
  41. """
  42. @property
  43. def order(self):
  44. """ The order of the Bessel-type function. """
  45. return self.args[0]
  46. @property
  47. def argument(self):
  48. """ The argument of the Bessel-type function. """
  49. return self.args[1]
  50. @classmethod
  51. def eval(cls, nu, z):
  52. return
  53. def fdiff(self, argindex=2):
  54. if argindex != 2:
  55. raise ArgumentIndexError(self, argindex)
  56. return (self._b/2 * self.__class__(self.order - 1, self.argument) -
  57. self._a/2 * self.__class__(self.order + 1, self.argument))
  58. def _eval_conjugate(self):
  59. z = self.argument
  60. if z.is_extended_negative is False:
  61. return self.__class__(self.order.conjugate(), z.conjugate())
  62. def _eval_is_meromorphic(self, x, a):
  63. nu, z = self.order, self.argument
  64. if nu.has(x):
  65. return False
  66. if not z._eval_is_meromorphic(x, a):
  67. return None
  68. z0 = z.subs(x, a)
  69. if nu.is_integer:
  70. if isinstance(self, (besselj, besseli, hn1, hn2, jn, yn)) or not nu.is_zero:
  71. return fuzzy_not(z0.is_infinite)
  72. return fuzzy_not(fuzzy_or([z0.is_zero, z0.is_infinite]))
  73. def _eval_expand_func(self, **hints):
  74. nu, z, f = self.order, self.argument, self.__class__
  75. if nu.is_real:
  76. if (nu - 1).is_positive:
  77. return (-self._a*self._b*f(nu - 2, z)._eval_expand_func() +
  78. 2*self._a*(nu - 1)*f(nu - 1, z)._eval_expand_func()/z)
  79. elif (nu + 1).is_negative:
  80. return (2*self._b*(nu + 1)*f(nu + 1, z)._eval_expand_func()/z -
  81. self._a*self._b*f(nu + 2, z)._eval_expand_func())
  82. return self
  83. def _eval_simplify(self, **kwargs):
  84. from sympy.simplify.simplify import besselsimp
  85. return besselsimp(self)
  86. class besselj(BesselBase):
  87. r"""
  88. Bessel function of the first kind.
  89. Explanation
  90. ===========
  91. The Bessel $J$ function of order $\nu$ is defined to be the function
  92. satisfying Bessel's differential equation
  93. .. math ::
  94. z^2 \frac{\mathrm{d}^2 w}{\mathrm{d}z^2}
  95. + z \frac{\mathrm{d}w}{\mathrm{d}z} + (z^2 - \nu^2) w = 0,
  96. with Laurent expansion
  97. .. math ::
  98. J_\nu(z) = z^\nu \left(\frac{1}{\Gamma(\nu + 1) 2^\nu} + O(z^2) \right),
  99. if $\nu$ is not a negative integer. If $\nu=-n \in \mathbb{Z}_{<0}$
  100. *is* a negative integer, then the definition is
  101. .. math ::
  102. J_{-n}(z) = (-1)^n J_n(z).
  103. Examples
  104. ========
  105. Create a Bessel function object:
  106. >>> from sympy import besselj, jn
  107. >>> from sympy.abc import z, n
  108. >>> b = besselj(n, z)
  109. Differentiate it:
  110. >>> b.diff(z)
  111. besselj(n - 1, z)/2 - besselj(n + 1, z)/2
  112. Rewrite in terms of spherical Bessel functions:
  113. >>> b.rewrite(jn)
  114. sqrt(2)*sqrt(z)*jn(n - 1/2, z)/sqrt(pi)
  115. Access the parameter and argument:
  116. >>> b.order
  117. n
  118. >>> b.argument
  119. z
  120. See Also
  121. ========
  122. bessely, besseli, besselk
  123. References
  124. ==========
  125. .. [1] Abramowitz, Milton; Stegun, Irene A., eds. (1965), "Chapter 9",
  126. Handbook of Mathematical Functions with Formulas, Graphs, and
  127. Mathematical Tables
  128. .. [2] Luke, Y. L. (1969), The Special Functions and Their
  129. Approximations, Volume 1
  130. .. [3] https://en.wikipedia.org/wiki/Bessel_function
  131. .. [4] https://functions.wolfram.com/Bessel-TypeFunctions/BesselJ/
  132. """
  133. _a = S.One
  134. _b = S.One
  135. @classmethod
  136. def eval(cls, nu, z):
  137. if z.is_zero:
  138. if nu.is_zero:
  139. return S.One
  140. elif (nu.is_integer and nu.is_zero is False) or re(nu).is_positive:
  141. return S.Zero
  142. elif re(nu).is_negative and not (nu.is_integer is True):
  143. return S.ComplexInfinity
  144. elif nu.is_imaginary:
  145. return S.NaN
  146. if z in (S.Infinity, S.NegativeInfinity):
  147. return S.Zero
  148. if z.could_extract_minus_sign():
  149. return (z)**nu*(-z)**(-nu)*besselj(nu, -z)
  150. if nu.is_integer:
  151. if nu.could_extract_minus_sign():
  152. return S.NegativeOne**(-nu)*besselj(-nu, z)
  153. newz = z.extract_multiplicatively(I)
  154. if newz: # NOTE we don't want to change the function if z==0
  155. return I**(nu)*besseli(nu, newz)
  156. # branch handling:
  157. if nu.is_integer:
  158. newz = unpolarify(z)
  159. if newz != z:
  160. return besselj(nu, newz)
  161. else:
  162. newz, n = z.extract_branch_factor()
  163. if n != 0:
  164. return exp(2*n*pi*nu*I)*besselj(nu, newz)
  165. nnu = unpolarify(nu)
  166. if nu != nnu:
  167. return besselj(nnu, z)
  168. def _eval_rewrite_as_besseli(self, nu, z, **kwargs):
  169. return exp(I*pi*nu/2)*besseli(nu, polar_lift(-I)*z)
  170. def _eval_rewrite_as_bessely(self, nu, z, **kwargs):
  171. if nu.is_integer is False:
  172. return csc(pi*nu)*bessely(-nu, z) - cot(pi*nu)*bessely(nu, z)
  173. def _eval_rewrite_as_jn(self, nu, z, **kwargs):
  174. return sqrt(2*z/pi)*jn(nu - S.Half, self.argument)
  175. def _eval_as_leading_term(self, x, logx=None, cdir=0):
  176. nu, z = self.args
  177. try:
  178. arg = z.as_leading_term(x)
  179. except NotImplementedError:
  180. return self
  181. c, e = arg.as_coeff_exponent(x)
  182. if e.is_positive:
  183. return arg**nu/(2**nu*gamma(nu + 1))
  184. elif e.is_negative:
  185. cdir = 1 if cdir == 0 else cdir
  186. sign = c*cdir**e
  187. if not sign.is_negative:
  188. # Refer Abramowitz and Stegun 1965, p. 364 for more information on
  189. # asymptotic approximation of besselj function.
  190. return sqrt(2)*cos(z - pi*(2*nu + 1)/4)/sqrt(pi*z)
  191. return self
  192. return super(besselj, self)._eval_as_leading_term(x, logx, cdir)
  193. def _eval_is_extended_real(self):
  194. nu, z = self.args
  195. if nu.is_integer and z.is_extended_real:
  196. return True
  197. def _eval_nseries(self, x, n, logx, cdir=0):
  198. # Refer https://functions.wolfram.com/Bessel-TypeFunctions/BesselJ/06/01/04/01/01/0003/
  199. # for more information on nseries expansion of besselj function.
  200. from sympy.series.order import Order
  201. nu, z = self.args
  202. # In case of powers less than 1, number of terms need to be computed
  203. # separately to avoid repeated callings of _eval_nseries with wrong n
  204. try:
  205. _, exp = z.leadterm(x)
  206. except (ValueError, NotImplementedError):
  207. return self
  208. if exp.is_positive:
  209. newn = ceiling(n/exp)
  210. o = Order(x**n, x)
  211. r = (z/2)._eval_nseries(x, n, logx, cdir).removeO()
  212. if r is S.Zero:
  213. return o
  214. t = (_mexpand(r**2) + o).removeO()
  215. term = r**nu/gamma(nu + 1)
  216. s = [term]
  217. for k in range(1, (newn + 1)//2):
  218. term *= -t/(k*(nu + k))
  219. term = (_mexpand(term) + o).removeO()
  220. s.append(term)
  221. return Add(*s) + o
  222. return super(besselj, self)._eval_nseries(x, n, logx, cdir)
  223. class bessely(BesselBase):
  224. r"""
  225. Bessel function of the second kind.
  226. Explanation
  227. ===========
  228. The Bessel $Y$ function of order $\nu$ is defined as
  229. .. math ::
  230. Y_\nu(z) = \lim_{\mu \to \nu} \frac{J_\mu(z) \cos(\pi \mu)
  231. - J_{-\mu}(z)}{\sin(\pi \mu)},
  232. where $J_\mu(z)$ is the Bessel function of the first kind.
  233. It is a solution to Bessel's equation, and linearly independent from
  234. $J_\nu$.
  235. Examples
  236. ========
  237. >>> from sympy import bessely, yn
  238. >>> from sympy.abc import z, n
  239. >>> b = bessely(n, z)
  240. >>> b.diff(z)
  241. bessely(n - 1, z)/2 - bessely(n + 1, z)/2
  242. >>> b.rewrite(yn)
  243. sqrt(2)*sqrt(z)*yn(n - 1/2, z)/sqrt(pi)
  244. See Also
  245. ========
  246. besselj, besseli, besselk
  247. References
  248. ==========
  249. .. [1] https://functions.wolfram.com/Bessel-TypeFunctions/BesselY/
  250. """
  251. _a = S.One
  252. _b = S.One
  253. @classmethod
  254. def eval(cls, nu, z):
  255. if z.is_zero:
  256. if nu.is_zero:
  257. return S.NegativeInfinity
  258. elif re(nu).is_zero is False:
  259. return S.ComplexInfinity
  260. elif re(nu).is_zero:
  261. return S.NaN
  262. if z in (S.Infinity, S.NegativeInfinity):
  263. return S.Zero
  264. if z == I*S.Infinity:
  265. return exp(I*pi*(nu + 1)/2) * S.Infinity
  266. if z == I*S.NegativeInfinity:
  267. return exp(-I*pi*(nu + 1)/2) * S.Infinity
  268. if nu.is_integer:
  269. if nu.could_extract_minus_sign():
  270. return S.NegativeOne**(-nu)*bessely(-nu, z)
  271. def _eval_rewrite_as_besselj(self, nu, z, **kwargs):
  272. if nu.is_integer is False:
  273. return csc(pi*nu)*(cos(pi*nu)*besselj(nu, z) - besselj(-nu, z))
  274. def _eval_rewrite_as_besseli(self, nu, z, **kwargs):
  275. aj = self._eval_rewrite_as_besselj(*self.args)
  276. if aj:
  277. return aj.rewrite(besseli)
  278. def _eval_rewrite_as_yn(self, nu, z, **kwargs):
  279. return sqrt(2*z/pi) * yn(nu - S.Half, self.argument)
  280. def _eval_as_leading_term(self, x, logx=None, cdir=0):
  281. nu, z = self.args
  282. try:
  283. arg = z.as_leading_term(x)
  284. except NotImplementedError:
  285. return self
  286. c, e = arg.as_coeff_exponent(x)
  287. if e.is_positive:
  288. term_one = ((2/pi)*log(z/2)*besselj(nu, z))
  289. term_two = -(z/2)**(-nu)*factorial(nu - 1)/pi if (nu).is_positive else S.Zero
  290. term_three = -(z/2)**nu/(pi*factorial(nu))*(digamma(nu + 1) - S.EulerGamma)
  291. arg = Add(*[term_one, term_two, term_three]).as_leading_term(x, logx=logx)
  292. return arg
  293. elif e.is_negative:
  294. cdir = 1 if cdir == 0 else cdir
  295. sign = c*cdir**e
  296. if not sign.is_negative:
  297. # Refer Abramowitz and Stegun 1965, p. 364 for more information on
  298. # asymptotic approximation of bessely function.
  299. return sqrt(2)*(-sin(pi*nu/2 - z + pi/4) + 3*cos(pi*nu/2 - z + pi/4)/(8*z))*sqrt(1/z)/sqrt(pi)
  300. return self
  301. return super(bessely, self)._eval_as_leading_term(x, logx, cdir)
  302. def _eval_is_extended_real(self):
  303. nu, z = self.args
  304. if nu.is_integer and z.is_positive:
  305. return True
  306. def _eval_nseries(self, x, n, logx, cdir=0):
  307. # Refer https://functions.wolfram.com/Bessel-TypeFunctions/BesselY/06/01/04/01/02/0008/
  308. # for more information on nseries expansion of bessely function.
  309. from sympy.series.order import Order
  310. nu, z = self.args
  311. # In case of powers less than 1, number of terms need to be computed
  312. # separately to avoid repeated callings of _eval_nseries with wrong n
  313. try:
  314. _, exp = z.leadterm(x)
  315. except (ValueError, NotImplementedError):
  316. return self
  317. if exp.is_positive and nu.is_integer:
  318. newn = ceiling(n/exp)
  319. bn = besselj(nu, z)
  320. a = ((2/pi)*log(z/2)*bn)._eval_nseries(x, n, logx, cdir)
  321. b, c = [], []
  322. o = Order(x**n, x)
  323. r = (z/2)._eval_nseries(x, n, logx, cdir).removeO()
  324. if r is S.Zero:
  325. return o
  326. t = (_mexpand(r**2) + o).removeO()
  327. if nu > S.Zero:
  328. term = r**(-nu)*factorial(nu - 1)/pi
  329. b.append(term)
  330. for k in range(1, nu):
  331. denom = (nu - k)*k
  332. if denom == S.Zero:
  333. term *= t/k
  334. else:
  335. term *= t/denom
  336. term = (_mexpand(term) + o).removeO()
  337. b.append(term)
  338. p = r**nu/(pi*factorial(nu))
  339. term = p*(digamma(nu + 1) - S.EulerGamma)
  340. c.append(term)
  341. for k in range(1, (newn + 1)//2):
  342. p *= -t/(k*(k + nu))
  343. p = (_mexpand(p) + o).removeO()
  344. term = p*(digamma(k + nu + 1) + digamma(k + 1))
  345. c.append(term)
  346. return a - Add(*b) - Add(*c) # Order term comes from a
  347. return super(bessely, self)._eval_nseries(x, n, logx, cdir)
  348. class besseli(BesselBase):
  349. r"""
  350. Modified Bessel function of the first kind.
  351. Explanation
  352. ===========
  353. The Bessel $I$ function is a solution to the modified Bessel equation
  354. .. math ::
  355. z^2 \frac{\mathrm{d}^2 w}{\mathrm{d}z^2}
  356. + z \frac{\mathrm{d}w}{\mathrm{d}z} + (z^2 + \nu^2)^2 w = 0.
  357. It can be defined as
  358. .. math ::
  359. I_\nu(z) = i^{-\nu} J_\nu(iz),
  360. where $J_\nu(z)$ is the Bessel function of the first kind.
  361. Examples
  362. ========
  363. >>> from sympy import besseli
  364. >>> from sympy.abc import z, n
  365. >>> besseli(n, z).diff(z)
  366. besseli(n - 1, z)/2 + besseli(n + 1, z)/2
  367. See Also
  368. ========
  369. besselj, bessely, besselk
  370. References
  371. ==========
  372. .. [1] https://functions.wolfram.com/Bessel-TypeFunctions/BesselI/
  373. """
  374. _a = -S.One
  375. _b = S.One
  376. @classmethod
  377. def eval(cls, nu, z):
  378. if z.is_zero:
  379. if nu.is_zero:
  380. return S.One
  381. elif (nu.is_integer and nu.is_zero is False) or re(nu).is_positive:
  382. return S.Zero
  383. elif re(nu).is_negative and not (nu.is_integer is True):
  384. return S.ComplexInfinity
  385. elif nu.is_imaginary:
  386. return S.NaN
  387. if im(z) in (S.Infinity, S.NegativeInfinity):
  388. return S.Zero
  389. if z is S.Infinity:
  390. return S.Infinity
  391. if z is S.NegativeInfinity:
  392. return (-1)**nu*S.Infinity
  393. if z.could_extract_minus_sign():
  394. return (z)**nu*(-z)**(-nu)*besseli(nu, -z)
  395. if nu.is_integer:
  396. if nu.could_extract_minus_sign():
  397. return besseli(-nu, z)
  398. newz = z.extract_multiplicatively(I)
  399. if newz: # NOTE we don't want to change the function if z==0
  400. return I**(-nu)*besselj(nu, -newz)
  401. # branch handling:
  402. if nu.is_integer:
  403. newz = unpolarify(z)
  404. if newz != z:
  405. return besseli(nu, newz)
  406. else:
  407. newz, n = z.extract_branch_factor()
  408. if n != 0:
  409. return exp(2*n*pi*nu*I)*besseli(nu, newz)
  410. nnu = unpolarify(nu)
  411. if nu != nnu:
  412. return besseli(nnu, z)
  413. def _eval_rewrite_as_besselj(self, nu, z, **kwargs):
  414. return exp(-I*pi*nu/2)*besselj(nu, polar_lift(I)*z)
  415. def _eval_rewrite_as_bessely(self, nu, z, **kwargs):
  416. aj = self._eval_rewrite_as_besselj(*self.args)
  417. if aj:
  418. return aj.rewrite(bessely)
  419. def _eval_rewrite_as_jn(self, nu, z, **kwargs):
  420. return self._eval_rewrite_as_besselj(*self.args).rewrite(jn)
  421. def _eval_is_extended_real(self):
  422. nu, z = self.args
  423. if nu.is_integer and z.is_extended_real:
  424. return True
  425. def _eval_as_leading_term(self, x, logx=None, cdir=0):
  426. nu, z = self.args
  427. try:
  428. arg = z.as_leading_term(x)
  429. except NotImplementedError:
  430. return self
  431. c, e = arg.as_coeff_exponent(x)
  432. if e.is_positive:
  433. return arg**nu/(2**nu*gamma(nu + 1))
  434. elif e.is_negative:
  435. cdir = 1 if cdir == 0 else cdir
  436. sign = c*cdir**e
  437. if not sign.is_negative:
  438. # Refer Abramowitz and Stegun 1965, p. 377 for more information on
  439. # asymptotic approximation of besseli function.
  440. return exp(z)/sqrt(2*pi*z)
  441. return self
  442. return super(besseli, self)._eval_as_leading_term(x, logx, cdir)
  443. def _eval_nseries(self, x, n, logx, cdir=0):
  444. # Refer https://functions.wolfram.com/Bessel-TypeFunctions/BesselI/06/01/04/01/01/0003/
  445. # for more information on nseries expansion of besseli function.
  446. from sympy.series.order import Order
  447. nu, z = self.args
  448. # In case of powers less than 1, number of terms need to be computed
  449. # separately to avoid repeated callings of _eval_nseries with wrong n
  450. try:
  451. _, exp = z.leadterm(x)
  452. except (ValueError, NotImplementedError):
  453. return self
  454. if exp.is_positive:
  455. newn = ceiling(n/exp)
  456. o = Order(x**n, x)
  457. r = (z/2)._eval_nseries(x, n, logx, cdir).removeO()
  458. if r is S.Zero:
  459. return o
  460. t = (_mexpand(r**2) + o).removeO()
  461. term = r**nu/gamma(nu + 1)
  462. s = [term]
  463. for k in range(1, (newn + 1)//2):
  464. term *= t/(k*(nu + k))
  465. term = (_mexpand(term) + o).removeO()
  466. s.append(term)
  467. return Add(*s) + o
  468. return super(besseli, self)._eval_nseries(x, n, logx, cdir)
  469. class besselk(BesselBase):
  470. r"""
  471. Modified Bessel function of the second kind.
  472. Explanation
  473. ===========
  474. The Bessel $K$ function of order $\nu$ is defined as
  475. .. math ::
  476. K_\nu(z) = \lim_{\mu \to \nu} \frac{\pi}{2}
  477. \frac{I_{-\mu}(z) -I_\mu(z)}{\sin(\pi \mu)},
  478. where $I_\mu(z)$ is the modified Bessel function of the first kind.
  479. It is a solution of the modified Bessel equation, and linearly independent
  480. from $Y_\nu$.
  481. Examples
  482. ========
  483. >>> from sympy import besselk
  484. >>> from sympy.abc import z, n
  485. >>> besselk(n, z).diff(z)
  486. -besselk(n - 1, z)/2 - besselk(n + 1, z)/2
  487. See Also
  488. ========
  489. besselj, besseli, bessely
  490. References
  491. ==========
  492. .. [1] https://functions.wolfram.com/Bessel-TypeFunctions/BesselK/
  493. """
  494. _a = S.One
  495. _b = -S.One
  496. @classmethod
  497. def eval(cls, nu, z):
  498. if z.is_zero:
  499. if nu.is_zero:
  500. return S.Infinity
  501. elif re(nu).is_zero is False:
  502. return S.ComplexInfinity
  503. elif re(nu).is_zero:
  504. return S.NaN
  505. if z in (S.Infinity, I*S.Infinity, I*S.NegativeInfinity):
  506. return S.Zero
  507. if nu.is_integer:
  508. if nu.could_extract_minus_sign():
  509. return besselk(-nu, z)
  510. def _eval_rewrite_as_besseli(self, nu, z, **kwargs):
  511. if nu.is_integer is False:
  512. return pi*csc(pi*nu)*(besseli(-nu, z) - besseli(nu, z))/2
  513. def _eval_rewrite_as_besselj(self, nu, z, **kwargs):
  514. ai = self._eval_rewrite_as_besseli(*self.args)
  515. if ai:
  516. return ai.rewrite(besselj)
  517. def _eval_rewrite_as_bessely(self, nu, z, **kwargs):
  518. aj = self._eval_rewrite_as_besselj(*self.args)
  519. if aj:
  520. return aj.rewrite(bessely)
  521. def _eval_rewrite_as_yn(self, nu, z, **kwargs):
  522. ay = self._eval_rewrite_as_bessely(*self.args)
  523. if ay:
  524. return ay.rewrite(yn)
  525. def _eval_is_extended_real(self):
  526. nu, z = self.args
  527. if nu.is_integer and z.is_positive:
  528. return True
  529. def _eval_as_leading_term(self, x, logx=None, cdir=0):
  530. nu, z = self.args
  531. try:
  532. arg = z.as_leading_term(x)
  533. except NotImplementedError:
  534. return self
  535. _, e = arg.as_coeff_exponent(x)
  536. if e.is_positive:
  537. term_one = ((-1)**(nu -1)*log(z/2)*besseli(nu, z))
  538. term_two = (z/2)**(-nu)*factorial(nu - 1)/2 if (nu).is_positive else S.Zero
  539. term_three = (-1)**nu*(z/2)**nu/(2*factorial(nu))*(digamma(nu + 1) - S.EulerGamma)
  540. arg = Add(*[term_one, term_two, term_three]).as_leading_term(x, logx=logx)
  541. return arg
  542. elif e.is_negative:
  543. # Refer Abramowitz and Stegun 1965, p. 378 for more information on
  544. # asymptotic approximation of besselk function.
  545. return sqrt(pi)*exp(-z)/sqrt(2*z)
  546. return super(besselk, self)._eval_as_leading_term(x, logx, cdir)
  547. def _eval_nseries(self, x, n, logx, cdir=0):
  548. # Refer https://functions.wolfram.com/Bessel-TypeFunctions/BesselK/06/01/04/01/02/0008/
  549. # for more information on nseries expansion of besselk function.
  550. from sympy.series.order import Order
  551. nu, z = self.args
  552. # In case of powers less than 1, number of terms need to be computed
  553. # separately to avoid repeated callings of _eval_nseries with wrong n
  554. try:
  555. _, exp = z.leadterm(x)
  556. except (ValueError, NotImplementedError):
  557. return self
  558. if exp.is_positive and nu.is_integer:
  559. newn = ceiling(n/exp)
  560. bn = besseli(nu, z)
  561. a = ((-1)**(nu - 1)*log(z/2)*bn)._eval_nseries(x, n, logx, cdir)
  562. b, c = [], []
  563. o = Order(x**n, x)
  564. r = (z/2)._eval_nseries(x, n, logx, cdir).removeO()
  565. if r is S.Zero:
  566. return o
  567. t = (_mexpand(r**2) + o).removeO()
  568. if nu > S.Zero:
  569. term = r**(-nu)*factorial(nu - 1)/2
  570. b.append(term)
  571. for k in range(1, nu):
  572. denom = (k - nu)*k
  573. if denom == S.Zero:
  574. term *= t/k
  575. else:
  576. term *= t/denom
  577. term = (_mexpand(term) + o).removeO()
  578. b.append(term)
  579. p = r**nu*(-1)**nu/(2*factorial(nu))
  580. term = p*(digamma(nu + 1) - S.EulerGamma)
  581. c.append(term)
  582. for k in range(1, (newn + 1)//2):
  583. p *= t/(k*(k + nu))
  584. p = (_mexpand(p) + o).removeO()
  585. term = p*(digamma(k + nu + 1) + digamma(k + 1))
  586. c.append(term)
  587. return a + Add(*b) + Add(*c) # Order term comes from a
  588. return super(besselk, self)._eval_nseries(x, n, logx, cdir)
  589. class hankel1(BesselBase):
  590. r"""
  591. Hankel function of the first kind.
  592. Explanation
  593. ===========
  594. This function is defined as
  595. .. math ::
  596. H_\nu^{(1)} = J_\nu(z) + iY_\nu(z),
  597. where $J_\nu(z)$ is the Bessel function of the first kind, and
  598. $Y_\nu(z)$ is the Bessel function of the second kind.
  599. It is a solution to Bessel's equation.
  600. Examples
  601. ========
  602. >>> from sympy import hankel1
  603. >>> from sympy.abc import z, n
  604. >>> hankel1(n, z).diff(z)
  605. hankel1(n - 1, z)/2 - hankel1(n + 1, z)/2
  606. See Also
  607. ========
  608. hankel2, besselj, bessely
  609. References
  610. ==========
  611. .. [1] https://functions.wolfram.com/Bessel-TypeFunctions/HankelH1/
  612. """
  613. _a = S.One
  614. _b = S.One
  615. def _eval_conjugate(self):
  616. z = self.argument
  617. if z.is_extended_negative is False:
  618. return hankel2(self.order.conjugate(), z.conjugate())
  619. class hankel2(BesselBase):
  620. r"""
  621. Hankel function of the second kind.
  622. Explanation
  623. ===========
  624. This function is defined as
  625. .. math ::
  626. H_\nu^{(2)} = J_\nu(z) - iY_\nu(z),
  627. where $J_\nu(z)$ is the Bessel function of the first kind, and
  628. $Y_\nu(z)$ is the Bessel function of the second kind.
  629. It is a solution to Bessel's equation, and linearly independent from
  630. $H_\nu^{(1)}$.
  631. Examples
  632. ========
  633. >>> from sympy import hankel2
  634. >>> from sympy.abc import z, n
  635. >>> hankel2(n, z).diff(z)
  636. hankel2(n - 1, z)/2 - hankel2(n + 1, z)/2
  637. See Also
  638. ========
  639. hankel1, besselj, bessely
  640. References
  641. ==========
  642. .. [1] https://functions.wolfram.com/Bessel-TypeFunctions/HankelH2/
  643. """
  644. _a = S.One
  645. _b = S.One
  646. def _eval_conjugate(self):
  647. z = self.argument
  648. if z.is_extended_negative is False:
  649. return hankel1(self.order.conjugate(), z.conjugate())
  650. def assume_integer_order(fn):
  651. @wraps(fn)
  652. def g(self, nu, z):
  653. if nu.is_integer:
  654. return fn(self, nu, z)
  655. return g
  656. class SphericalBesselBase(BesselBase):
  657. """
  658. Base class for spherical Bessel functions.
  659. These are thin wrappers around ordinary Bessel functions,
  660. since spherical Bessel functions differ from the ordinary
  661. ones just by a slight change in order.
  662. To use this class, define the ``_eval_evalf()`` and ``_expand()`` methods.
  663. """
  664. def _expand(self, **hints):
  665. """ Expand self into a polynomial. Nu is guaranteed to be Integer. """
  666. raise NotImplementedError('expansion')
  667. def _eval_expand_func(self, **hints):
  668. if self.order.is_Integer:
  669. return self._expand(**hints)
  670. return self
  671. def fdiff(self, argindex=2):
  672. if argindex != 2:
  673. raise ArgumentIndexError(self, argindex)
  674. return self.__class__(self.order - 1, self.argument) - \
  675. self * (self.order + 1)/self.argument
  676. def _jn(n, z):
  677. return (spherical_bessel_fn(n, z)*sin(z) +
  678. S.NegativeOne**(n + 1)*spherical_bessel_fn(-n - 1, z)*cos(z))
  679. def _yn(n, z):
  680. # (-1)**(n + 1) * _jn(-n - 1, z)
  681. return (S.NegativeOne**(n + 1) * spherical_bessel_fn(-n - 1, z)*sin(z) -
  682. spherical_bessel_fn(n, z)*cos(z))
  683. class jn(SphericalBesselBase):
  684. r"""
  685. Spherical Bessel function of the first kind.
  686. Explanation
  687. ===========
  688. This function is a solution to the spherical Bessel equation
  689. .. math ::
  690. z^2 \frac{\mathrm{d}^2 w}{\mathrm{d}z^2}
  691. + 2z \frac{\mathrm{d}w}{\mathrm{d}z} + (z^2 - \nu(\nu + 1)) w = 0.
  692. It can be defined as
  693. .. math ::
  694. j_\nu(z) = \sqrt{\frac{\pi}{2z}} J_{\nu + \frac{1}{2}}(z),
  695. where $J_\nu(z)$ is the Bessel function of the first kind.
  696. The spherical Bessel functions of integral order are
  697. calculated using the formula:
  698. .. math:: j_n(z) = f_n(z) \sin{z} + (-1)^{n+1} f_{-n-1}(z) \cos{z},
  699. where the coefficients $f_n(z)$ are available as
  700. :func:`sympy.polys.orthopolys.spherical_bessel_fn`.
  701. Examples
  702. ========
  703. >>> from sympy import Symbol, jn, sin, cos, expand_func, besselj, bessely
  704. >>> z = Symbol("z")
  705. >>> nu = Symbol("nu", integer=True)
  706. >>> print(expand_func(jn(0, z)))
  707. sin(z)/z
  708. >>> expand_func(jn(1, z)) == sin(z)/z**2 - cos(z)/z
  709. True
  710. >>> expand_func(jn(3, z))
  711. (-6/z**2 + 15/z**4)*sin(z) + (1/z - 15/z**3)*cos(z)
  712. >>> jn(nu, z).rewrite(besselj)
  713. sqrt(2)*sqrt(pi)*sqrt(1/z)*besselj(nu + 1/2, z)/2
  714. >>> jn(nu, z).rewrite(bessely)
  715. (-1)**nu*sqrt(2)*sqrt(pi)*sqrt(1/z)*bessely(-nu - 1/2, z)/2
  716. >>> jn(2, 5.2+0.3j).evalf(20)
  717. 0.099419756723640344491 - 0.054525080242173562897*I
  718. See Also
  719. ========
  720. besselj, bessely, besselk, yn
  721. References
  722. ==========
  723. .. [1] https://dlmf.nist.gov/10.47
  724. """
  725. @classmethod
  726. def eval(cls, nu, z):
  727. if z.is_zero:
  728. if nu.is_zero:
  729. return S.One
  730. elif nu.is_integer:
  731. if nu.is_positive:
  732. return S.Zero
  733. else:
  734. return S.ComplexInfinity
  735. if z in (S.NegativeInfinity, S.Infinity):
  736. return S.Zero
  737. def _eval_rewrite_as_besselj(self, nu, z, **kwargs):
  738. return sqrt(pi/(2*z)) * besselj(nu + S.Half, z)
  739. def _eval_rewrite_as_bessely(self, nu, z, **kwargs):
  740. return S.NegativeOne**nu * sqrt(pi/(2*z)) * bessely(-nu - S.Half, z)
  741. def _eval_rewrite_as_yn(self, nu, z, **kwargs):
  742. return S.NegativeOne**(nu) * yn(-nu - 1, z)
  743. def _expand(self, **hints):
  744. return _jn(self.order, self.argument)
  745. def _eval_evalf(self, prec):
  746. if self.order.is_Integer:
  747. return self.rewrite(besselj)._eval_evalf(prec)
  748. class yn(SphericalBesselBase):
  749. r"""
  750. Spherical Bessel function of the second kind.
  751. Explanation
  752. ===========
  753. This function is another solution to the spherical Bessel equation, and
  754. linearly independent from $j_n$. It can be defined as
  755. .. math ::
  756. y_\nu(z) = \sqrt{\frac{\pi}{2z}} Y_{\nu + \frac{1}{2}}(z),
  757. where $Y_\nu(z)$ is the Bessel function of the second kind.
  758. For integral orders $n$, $y_n$ is calculated using the formula:
  759. .. math:: y_n(z) = (-1)^{n+1} j_{-n-1}(z)
  760. Examples
  761. ========
  762. >>> from sympy import Symbol, yn, sin, cos, expand_func, besselj, bessely
  763. >>> z = Symbol("z")
  764. >>> nu = Symbol("nu", integer=True)
  765. >>> print(expand_func(yn(0, z)))
  766. -cos(z)/z
  767. >>> expand_func(yn(1, z)) == -cos(z)/z**2-sin(z)/z
  768. True
  769. >>> yn(nu, z).rewrite(besselj)
  770. (-1)**(nu + 1)*sqrt(2)*sqrt(pi)*sqrt(1/z)*besselj(-nu - 1/2, z)/2
  771. >>> yn(nu, z).rewrite(bessely)
  772. sqrt(2)*sqrt(pi)*sqrt(1/z)*bessely(nu + 1/2, z)/2
  773. >>> yn(2, 5.2+0.3j).evalf(20)
  774. 0.18525034196069722536 + 0.014895573969924817587*I
  775. See Also
  776. ========
  777. besselj, bessely, besselk, jn
  778. References
  779. ==========
  780. .. [1] https://dlmf.nist.gov/10.47
  781. """
  782. @assume_integer_order
  783. def _eval_rewrite_as_besselj(self, nu, z, **kwargs):
  784. return S.NegativeOne**(nu+1) * sqrt(pi/(2*z)) * besselj(-nu - S.Half, z)
  785. @assume_integer_order
  786. def _eval_rewrite_as_bessely(self, nu, z, **kwargs):
  787. return sqrt(pi/(2*z)) * bessely(nu + S.Half, z)
  788. def _eval_rewrite_as_jn(self, nu, z, **kwargs):
  789. return S.NegativeOne**(nu + 1) * jn(-nu - 1, z)
  790. def _expand(self, **hints):
  791. return _yn(self.order, self.argument)
  792. def _eval_evalf(self, prec):
  793. if self.order.is_Integer:
  794. return self.rewrite(bessely)._eval_evalf(prec)
  795. class SphericalHankelBase(SphericalBesselBase):
  796. @assume_integer_order
  797. def _eval_rewrite_as_besselj(self, nu, z, **kwargs):
  798. # jn +- I*yn
  799. # jn as beeselj: sqrt(pi/(2*z)) * besselj(nu + S.Half, z)
  800. # yn as besselj: (-1)**(nu+1) * sqrt(pi/(2*z)) * besselj(-nu - S.Half, z)
  801. hks = self._hankel_kind_sign
  802. return sqrt(pi/(2*z))*(besselj(nu + S.Half, z) +
  803. hks*I*S.NegativeOne**(nu+1)*besselj(-nu - S.Half, z))
  804. @assume_integer_order
  805. def _eval_rewrite_as_bessely(self, nu, z, **kwargs):
  806. # jn +- I*yn
  807. # jn as bessely: (-1)**nu * sqrt(pi/(2*z)) * bessely(-nu - S.Half, z)
  808. # yn as bessely: sqrt(pi/(2*z)) * bessely(nu + S.Half, z)
  809. hks = self._hankel_kind_sign
  810. return sqrt(pi/(2*z))*(S.NegativeOne**nu*bessely(-nu - S.Half, z) +
  811. hks*I*bessely(nu + S.Half, z))
  812. def _eval_rewrite_as_yn(self, nu, z, **kwargs):
  813. hks = self._hankel_kind_sign
  814. return jn(nu, z).rewrite(yn) + hks*I*yn(nu, z)
  815. def _eval_rewrite_as_jn(self, nu, z, **kwargs):
  816. hks = self._hankel_kind_sign
  817. return jn(nu, z) + hks*I*yn(nu, z).rewrite(jn)
  818. def _eval_expand_func(self, **hints):
  819. if self.order.is_Integer:
  820. return self._expand(**hints)
  821. else:
  822. nu = self.order
  823. z = self.argument
  824. hks = self._hankel_kind_sign
  825. return jn(nu, z) + hks*I*yn(nu, z)
  826. def _expand(self, **hints):
  827. n = self.order
  828. z = self.argument
  829. hks = self._hankel_kind_sign
  830. # fully expanded version
  831. # return ((fn(n, z) * sin(z) +
  832. # (-1)**(n + 1) * fn(-n - 1, z) * cos(z)) + # jn
  833. # (hks * I * (-1)**(n + 1) *
  834. # (fn(-n - 1, z) * hk * I * sin(z) +
  835. # (-1)**(-n) * fn(n, z) * I * cos(z))) # +-I*yn
  836. # )
  837. return (_jn(n, z) + hks*I*_yn(n, z)).expand()
  838. def _eval_evalf(self, prec):
  839. if self.order.is_Integer:
  840. return self.rewrite(besselj)._eval_evalf(prec)
  841. class hn1(SphericalHankelBase):
  842. r"""
  843. Spherical Hankel function of the first kind.
  844. Explanation
  845. ===========
  846. This function is defined as
  847. .. math:: h_\nu^(1)(z) = j_\nu(z) + i y_\nu(z),
  848. where $j_\nu(z)$ and $y_\nu(z)$ are the spherical
  849. Bessel function of the first and second kinds.
  850. For integral orders $n$, $h_n^(1)$ is calculated using the formula:
  851. .. math:: h_n^(1)(z) = j_{n}(z) + i (-1)^{n+1} j_{-n-1}(z)
  852. Examples
  853. ========
  854. >>> from sympy import Symbol, hn1, hankel1, expand_func, yn, jn
  855. >>> z = Symbol("z")
  856. >>> nu = Symbol("nu", integer=True)
  857. >>> print(expand_func(hn1(nu, z)))
  858. jn(nu, z) + I*yn(nu, z)
  859. >>> print(expand_func(hn1(0, z)))
  860. sin(z)/z - I*cos(z)/z
  861. >>> print(expand_func(hn1(1, z)))
  862. -I*sin(z)/z - cos(z)/z + sin(z)/z**2 - I*cos(z)/z**2
  863. >>> hn1(nu, z).rewrite(jn)
  864. (-1)**(nu + 1)*I*jn(-nu - 1, z) + jn(nu, z)
  865. >>> hn1(nu, z).rewrite(yn)
  866. (-1)**nu*yn(-nu - 1, z) + I*yn(nu, z)
  867. >>> hn1(nu, z).rewrite(hankel1)
  868. sqrt(2)*sqrt(pi)*sqrt(1/z)*hankel1(nu, z)/2
  869. See Also
  870. ========
  871. hn2, jn, yn, hankel1, hankel2
  872. References
  873. ==========
  874. .. [1] https://dlmf.nist.gov/10.47
  875. """
  876. _hankel_kind_sign = S.One
  877. @assume_integer_order
  878. def _eval_rewrite_as_hankel1(self, nu, z, **kwargs):
  879. return sqrt(pi/(2*z))*hankel1(nu, z)
  880. class hn2(SphericalHankelBase):
  881. r"""
  882. Spherical Hankel function of the second kind.
  883. Explanation
  884. ===========
  885. This function is defined as
  886. .. math:: h_\nu^(2)(z) = j_\nu(z) - i y_\nu(z),
  887. where $j_\nu(z)$ and $y_\nu(z)$ are the spherical
  888. Bessel function of the first and second kinds.
  889. For integral orders $n$, $h_n^(2)$ is calculated using the formula:
  890. .. math:: h_n^(2)(z) = j_{n} - i (-1)^{n+1} j_{-n-1}(z)
  891. Examples
  892. ========
  893. >>> from sympy import Symbol, hn2, hankel2, expand_func, jn, yn
  894. >>> z = Symbol("z")
  895. >>> nu = Symbol("nu", integer=True)
  896. >>> print(expand_func(hn2(nu, z)))
  897. jn(nu, z) - I*yn(nu, z)
  898. >>> print(expand_func(hn2(0, z)))
  899. sin(z)/z + I*cos(z)/z
  900. >>> print(expand_func(hn2(1, z)))
  901. I*sin(z)/z - cos(z)/z + sin(z)/z**2 + I*cos(z)/z**2
  902. >>> hn2(nu, z).rewrite(hankel2)
  903. sqrt(2)*sqrt(pi)*sqrt(1/z)*hankel2(nu, z)/2
  904. >>> hn2(nu, z).rewrite(jn)
  905. -(-1)**(nu + 1)*I*jn(-nu - 1, z) + jn(nu, z)
  906. >>> hn2(nu, z).rewrite(yn)
  907. (-1)**nu*yn(-nu - 1, z) - I*yn(nu, z)
  908. See Also
  909. ========
  910. hn1, jn, yn, hankel1, hankel2
  911. References
  912. ==========
  913. .. [1] https://dlmf.nist.gov/10.47
  914. """
  915. _hankel_kind_sign = -S.One
  916. @assume_integer_order
  917. def _eval_rewrite_as_hankel2(self, nu, z, **kwargs):
  918. return sqrt(pi/(2*z))*hankel2(nu, z)
  919. def jn_zeros(n, k, method="sympy", dps=15):
  920. """
  921. Zeros of the spherical Bessel function of the first kind.
  922. Explanation
  923. ===========
  924. This returns an array of zeros of $jn$ up to the $k$-th zero.
  925. * method = "sympy": uses `mpmath.besseljzero
  926. <https://mpmath.org/doc/current/functions/bessel.html#mpmath.besseljzero>`_
  927. * method = "scipy": uses the
  928. `SciPy's sph_jn <https://docs.scipy.org/doc/scipy/reference/generated/scipy.special.jn_zeros.html>`_
  929. and
  930. `newton <https://docs.scipy.org/doc/scipy/reference/generated/scipy.optimize.newton.html>`_
  931. to find all
  932. roots, which is faster than computing the zeros using a general
  933. numerical solver, but it requires SciPy and only works with low
  934. precision floating point numbers. (The function used with
  935. method="sympy" is a recent addition to mpmath; before that a general
  936. solver was used.)
  937. Examples
  938. ========
  939. >>> from sympy import jn_zeros
  940. >>> jn_zeros(2, 4, dps=5)
  941. [5.7635, 9.095, 12.323, 15.515]
  942. See Also
  943. ========
  944. jn, yn, besselj, besselk, bessely
  945. Parameters
  946. ==========
  947. n : integer
  948. order of Bessel function
  949. k : integer
  950. number of zeros to return
  951. """
  952. from math import pi as math_pi
  953. if method == "sympy":
  954. from mpmath import besseljzero
  955. from mpmath.libmp.libmpf import dps_to_prec
  956. prec = dps_to_prec(dps)
  957. return [Expr._from_mpmath(besseljzero(S(n + 0.5)._to_mpmath(prec),
  958. int(l)), prec)
  959. for l in range(1, k + 1)]
  960. elif method == "scipy":
  961. from scipy.optimize import newton
  962. try:
  963. from scipy.special import spherical_jn
  964. f = lambda x: spherical_jn(n, x)
  965. except ImportError:
  966. from scipy.special import sph_jn
  967. f = lambda x: sph_jn(n, x)[0][-1]
  968. else:
  969. raise NotImplementedError("Unknown method.")
  970. def solver(f, x):
  971. if method == "scipy":
  972. root = newton(f, x)
  973. else:
  974. raise NotImplementedError("Unknown method.")
  975. return root
  976. # we need to approximate the position of the first root:
  977. root = n + math_pi
  978. # determine the first root exactly:
  979. root = solver(f, root)
  980. roots = [root]
  981. for i in range(k - 1):
  982. # estimate the position of the next root using the last root + pi:
  983. root = solver(f, root + math_pi)
  984. roots.append(root)
  985. return roots
  986. class AiryBase(Function):
  987. """
  988. Abstract base class for Airy functions.
  989. This class is meant to reduce code duplication.
  990. """
  991. def _eval_conjugate(self):
  992. return self.func(self.args[0].conjugate())
  993. def _eval_is_extended_real(self):
  994. return self.args[0].is_extended_real
  995. def as_real_imag(self, deep=True, **hints):
  996. z = self.args[0]
  997. zc = z.conjugate()
  998. f = self.func
  999. u = (f(z)+f(zc))/2
  1000. v = I*(f(zc)-f(z))/2
  1001. return u, v
  1002. def _eval_expand_complex(self, deep=True, **hints):
  1003. re_part, im_part = self.as_real_imag(deep=deep, **hints)
  1004. return re_part + im_part*I
  1005. class airyai(AiryBase):
  1006. r"""
  1007. The Airy function $\operatorname{Ai}$ of the first kind.
  1008. Explanation
  1009. ===========
  1010. The Airy function $\operatorname{Ai}(z)$ is defined to be the function
  1011. satisfying Airy's differential equation
  1012. .. math::
  1013. \frac{\mathrm{d}^2 w(z)}{\mathrm{d}z^2} - z w(z) = 0.
  1014. Equivalently, for real $z$
  1015. .. math::
  1016. \operatorname{Ai}(z) := \frac{1}{\pi}
  1017. \int_0^\infty \cos\left(\frac{t^3}{3} + z t\right) \mathrm{d}t.
  1018. Examples
  1019. ========
  1020. Create an Airy function object:
  1021. >>> from sympy import airyai
  1022. >>> from sympy.abc import z
  1023. >>> airyai(z)
  1024. airyai(z)
  1025. Several special values are known:
  1026. >>> airyai(0)
  1027. 3**(1/3)/(3*gamma(2/3))
  1028. >>> from sympy import oo
  1029. >>> airyai(oo)
  1030. 0
  1031. >>> airyai(-oo)
  1032. 0
  1033. The Airy function obeys the mirror symmetry:
  1034. >>> from sympy import conjugate
  1035. >>> conjugate(airyai(z))
  1036. airyai(conjugate(z))
  1037. Differentiation with respect to $z$ is supported:
  1038. >>> from sympy import diff
  1039. >>> diff(airyai(z), z)
  1040. airyaiprime(z)
  1041. >>> diff(airyai(z), z, 2)
  1042. z*airyai(z)
  1043. Series expansion is also supported:
  1044. >>> from sympy import series
  1045. >>> series(airyai(z), z, 0, 3)
  1046. 3**(5/6)*gamma(1/3)/(6*pi) - 3**(1/6)*z*gamma(2/3)/(2*pi) + O(z**3)
  1047. We can numerically evaluate the Airy function to arbitrary precision
  1048. on the whole complex plane:
  1049. >>> airyai(-2).evalf(50)
  1050. 0.22740742820168557599192443603787379946077222541710
  1051. Rewrite $\operatorname{Ai}(z)$ in terms of hypergeometric functions:
  1052. >>> from sympy import hyper
  1053. >>> airyai(z).rewrite(hyper)
  1054. -3**(2/3)*z*hyper((), (4/3,), z**3/9)/(3*gamma(1/3)) + 3**(1/3)*hyper((), (2/3,), z**3/9)/(3*gamma(2/3))
  1055. See Also
  1056. ========
  1057. airybi: Airy function of the second kind.
  1058. airyaiprime: Derivative of the Airy function of the first kind.
  1059. airybiprime: Derivative of the Airy function of the second kind.
  1060. References
  1061. ==========
  1062. .. [1] https://en.wikipedia.org/wiki/Airy_function
  1063. .. [2] https://dlmf.nist.gov/9
  1064. .. [3] https://encyclopediaofmath.org/wiki/Airy_functions
  1065. .. [4] https://mathworld.wolfram.com/AiryFunctions.html
  1066. """
  1067. nargs = 1
  1068. unbranched = True
  1069. @classmethod
  1070. def eval(cls, arg):
  1071. if arg.is_Number:
  1072. if arg is S.NaN:
  1073. return S.NaN
  1074. elif arg is S.Infinity:
  1075. return S.Zero
  1076. elif arg is S.NegativeInfinity:
  1077. return S.Zero
  1078. elif arg.is_zero:
  1079. return S.One / (3**Rational(2, 3) * gamma(Rational(2, 3)))
  1080. if arg.is_zero:
  1081. return S.One / (3**Rational(2, 3) * gamma(Rational(2, 3)))
  1082. def fdiff(self, argindex=1):
  1083. if argindex == 1:
  1084. return airyaiprime(self.args[0])
  1085. else:
  1086. raise ArgumentIndexError(self, argindex)
  1087. @staticmethod
  1088. @cacheit
  1089. def taylor_term(n, x, *previous_terms):
  1090. if n < 0:
  1091. return S.Zero
  1092. else:
  1093. x = sympify(x)
  1094. if len(previous_terms) > 1:
  1095. p = previous_terms[-1]
  1096. return ((cbrt(3)*x)**(-n)*(cbrt(3)*x)**(n + 1)*sin(pi*(n*Rational(2, 3) + Rational(4, 3)))*factorial(n) *
  1097. gamma(n/3 + Rational(2, 3))/(sin(pi*(n*Rational(2, 3) + Rational(2, 3)))*factorial(n + 1)*gamma(n/3 + Rational(1, 3))) * p)
  1098. else:
  1099. return (S.One/(3**Rational(2, 3)*pi) * gamma((n+S.One)/S(3)) * sin(Rational(2, 3)*pi*(n+S.One)) /
  1100. factorial(n) * (cbrt(3)*x)**n)
  1101. def _eval_rewrite_as_besselj(self, z, **kwargs):
  1102. ot = Rational(1, 3)
  1103. tt = Rational(2, 3)
  1104. a = Pow(-z, Rational(3, 2))
  1105. if re(z).is_negative:
  1106. return ot*sqrt(-z) * (besselj(-ot, tt*a) + besselj(ot, tt*a))
  1107. def _eval_rewrite_as_besseli(self, z, **kwargs):
  1108. ot = Rational(1, 3)
  1109. tt = Rational(2, 3)
  1110. a = Pow(z, Rational(3, 2))
  1111. if re(z).is_positive:
  1112. return ot*sqrt(z) * (besseli(-ot, tt*a) - besseli(ot, tt*a))
  1113. else:
  1114. return ot*(Pow(a, ot)*besseli(-ot, tt*a) - z*Pow(a, -ot)*besseli(ot, tt*a))
  1115. def _eval_rewrite_as_hyper(self, z, **kwargs):
  1116. pf1 = S.One / (3**Rational(2, 3)*gamma(Rational(2, 3)))
  1117. pf2 = z / (root(3, 3)*gamma(Rational(1, 3)))
  1118. return pf1 * hyper([], [Rational(2, 3)], z**3/9) - pf2 * hyper([], [Rational(4, 3)], z**3/9)
  1119. def _eval_expand_func(self, **hints):
  1120. arg = self.args[0]
  1121. symbs = arg.free_symbols
  1122. if len(symbs) == 1:
  1123. z = symbs.pop()
  1124. c = Wild("c", exclude=[z])
  1125. d = Wild("d", exclude=[z])
  1126. m = Wild("m", exclude=[z])
  1127. n = Wild("n", exclude=[z])
  1128. M = arg.match(c*(d*z**n)**m)
  1129. if M is not None:
  1130. m = M[m]
  1131. # The transformation is given by 03.05.16.0001.01
  1132. # https://functions.wolfram.com/Bessel-TypeFunctions/AiryAi/16/01/01/0001/
  1133. if (3*m).is_integer:
  1134. c = M[c]
  1135. d = M[d]
  1136. n = M[n]
  1137. pf = (d * z**n)**m / (d**m * z**(m*n))
  1138. newarg = c * d**m * z**(m*n)
  1139. return S.Half * ((pf + S.One)*airyai(newarg) - (pf - S.One)/sqrt(3)*airybi(newarg))
  1140. class airybi(AiryBase):
  1141. r"""
  1142. The Airy function $\operatorname{Bi}$ of the second kind.
  1143. Explanation
  1144. ===========
  1145. The Airy function $\operatorname{Bi}(z)$ is defined to be the function
  1146. satisfying Airy's differential equation
  1147. .. math::
  1148. \frac{\mathrm{d}^2 w(z)}{\mathrm{d}z^2} - z w(z) = 0.
  1149. Equivalently, for real $z$
  1150. .. math::
  1151. \operatorname{Bi}(z) := \frac{1}{\pi}
  1152. \int_0^\infty
  1153. \exp\left(-\frac{t^3}{3} + z t\right)
  1154. + \sin\left(\frac{t^3}{3} + z t\right) \mathrm{d}t.
  1155. Examples
  1156. ========
  1157. Create an Airy function object:
  1158. >>> from sympy import airybi
  1159. >>> from sympy.abc import z
  1160. >>> airybi(z)
  1161. airybi(z)
  1162. Several special values are known:
  1163. >>> airybi(0)
  1164. 3**(5/6)/(3*gamma(2/3))
  1165. >>> from sympy import oo
  1166. >>> airybi(oo)
  1167. oo
  1168. >>> airybi(-oo)
  1169. 0
  1170. The Airy function obeys the mirror symmetry:
  1171. >>> from sympy import conjugate
  1172. >>> conjugate(airybi(z))
  1173. airybi(conjugate(z))
  1174. Differentiation with respect to $z$ is supported:
  1175. >>> from sympy import diff
  1176. >>> diff(airybi(z), z)
  1177. airybiprime(z)
  1178. >>> diff(airybi(z), z, 2)
  1179. z*airybi(z)
  1180. Series expansion is also supported:
  1181. >>> from sympy import series
  1182. >>> series(airybi(z), z, 0, 3)
  1183. 3**(1/3)*gamma(1/3)/(2*pi) + 3**(2/3)*z*gamma(2/3)/(2*pi) + O(z**3)
  1184. We can numerically evaluate the Airy function to arbitrary precision
  1185. on the whole complex plane:
  1186. >>> airybi(-2).evalf(50)
  1187. -0.41230258795639848808323405461146104203453483447240
  1188. Rewrite $\operatorname{Bi}(z)$ in terms of hypergeometric functions:
  1189. >>> from sympy import hyper
  1190. >>> airybi(z).rewrite(hyper)
  1191. 3**(1/6)*z*hyper((), (4/3,), z**3/9)/gamma(1/3) + 3**(5/6)*hyper((), (2/3,), z**3/9)/(3*gamma(2/3))
  1192. See Also
  1193. ========
  1194. airyai: Airy function of the first kind.
  1195. airyaiprime: Derivative of the Airy function of the first kind.
  1196. airybiprime: Derivative of the Airy function of the second kind.
  1197. References
  1198. ==========
  1199. .. [1] https://en.wikipedia.org/wiki/Airy_function
  1200. .. [2] https://dlmf.nist.gov/9
  1201. .. [3] https://encyclopediaofmath.org/wiki/Airy_functions
  1202. .. [4] https://mathworld.wolfram.com/AiryFunctions.html
  1203. """
  1204. nargs = 1
  1205. unbranched = True
  1206. @classmethod
  1207. def eval(cls, arg):
  1208. if arg.is_Number:
  1209. if arg is S.NaN:
  1210. return S.NaN
  1211. elif arg is S.Infinity:
  1212. return S.Infinity
  1213. elif arg is S.NegativeInfinity:
  1214. return S.Zero
  1215. elif arg.is_zero:
  1216. return S.One / (3**Rational(1, 6) * gamma(Rational(2, 3)))
  1217. if arg.is_zero:
  1218. return S.One / (3**Rational(1, 6) * gamma(Rational(2, 3)))
  1219. def fdiff(self, argindex=1):
  1220. if argindex == 1:
  1221. return airybiprime(self.args[0])
  1222. else:
  1223. raise ArgumentIndexError(self, argindex)
  1224. @staticmethod
  1225. @cacheit
  1226. def taylor_term(n, x, *previous_terms):
  1227. if n < 0:
  1228. return S.Zero
  1229. else:
  1230. x = sympify(x)
  1231. if len(previous_terms) > 1:
  1232. p = previous_terms[-1]
  1233. return (cbrt(3)*x * Abs(sin(Rational(2, 3)*pi*(n + S.One))) * factorial((n - S.One)/S(3)) /
  1234. ((n + S.One) * Abs(cos(Rational(2, 3)*pi*(n + S.Half))) * factorial((n - 2)/S(3))) * p)
  1235. else:
  1236. return (S.One/(root(3, 6)*pi) * gamma((n + S.One)/S(3)) * Abs(sin(Rational(2, 3)*pi*(n + S.One))) /
  1237. factorial(n) * (cbrt(3)*x)**n)
  1238. def _eval_rewrite_as_besselj(self, z, **kwargs):
  1239. ot = Rational(1, 3)
  1240. tt = Rational(2, 3)
  1241. a = Pow(-z, Rational(3, 2))
  1242. if re(z).is_negative:
  1243. return sqrt(-z/3) * (besselj(-ot, tt*a) - besselj(ot, tt*a))
  1244. def _eval_rewrite_as_besseli(self, z, **kwargs):
  1245. ot = Rational(1, 3)
  1246. tt = Rational(2, 3)
  1247. a = Pow(z, Rational(3, 2))
  1248. if re(z).is_positive:
  1249. return sqrt(z)/sqrt(3) * (besseli(-ot, tt*a) + besseli(ot, tt*a))
  1250. else:
  1251. b = Pow(a, ot)
  1252. c = Pow(a, -ot)
  1253. return sqrt(ot)*(b*besseli(-ot, tt*a) + z*c*besseli(ot, tt*a))
  1254. def _eval_rewrite_as_hyper(self, z, **kwargs):
  1255. pf1 = S.One / (root(3, 6)*gamma(Rational(2, 3)))
  1256. pf2 = z*root(3, 6) / gamma(Rational(1, 3))
  1257. return pf1 * hyper([], [Rational(2, 3)], z**3/9) + pf2 * hyper([], [Rational(4, 3)], z**3/9)
  1258. def _eval_expand_func(self, **hints):
  1259. arg = self.args[0]
  1260. symbs = arg.free_symbols
  1261. if len(symbs) == 1:
  1262. z = symbs.pop()
  1263. c = Wild("c", exclude=[z])
  1264. d = Wild("d", exclude=[z])
  1265. m = Wild("m", exclude=[z])
  1266. n = Wild("n", exclude=[z])
  1267. M = arg.match(c*(d*z**n)**m)
  1268. if M is not None:
  1269. m = M[m]
  1270. # The transformation is given by 03.06.16.0001.01
  1271. # https://functions.wolfram.com/Bessel-TypeFunctions/AiryBi/16/01/01/0001/
  1272. if (3*m).is_integer:
  1273. c = M[c]
  1274. d = M[d]
  1275. n = M[n]
  1276. pf = (d * z**n)**m / (d**m * z**(m*n))
  1277. newarg = c * d**m * z**(m*n)
  1278. return S.Half * (sqrt(3)*(S.One - pf)*airyai(newarg) + (S.One + pf)*airybi(newarg))
  1279. class airyaiprime(AiryBase):
  1280. r"""
  1281. The derivative $\operatorname{Ai}^\prime$ of the Airy function of the first
  1282. kind.
  1283. Explanation
  1284. ===========
  1285. The Airy function $\operatorname{Ai}^\prime(z)$ is defined to be the
  1286. function
  1287. .. math::
  1288. \operatorname{Ai}^\prime(z) := \frac{\mathrm{d} \operatorname{Ai}(z)}{\mathrm{d} z}.
  1289. Examples
  1290. ========
  1291. Create an Airy function object:
  1292. >>> from sympy import airyaiprime
  1293. >>> from sympy.abc import z
  1294. >>> airyaiprime(z)
  1295. airyaiprime(z)
  1296. Several special values are known:
  1297. >>> airyaiprime(0)
  1298. -3**(2/3)/(3*gamma(1/3))
  1299. >>> from sympy import oo
  1300. >>> airyaiprime(oo)
  1301. 0
  1302. The Airy function obeys the mirror symmetry:
  1303. >>> from sympy import conjugate
  1304. >>> conjugate(airyaiprime(z))
  1305. airyaiprime(conjugate(z))
  1306. Differentiation with respect to $z$ is supported:
  1307. >>> from sympy import diff
  1308. >>> diff(airyaiprime(z), z)
  1309. z*airyai(z)
  1310. >>> diff(airyaiprime(z), z, 2)
  1311. z*airyaiprime(z) + airyai(z)
  1312. Series expansion is also supported:
  1313. >>> from sympy import series
  1314. >>> series(airyaiprime(z), z, 0, 3)
  1315. -3**(2/3)/(3*gamma(1/3)) + 3**(1/3)*z**2/(6*gamma(2/3)) + O(z**3)
  1316. We can numerically evaluate the Airy function to arbitrary precision
  1317. on the whole complex plane:
  1318. >>> airyaiprime(-2).evalf(50)
  1319. 0.61825902074169104140626429133247528291577794512415
  1320. Rewrite $\operatorname{Ai}^\prime(z)$ in terms of hypergeometric functions:
  1321. >>> from sympy import hyper
  1322. >>> airyaiprime(z).rewrite(hyper)
  1323. 3**(1/3)*z**2*hyper((), (5/3,), z**3/9)/(6*gamma(2/3)) - 3**(2/3)*hyper((), (1/3,), z**3/9)/(3*gamma(1/3))
  1324. See Also
  1325. ========
  1326. airyai: Airy function of the first kind.
  1327. airybi: Airy function of the second kind.
  1328. airybiprime: Derivative of the Airy function of the second kind.
  1329. References
  1330. ==========
  1331. .. [1] https://en.wikipedia.org/wiki/Airy_function
  1332. .. [2] https://dlmf.nist.gov/9
  1333. .. [3] https://encyclopediaofmath.org/wiki/Airy_functions
  1334. .. [4] https://mathworld.wolfram.com/AiryFunctions.html
  1335. """
  1336. nargs = 1
  1337. unbranched = True
  1338. @classmethod
  1339. def eval(cls, arg):
  1340. if arg.is_Number:
  1341. if arg is S.NaN:
  1342. return S.NaN
  1343. elif arg is S.Infinity:
  1344. return S.Zero
  1345. if arg.is_zero:
  1346. return S.NegativeOne / (3**Rational(1, 3) * gamma(Rational(1, 3)))
  1347. def fdiff(self, argindex=1):
  1348. if argindex == 1:
  1349. return self.args[0]*airyai(self.args[0])
  1350. else:
  1351. raise ArgumentIndexError(self, argindex)
  1352. def _eval_evalf(self, prec):
  1353. z = self.args[0]._to_mpmath(prec)
  1354. with workprec(prec):
  1355. res = mp.airyai(z, derivative=1)
  1356. return Expr._from_mpmath(res, prec)
  1357. def _eval_rewrite_as_besselj(self, z, **kwargs):
  1358. tt = Rational(2, 3)
  1359. a = Pow(-z, Rational(3, 2))
  1360. if re(z).is_negative:
  1361. return z/3 * (besselj(-tt, tt*a) - besselj(tt, tt*a))
  1362. def _eval_rewrite_as_besseli(self, z, **kwargs):
  1363. ot = Rational(1, 3)
  1364. tt = Rational(2, 3)
  1365. a = tt * Pow(z, Rational(3, 2))
  1366. if re(z).is_positive:
  1367. return z/3 * (besseli(tt, a) - besseli(-tt, a))
  1368. else:
  1369. a = Pow(z, Rational(3, 2))
  1370. b = Pow(a, tt)
  1371. c = Pow(a, -tt)
  1372. return ot * (z**2*c*besseli(tt, tt*a) - b*besseli(-ot, tt*a))
  1373. def _eval_rewrite_as_hyper(self, z, **kwargs):
  1374. pf1 = z**2 / (2*3**Rational(2, 3)*gamma(Rational(2, 3)))
  1375. pf2 = 1 / (root(3, 3)*gamma(Rational(1, 3)))
  1376. return pf1 * hyper([], [Rational(5, 3)], z**3/9) - pf2 * hyper([], [Rational(1, 3)], z**3/9)
  1377. def _eval_expand_func(self, **hints):
  1378. arg = self.args[0]
  1379. symbs = arg.free_symbols
  1380. if len(symbs) == 1:
  1381. z = symbs.pop()
  1382. c = Wild("c", exclude=[z])
  1383. d = Wild("d", exclude=[z])
  1384. m = Wild("m", exclude=[z])
  1385. n = Wild("n", exclude=[z])
  1386. M = arg.match(c*(d*z**n)**m)
  1387. if M is not None:
  1388. m = M[m]
  1389. # The transformation is in principle
  1390. # given by 03.07.16.0001.01 but note
  1391. # that there is an error in this formula.
  1392. # https://functions.wolfram.com/Bessel-TypeFunctions/AiryAiPrime/16/01/01/0001/
  1393. if (3*m).is_integer:
  1394. c = M[c]
  1395. d = M[d]
  1396. n = M[n]
  1397. pf = (d**m * z**(n*m)) / (d * z**n)**m
  1398. newarg = c * d**m * z**(n*m)
  1399. return S.Half * ((pf + S.One)*airyaiprime(newarg) + (pf - S.One)/sqrt(3)*airybiprime(newarg))
  1400. class airybiprime(AiryBase):
  1401. r"""
  1402. The derivative $\operatorname{Bi}^\prime$ of the Airy function of the first
  1403. kind.
  1404. Explanation
  1405. ===========
  1406. The Airy function $\operatorname{Bi}^\prime(z)$ is defined to be the
  1407. function
  1408. .. math::
  1409. \operatorname{Bi}^\prime(z) := \frac{\mathrm{d} \operatorname{Bi}(z)}{\mathrm{d} z}.
  1410. Examples
  1411. ========
  1412. Create an Airy function object:
  1413. >>> from sympy import airybiprime
  1414. >>> from sympy.abc import z
  1415. >>> airybiprime(z)
  1416. airybiprime(z)
  1417. Several special values are known:
  1418. >>> airybiprime(0)
  1419. 3**(1/6)/gamma(1/3)
  1420. >>> from sympy import oo
  1421. >>> airybiprime(oo)
  1422. oo
  1423. >>> airybiprime(-oo)
  1424. 0
  1425. The Airy function obeys the mirror symmetry:
  1426. >>> from sympy import conjugate
  1427. >>> conjugate(airybiprime(z))
  1428. airybiprime(conjugate(z))
  1429. Differentiation with respect to $z$ is supported:
  1430. >>> from sympy import diff
  1431. >>> diff(airybiprime(z), z)
  1432. z*airybi(z)
  1433. >>> diff(airybiprime(z), z, 2)
  1434. z*airybiprime(z) + airybi(z)
  1435. Series expansion is also supported:
  1436. >>> from sympy import series
  1437. >>> series(airybiprime(z), z, 0, 3)
  1438. 3**(1/6)/gamma(1/3) + 3**(5/6)*z**2/(6*gamma(2/3)) + O(z**3)
  1439. We can numerically evaluate the Airy function to arbitrary precision
  1440. on the whole complex plane:
  1441. >>> airybiprime(-2).evalf(50)
  1442. 0.27879516692116952268509756941098324140300059345163
  1443. Rewrite $\operatorname{Bi}^\prime(z)$ in terms of hypergeometric functions:
  1444. >>> from sympy import hyper
  1445. >>> airybiprime(z).rewrite(hyper)
  1446. 3**(5/6)*z**2*hyper((), (5/3,), z**3/9)/(6*gamma(2/3)) + 3**(1/6)*hyper((), (1/3,), z**3/9)/gamma(1/3)
  1447. See Also
  1448. ========
  1449. airyai: Airy function of the first kind.
  1450. airybi: Airy function of the second kind.
  1451. airyaiprime: Derivative of the Airy function of the first kind.
  1452. References
  1453. ==========
  1454. .. [1] https://en.wikipedia.org/wiki/Airy_function
  1455. .. [2] https://dlmf.nist.gov/9
  1456. .. [3] https://encyclopediaofmath.org/wiki/Airy_functions
  1457. .. [4] https://mathworld.wolfram.com/AiryFunctions.html
  1458. """
  1459. nargs = 1
  1460. unbranched = True
  1461. @classmethod
  1462. def eval(cls, arg):
  1463. if arg.is_Number:
  1464. if arg is S.NaN:
  1465. return S.NaN
  1466. elif arg is S.Infinity:
  1467. return S.Infinity
  1468. elif arg is S.NegativeInfinity:
  1469. return S.Zero
  1470. elif arg.is_zero:
  1471. return 3**Rational(1, 6) / gamma(Rational(1, 3))
  1472. if arg.is_zero:
  1473. return 3**Rational(1, 6) / gamma(Rational(1, 3))
  1474. def fdiff(self, argindex=1):
  1475. if argindex == 1:
  1476. return self.args[0]*airybi(self.args[0])
  1477. else:
  1478. raise ArgumentIndexError(self, argindex)
  1479. def _eval_evalf(self, prec):
  1480. z = self.args[0]._to_mpmath(prec)
  1481. with workprec(prec):
  1482. res = mp.airybi(z, derivative=1)
  1483. return Expr._from_mpmath(res, prec)
  1484. def _eval_rewrite_as_besselj(self, z, **kwargs):
  1485. tt = Rational(2, 3)
  1486. a = tt * Pow(-z, Rational(3, 2))
  1487. if re(z).is_negative:
  1488. return -z/sqrt(3) * (besselj(-tt, a) + besselj(tt, a))
  1489. def _eval_rewrite_as_besseli(self, z, **kwargs):
  1490. ot = Rational(1, 3)
  1491. tt = Rational(2, 3)
  1492. a = tt * Pow(z, Rational(3, 2))
  1493. if re(z).is_positive:
  1494. return z/sqrt(3) * (besseli(-tt, a) + besseli(tt, a))
  1495. else:
  1496. a = Pow(z, Rational(3, 2))
  1497. b = Pow(a, tt)
  1498. c = Pow(a, -tt)
  1499. return sqrt(ot) * (b*besseli(-tt, tt*a) + z**2*c*besseli(tt, tt*a))
  1500. def _eval_rewrite_as_hyper(self, z, **kwargs):
  1501. pf1 = z**2 / (2*root(3, 6)*gamma(Rational(2, 3)))
  1502. pf2 = root(3, 6) / gamma(Rational(1, 3))
  1503. return pf1 * hyper([], [Rational(5, 3)], z**3/9) + pf2 * hyper([], [Rational(1, 3)], z**3/9)
  1504. def _eval_expand_func(self, **hints):
  1505. arg = self.args[0]
  1506. symbs = arg.free_symbols
  1507. if len(symbs) == 1:
  1508. z = symbs.pop()
  1509. c = Wild("c", exclude=[z])
  1510. d = Wild("d", exclude=[z])
  1511. m = Wild("m", exclude=[z])
  1512. n = Wild("n", exclude=[z])
  1513. M = arg.match(c*(d*z**n)**m)
  1514. if M is not None:
  1515. m = M[m]
  1516. # The transformation is in principle
  1517. # given by 03.08.16.0001.01 but note
  1518. # that there is an error in this formula.
  1519. # https://functions.wolfram.com/Bessel-TypeFunctions/AiryBiPrime/16/01/01/0001/
  1520. if (3*m).is_integer:
  1521. c = M[c]
  1522. d = M[d]
  1523. n = M[n]
  1524. pf = (d**m * z**(n*m)) / (d * z**n)**m
  1525. newarg = c * d**m * z**(n*m)
  1526. return S.Half * (sqrt(3)*(pf - S.One)*airyaiprime(newarg) + (pf + S.One)*airybiprime(newarg))
  1527. class marcumq(Function):
  1528. r"""
  1529. The Marcum Q-function.
  1530. Explanation
  1531. ===========
  1532. The Marcum Q-function is defined by the meromorphic continuation of
  1533. .. math::
  1534. Q_m(a, b) = a^{- m + 1} \int_{b}^{\infty} x^{m} e^{- \frac{a^{2}}{2} - \frac{x^{2}}{2}} I_{m - 1}\left(a x\right)\, dx
  1535. Examples
  1536. ========
  1537. >>> from sympy import marcumq
  1538. >>> from sympy.abc import m, a, b
  1539. >>> marcumq(m, a, b)
  1540. marcumq(m, a, b)
  1541. Special values:
  1542. >>> marcumq(m, 0, b)
  1543. uppergamma(m, b**2/2)/gamma(m)
  1544. >>> marcumq(0, 0, 0)
  1545. 0
  1546. >>> marcumq(0, a, 0)
  1547. 1 - exp(-a**2/2)
  1548. >>> marcumq(1, a, a)
  1549. 1/2 + exp(-a**2)*besseli(0, a**2)/2
  1550. >>> marcumq(2, a, a)
  1551. 1/2 + exp(-a**2)*besseli(0, a**2)/2 + exp(-a**2)*besseli(1, a**2)
  1552. Differentiation with respect to $a$ and $b$ is supported:
  1553. >>> from sympy import diff
  1554. >>> diff(marcumq(m, a, b), a)
  1555. a*(-marcumq(m, a, b) + marcumq(m + 1, a, b))
  1556. >>> diff(marcumq(m, a, b), b)
  1557. -a**(1 - m)*b**m*exp(-a**2/2 - b**2/2)*besseli(m - 1, a*b)
  1558. References
  1559. ==========
  1560. .. [1] https://en.wikipedia.org/wiki/Marcum_Q-function
  1561. .. [2] https://mathworld.wolfram.com/MarcumQ-Function.html
  1562. """
  1563. @classmethod
  1564. def eval(cls, m, a, b):
  1565. if a is S.Zero:
  1566. if m is S.Zero and b is S.Zero:
  1567. return S.Zero
  1568. return uppergamma(m, b**2 * S.Half) / gamma(m)
  1569. if m is S.Zero and b is S.Zero:
  1570. return 1 - 1 / exp(a**2 * S.Half)
  1571. if a == b:
  1572. if m is S.One:
  1573. return (1 + exp(-a**2) * besseli(0, a**2))*S.Half
  1574. if m == 2:
  1575. return S.Half + S.Half * exp(-a**2) * besseli(0, a**2) + exp(-a**2) * besseli(1, a**2)
  1576. if a.is_zero:
  1577. if m.is_zero and b.is_zero:
  1578. return S.Zero
  1579. return uppergamma(m, b**2*S.Half) / gamma(m)
  1580. if m.is_zero and b.is_zero:
  1581. return 1 - 1 / exp(a**2*S.Half)
  1582. def fdiff(self, argindex=2):
  1583. m, a, b = self.args
  1584. if argindex == 2:
  1585. return a * (-marcumq(m, a, b) + marcumq(1+m, a, b))
  1586. elif argindex == 3:
  1587. return (-b**m / a**(m-1)) * exp(-(a**2 + b**2)/2) * besseli(m-1, a*b)
  1588. else:
  1589. raise ArgumentIndexError(self, argindex)
  1590. def _eval_rewrite_as_Integral(self, m, a, b, **kwargs):
  1591. from sympy.integrals.integrals import Integral
  1592. x = kwargs.get('x', Dummy('x'))
  1593. return a ** (1 - m) * \
  1594. Integral(x**m * exp(-(x**2 + a**2)/2) * besseli(m-1, a*x), [x, b, S.Infinity])
  1595. def _eval_rewrite_as_Sum(self, m, a, b, **kwargs):
  1596. from sympy.concrete.summations import Sum
  1597. k = kwargs.get('k', Dummy('k'))
  1598. return exp(-(a**2 + b**2) / 2) * Sum((a/b)**k * besseli(k, a*b), [k, 1-m, S.Infinity])
  1599. def _eval_rewrite_as_besseli(self, m, a, b, **kwargs):
  1600. if a == b:
  1601. if m == 1:
  1602. return (1 + exp(-a**2) * besseli(0, a**2)) / 2
  1603. if m.is_Integer and m >= 2:
  1604. s = sum([besseli(i, a**2) for i in range(1, m)])
  1605. return S.Half + exp(-a**2) * besseli(0, a**2) / 2 + exp(-a**2) * s
  1606. def _eval_is_zero(self):
  1607. if all(arg.is_zero for arg in self.args):
  1608. return True