bessel.py 67 KB

123456789101112131415161718192021222324252627282930313233343536373839404142434445464748495051525354555657585960616263646566676869707172737475767778798081828384858687888990919293949596979899100101102103104105106107108109110111112113114115116117118119120121122123124125126127128129130131132133134135136137138139140141142143144145146147148149150151152153154155156157158159160161162163164165166167168169170171172173174175176177178179180181182183184185186187188189190191192193194195196197198199200201202203204205206207208209210211212213214215216217218219220221222223224225226227228229230231232233234235236237238239240241242243244245246247248249250251252253254255256257258259260261262263264265266267268269270271272273274275276277278279280281282283284285286287288289290291292293294295296297298299300301302303304305306307308309310311312313314315316317318319320321322323324325326327328329330331332333334335336337338339340341342343344345346347348349350351352353354355356357358359360361362363364365366367368369370371372373374375376377378379380381382383384385386387388389390391392393394395396397398399400401402403404405406407408409410411412413414415416417418419420421422423424425426427428429430431432433434435436437438439440441442443444445446447448449450451452453454455456457458459460461462463464465466467468469470471472473474475476477478479480481482483484485486487488489490491492493494495496497498499500501502503504505506507508509510511512513514515516517518519520521522523524525526527528529530531532533534535536537538539540541542543544545546547548549550551552553554555556557558559560561562563564565566567568569570571572573574575576577578579580581582583584585586587588589590591592593594595596597598599600601602603604605606607608609610611612613614615616617618619620621622623624625626627628629630631632633634635636637638639640641642643644645646647648649650651652653654655656657658659660661662663664665666667668669670671672673674675676677678679680681682683684685686687688689690691692693694695696697698699700701702703704705706707708709710711712713714715716717718719720721722723724725726727728729730731732733734735736737738739740741742743744745746747748749750751752753754755756757758759760761762763764765766767768769770771772773774775776777778779780781782783784785786787788789790791792793794795796797798799800801802803804805806807808809810811812813814815816817818819820821822823824825826827828829830831832833834835836837838839840841842843844845846847848849850851852853854855856857858859860861862863864865866867868869870871872873874875876877878879880881882883884885886887888889890891892893894895896897898899900901902903904905906907908909910911912913914915916917918919920921922923924925926927928929930931932933934935936937938939940941942943944945946947948949950951952953954955956957958959960961962963964965966967968969970971972973974975976977978979980981982983984985986987988989990991992993994995996997998999100010011002100310041005100610071008100910101011101210131014101510161017101810191020102110221023102410251026102710281029103010311032103310341035103610371038103910401041104210431044104510461047104810491050105110521053105410551056105710581059106010611062106310641065106610671068106910701071107210731074107510761077107810791080108110821083108410851086108710881089109010911092109310941095109610971098109911001101110211031104110511061107110811091110111111121113111411151116111711181119112011211122112311241125112611271128112911301131113211331134113511361137113811391140114111421143114411451146114711481149115011511152115311541155115611571158115911601161116211631164116511661167116811691170117111721173117411751176117711781179118011811182118311841185118611871188118911901191119211931194119511961197119811991200120112021203120412051206120712081209121012111212121312141215121612171218121912201221122212231224122512261227122812291230123112321233123412351236123712381239124012411242124312441245124612471248124912501251125212531254125512561257125812591260126112621263126412651266126712681269127012711272127312741275127612771278127912801281128212831284128512861287128812891290129112921293129412951296129712981299130013011302130313041305130613071308130913101311131213131314131513161317131813191320132113221323132413251326132713281329133013311332133313341335133613371338133913401341134213431344134513461347134813491350135113521353135413551356135713581359136013611362136313641365136613671368136913701371137213731374137513761377137813791380138113821383138413851386138713881389139013911392139313941395139613971398139914001401140214031404140514061407140814091410141114121413141414151416141714181419142014211422142314241425142614271428142914301431143214331434143514361437143814391440144114421443144414451446144714481449145014511452145314541455145614571458145914601461146214631464146514661467146814691470147114721473147414751476147714781479148014811482148314841485148614871488148914901491149214931494149514961497149814991500150115021503150415051506150715081509151015111512151315141515151615171518151915201521152215231524152515261527152815291530153115321533153415351536153715381539154015411542154315441545154615471548154915501551155215531554155515561557155815591560156115621563156415651566156715681569157015711572157315741575157615771578157915801581158215831584158515861587158815891590159115921593159415951596159715981599160016011602160316041605160616071608160916101611161216131614161516161617161816191620162116221623162416251626162716281629163016311632163316341635163616371638163916401641164216431644164516461647164816491650165116521653165416551656165716581659166016611662166316641665166616671668166916701671167216731674167516761677167816791680168116821683168416851686168716881689169016911692169316941695169616971698169917001701170217031704170517061707170817091710171117121713171417151716171717181719172017211722172317241725172617271728172917301731173217331734173517361737173817391740174117421743174417451746174717481749175017511752175317541755175617571758175917601761176217631764176517661767176817691770177117721773177417751776177717781779178017811782178317841785178617871788178917901791179217931794179517961797179817991800180118021803180418051806180718081809181018111812181318141815181618171818181918201821182218231824182518261827182818291830183118321833183418351836183718381839184018411842184318441845184618471848184918501851185218531854185518561857185818591860186118621863186418651866186718681869187018711872187318741875187618771878187918801881188218831884188518861887188818891890189118921893189418951896189718981899190019011902190319041905190619071908190919101911191219131914191519161917191819191920192119221923192419251926192719281929193019311932193319341935193619371938193919401941194219431944194519461947194819491950195119521953195419551956195719581959196019611962196319641965196619671968196919701971197219731974197519761977197819791980198119821983198419851986198719881989199019911992199319941995199619971998199920002001200220032004200520062007200820092010201120122013201420152016201720182019202020212022202320242025202620272028202920302031203220332034203520362037203820392040204120422043204420452046204720482049205020512052205320542055205620572058205920602061206220632064206520662067206820692070207120722073207420752076207720782079208020812082208320842085208620872088208920902091209220932094209520962097209820992100210121022103210421052106210721082109211021112112211321142115211621172118211921202121212221232124212521262127212821292130213121322133213421352136213721382139214021412142214321442145214621472148214921502151215221532154215521562157215821592160216121622163216421652166216721682169217021712172217321742175217621772178217921802181218221832184218521862187218821892190219121922193219421952196219721982199220022012202220322042205220622072208
  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 DefinedFunction, 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, uniquely_named_symbol, Wild
  11. from sympy.core.sympify import sympify
  12. from sympy.functions.combinatorial.factorials import factorial, RisingFactorial
  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(DefinedFunction):
  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, cdir):
  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=logx, cdir=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, cdir):
  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=logx, cdir=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_tractable(self, nu, z, limitvar=None, **kwargs):
  414. if z.is_extended_real:
  415. return exp(z)*_besseli(nu, z)
  416. def _eval_rewrite_as_besselj(self, nu, z, **kwargs):
  417. return exp(-I*pi*nu/2)*besselj(nu, polar_lift(I)*z)
  418. def _eval_rewrite_as_bessely(self, nu, z, **kwargs):
  419. aj = self._eval_rewrite_as_besselj(*self.args)
  420. if aj:
  421. return aj.rewrite(bessely)
  422. def _eval_rewrite_as_jn(self, nu, z, **kwargs):
  423. return self._eval_rewrite_as_besselj(*self.args).rewrite(jn)
  424. def _eval_is_extended_real(self):
  425. nu, z = self.args
  426. if nu.is_integer and z.is_extended_real:
  427. return True
  428. def _eval_as_leading_term(self, x, logx, cdir):
  429. nu, z = self.args
  430. try:
  431. arg = z.as_leading_term(x)
  432. except NotImplementedError:
  433. return self
  434. c, e = arg.as_coeff_exponent(x)
  435. if e.is_positive:
  436. return arg**nu/(2**nu*gamma(nu + 1))
  437. elif e.is_negative:
  438. cdir = 1 if cdir == 0 else cdir
  439. sign = c*cdir**e
  440. if not sign.is_negative:
  441. # Refer Abramowitz and Stegun 1965, p. 377 for more information on
  442. # asymptotic approximation of besseli function.
  443. return exp(z)/sqrt(2*pi*z)
  444. return self
  445. return super(besseli, self)._eval_as_leading_term(x, logx=logx, cdir=cdir)
  446. def _eval_nseries(self, x, n, logx, cdir=0):
  447. # Refer https://functions.wolfram.com/Bessel-TypeFunctions/BesselI/06/01/04/01/01/0003/
  448. # for more information on nseries expansion of besseli function.
  449. from sympy.series.order import Order
  450. nu, z = self.args
  451. # In case of powers less than 1, number of terms need to be computed
  452. # separately to avoid repeated callings of _eval_nseries with wrong n
  453. try:
  454. _, exp = z.leadterm(x)
  455. except (ValueError, NotImplementedError):
  456. return self
  457. if exp.is_positive:
  458. newn = ceiling(n/exp)
  459. o = Order(x**n, x)
  460. r = (z/2)._eval_nseries(x, n, logx, cdir).removeO()
  461. if r is S.Zero:
  462. return o
  463. t = (_mexpand(r**2) + o).removeO()
  464. term = r**nu/gamma(nu + 1)
  465. s = [term]
  466. for k in range(1, (newn + 1)//2):
  467. term *= t/(k*(nu + k))
  468. term = (_mexpand(term) + o).removeO()
  469. s.append(term)
  470. return Add(*s) + o
  471. return super(besseli, self)._eval_nseries(x, n, logx, cdir)
  472. def _eval_aseries(self, n, args0, x, logx):
  473. from sympy.functions.combinatorial.factorials import RisingFactorial
  474. from sympy.series.order import Order
  475. point = args0[1]
  476. if point in [S.Infinity, S.NegativeInfinity]:
  477. nu, z = self.args
  478. s = [(RisingFactorial(Rational(2*nu - 1, 2), k)*RisingFactorial(Rational(2*nu + 1, 2), k))/\
  479. ((2)**(k)*z**(Rational(2*k + 1, 2))*factorial(k)) for k in range(n)] + [Order(1/z**(Rational(2*n + 1, 2)), x)]
  480. return exp(z)/sqrt(2*pi) * (Add(*s))
  481. return super()._eval_aseries(n, args0, x, logx)
  482. class besselk(BesselBase):
  483. r"""
  484. Modified Bessel function of the second kind.
  485. Explanation
  486. ===========
  487. The Bessel $K$ function of order $\nu$ is defined as
  488. .. math ::
  489. K_\nu(z) = \lim_{\mu \to \nu} \frac{\pi}{2}
  490. \frac{I_{-\mu}(z) -I_\mu(z)}{\sin(\pi \mu)},
  491. where $I_\mu(z)$ is the modified Bessel function of the first kind.
  492. It is a solution of the modified Bessel equation, and linearly independent
  493. from $Y_\nu$.
  494. Examples
  495. ========
  496. >>> from sympy import besselk
  497. >>> from sympy.abc import z, n
  498. >>> besselk(n, z).diff(z)
  499. -besselk(n - 1, z)/2 - besselk(n + 1, z)/2
  500. See Also
  501. ========
  502. besselj, besseli, bessely
  503. References
  504. ==========
  505. .. [1] https://functions.wolfram.com/Bessel-TypeFunctions/BesselK/
  506. """
  507. _a = S.One
  508. _b = -S.One
  509. @classmethod
  510. def eval(cls, nu, z):
  511. if z.is_zero:
  512. if nu.is_zero:
  513. return S.Infinity
  514. elif re(nu).is_zero is False:
  515. return S.ComplexInfinity
  516. elif re(nu).is_zero:
  517. return S.NaN
  518. if z in (S.Infinity, I*S.Infinity, I*S.NegativeInfinity):
  519. return S.Zero
  520. if nu.is_integer:
  521. if nu.could_extract_minus_sign():
  522. return besselk(-nu, z)
  523. def _eval_rewrite_as_besseli(self, nu, z, **kwargs):
  524. if nu.is_integer is False:
  525. return pi*csc(pi*nu)*(besseli(-nu, z) - besseli(nu, z))/2
  526. def _eval_rewrite_as_besselj(self, nu, z, **kwargs):
  527. ai = self._eval_rewrite_as_besseli(*self.args)
  528. if ai:
  529. return ai.rewrite(besselj)
  530. def _eval_rewrite_as_bessely(self, nu, z, **kwargs):
  531. aj = self._eval_rewrite_as_besselj(*self.args)
  532. if aj:
  533. return aj.rewrite(bessely)
  534. def _eval_rewrite_as_yn(self, nu, z, **kwargs):
  535. ay = self._eval_rewrite_as_bessely(*self.args)
  536. if ay:
  537. return ay.rewrite(yn)
  538. def _eval_is_extended_real(self):
  539. nu, z = self.args
  540. if nu.is_integer and z.is_positive:
  541. return True
  542. def _eval_rewrite_as_tractable(self, nu, z, limitvar=None, **kwargs):
  543. if z.is_extended_real:
  544. return exp(-z)*_besselk(nu, z)
  545. def _eval_as_leading_term(self, x, logx, cdir):
  546. nu, z = self.args
  547. try:
  548. arg = z.as_leading_term(x)
  549. except NotImplementedError:
  550. return self
  551. _, e = arg.as_coeff_exponent(x)
  552. if e.is_positive:
  553. if nu.is_zero:
  554. # Equation 9.6.8 of Abramowitz and Stegun (10th ed, 1972).
  555. term = -log(z) - S.EulerGamma + log(2)
  556. elif nu.is_nonzero:
  557. # Equation 9.6.9 of Abramowitz and Stegun (10th ed, 1972).
  558. term = gamma(Abs(nu))*(z/2)**(-Abs(nu))/2
  559. else:
  560. raise NotImplementedError(f"Cannot proceed without knowing if {nu} is zero or not.")
  561. return term.as_leading_term(x, logx=logx)
  562. elif e.is_negative:
  563. # Equation 9.7.2 of Abramowitz and Stegun (10th ed, 1972).
  564. return sqrt(pi)*exp(-arg)/sqrt(2*arg)
  565. else:
  566. return self.func(nu, arg)
  567. def _eval_nseries(self, x, n, logx, cdir=0):
  568. from sympy.series.order import Order
  569. nu, z = self.args
  570. try:
  571. _, exp = z.leadterm(x)
  572. except (ValueError, NotImplementedError):
  573. return self
  574. # In case of powers less than 1, number of terms need to be computed
  575. # separately to avoid repeated callings of _eval_nseries with wrong n
  576. if exp.is_positive:
  577. r = (z/2)._eval_nseries(x, n, logx, cdir).removeO()
  578. if r is S.Zero:
  579. return Order(z**(-nu) + z**nu, x)
  580. o = Order(x**n, x)
  581. if nu.is_integer:
  582. # Reference: https://functions.wolfram.com/Bessel-TypeFunctions/BesselK/06/01/04/01/02/0008/ (only for integer order)
  583. newn = ceiling(n/exp)
  584. bn = besseli(nu, z)
  585. a = ((-1)**(nu - 1)*log(z/2)*bn)._eval_nseries(x, n, logx, cdir)
  586. b, c = [], []
  587. t = _mexpand(r**2)
  588. if nu > S.Zero:
  589. term = r**(-nu)*factorial(nu - 1)/2
  590. b.append(term)
  591. for k in range(1, nu):
  592. term *= t/((k - nu)*k)
  593. term = (_mexpand(term) + o).removeO()
  594. b.append(term)
  595. p = r**nu*(-1)**nu/(2*factorial(nu))
  596. term = p*(digamma(nu + 1) - S.EulerGamma)
  597. c.append(term)
  598. for k in range(1, (newn + 1)//2):
  599. p *= t/(k*(k + nu))
  600. p = (_mexpand(p) + o).removeO()
  601. term = p*(digamma(k + nu + 1) + digamma(k + 1))
  602. c.append(term)
  603. return a + Add(*b) + Add(*c) + o
  604. elif nu.is_noninteger:
  605. # Reference: https://functions.wolfram.com/Bessel-TypeFunctions/BesselK/06/01/04/01/01/0003/
  606. # (only for non-integer order).
  607. # While the expression in the reference above seems correct
  608. # for non-real order as well, it would need some manipulation
  609. # (not implemented) to be written as a power series in x with
  610. # real exponents [e.g. Dunster 1990. "Bessel functions
  611. # of purely imaginary order, with an application to second-order
  612. # linear differential equations having a large parameter".
  613. # SIAM J. Math. Anal. Vol 21, No. 4, pp 995-1018.].
  614. newn_a = ceiling((n+nu)/exp)
  615. newn_b = ceiling((n-nu)/exp)
  616. a, b = [], []
  617. for k in range((newn_a+1)//2):
  618. term = gamma(nu)*r**(2*k-nu)/(2*RisingFactorial(1-nu, k)*factorial(k))
  619. a.append(_mexpand(term))
  620. for k in range((newn_b+1)//2):
  621. term = gamma(-nu)*r**(2*k+nu)/(2*RisingFactorial(nu+1, k)*factorial(k))
  622. b.append(_mexpand(term))
  623. return Add(*a) + Add(*b) + o
  624. else:
  625. raise NotImplementedError("besselk expansion is only implemented for real order")
  626. return super(besselk, self)._eval_nseries(x, n, logx, cdir)
  627. def _eval_aseries(self, n, args0, x, logx):
  628. from sympy.functions.combinatorial.factorials import RisingFactorial
  629. from sympy.series.order import Order
  630. point = args0[1]
  631. if point in [S.Infinity, S.NegativeInfinity]:
  632. nu, z = self.args
  633. s = [(RisingFactorial(Rational(2*nu - 1, 2), k)*RisingFactorial(Rational(2*nu + 1, 2), k))/\
  634. ((-2)**(k)*z**(Rational(2*k + 1, 2))*factorial(k)) for k in range(n)] +[Order(1/z**(Rational(2*n + 1, 2)), x)]
  635. return (exp(-z)*sqrt(pi/2))*Add(*s)
  636. return super()._eval_aseries(n, args0, x, logx)
  637. class hankel1(BesselBase):
  638. r"""
  639. Hankel function of the first kind.
  640. Explanation
  641. ===========
  642. This function is defined as
  643. .. math ::
  644. H_\nu^{(1)} = J_\nu(z) + iY_\nu(z),
  645. where $J_\nu(z)$ is the Bessel function of the first kind, and
  646. $Y_\nu(z)$ is the Bessel function of the second kind.
  647. It is a solution to Bessel's equation.
  648. Examples
  649. ========
  650. >>> from sympy import hankel1
  651. >>> from sympy.abc import z, n
  652. >>> hankel1(n, z).diff(z)
  653. hankel1(n - 1, z)/2 - hankel1(n + 1, z)/2
  654. See Also
  655. ========
  656. hankel2, besselj, bessely
  657. References
  658. ==========
  659. .. [1] https://functions.wolfram.com/Bessel-TypeFunctions/HankelH1/
  660. """
  661. _a = S.One
  662. _b = S.One
  663. def _eval_conjugate(self):
  664. z = self.argument
  665. if z.is_extended_negative is False:
  666. return hankel2(self.order.conjugate(), z.conjugate())
  667. class hankel2(BesselBase):
  668. r"""
  669. Hankel function of the second kind.
  670. Explanation
  671. ===========
  672. This function is defined as
  673. .. math ::
  674. H_\nu^{(2)} = J_\nu(z) - iY_\nu(z),
  675. where $J_\nu(z)$ is the Bessel function of the first kind, and
  676. $Y_\nu(z)$ is the Bessel function of the second kind.
  677. It is a solution to Bessel's equation, and linearly independent from
  678. $H_\nu^{(1)}$.
  679. Examples
  680. ========
  681. >>> from sympy import hankel2
  682. >>> from sympy.abc import z, n
  683. >>> hankel2(n, z).diff(z)
  684. hankel2(n - 1, z)/2 - hankel2(n + 1, z)/2
  685. See Also
  686. ========
  687. hankel1, besselj, bessely
  688. References
  689. ==========
  690. .. [1] https://functions.wolfram.com/Bessel-TypeFunctions/HankelH2/
  691. """
  692. _a = S.One
  693. _b = S.One
  694. def _eval_conjugate(self):
  695. z = self.argument
  696. if z.is_extended_negative is False:
  697. return hankel1(self.order.conjugate(), z.conjugate())
  698. def assume_integer_order(fn):
  699. @wraps(fn)
  700. def g(self, nu, z):
  701. if nu.is_integer:
  702. return fn(self, nu, z)
  703. return g
  704. class SphericalBesselBase(BesselBase):
  705. """
  706. Base class for spherical Bessel functions.
  707. These are thin wrappers around ordinary Bessel functions,
  708. since spherical Bessel functions differ from the ordinary
  709. ones just by a slight change in order.
  710. To use this class, define the ``_eval_evalf()`` and ``_expand()`` methods.
  711. """
  712. def _expand(self, **hints):
  713. """ Expand self into a polynomial. Nu is guaranteed to be Integer. """
  714. raise NotImplementedError('expansion')
  715. def _eval_expand_func(self, **hints):
  716. if self.order.is_Integer:
  717. return self._expand(**hints)
  718. return self
  719. def fdiff(self, argindex=2):
  720. if argindex != 2:
  721. raise ArgumentIndexError(self, argindex)
  722. return self.__class__(self.order - 1, self.argument) - \
  723. self * (self.order + 1)/self.argument
  724. def _jn(n, z):
  725. return (spherical_bessel_fn(n, z)*sin(z) +
  726. S.NegativeOne**(n + 1)*spherical_bessel_fn(-n - 1, z)*cos(z))
  727. def _yn(n, z):
  728. # (-1)**(n + 1) * _jn(-n - 1, z)
  729. return (S.NegativeOne**(n + 1) * spherical_bessel_fn(-n - 1, z)*sin(z) -
  730. spherical_bessel_fn(n, z)*cos(z))
  731. class jn(SphericalBesselBase):
  732. r"""
  733. Spherical Bessel function of the first kind.
  734. Explanation
  735. ===========
  736. This function is a solution to the spherical Bessel equation
  737. .. math ::
  738. z^2 \frac{\mathrm{d}^2 w}{\mathrm{d}z^2}
  739. + 2z \frac{\mathrm{d}w}{\mathrm{d}z} + (z^2 - \nu(\nu + 1)) w = 0.
  740. It can be defined as
  741. .. math ::
  742. j_\nu(z) = \sqrt{\frac{\pi}{2z}} J_{\nu + \frac{1}{2}}(z),
  743. where $J_\nu(z)$ is the Bessel function of the first kind.
  744. The spherical Bessel functions of integral order are
  745. calculated using the formula:
  746. .. math:: j_n(z) = f_n(z) \sin{z} + (-1)^{n+1} f_{-n-1}(z) \cos{z},
  747. where the coefficients $f_n(z)$ are available as
  748. :func:`sympy.polys.orthopolys.spherical_bessel_fn`.
  749. Examples
  750. ========
  751. >>> from sympy import Symbol, jn, sin, cos, expand_func, besselj, bessely
  752. >>> z = Symbol("z")
  753. >>> nu = Symbol("nu", integer=True)
  754. >>> print(expand_func(jn(0, z)))
  755. sin(z)/z
  756. >>> expand_func(jn(1, z)) == sin(z)/z**2 - cos(z)/z
  757. True
  758. >>> expand_func(jn(3, z))
  759. (-6/z**2 + 15/z**4)*sin(z) + (1/z - 15/z**3)*cos(z)
  760. >>> jn(nu, z).rewrite(besselj)
  761. sqrt(2)*sqrt(pi)*sqrt(1/z)*besselj(nu + 1/2, z)/2
  762. >>> jn(nu, z).rewrite(bessely)
  763. (-1)**nu*sqrt(2)*sqrt(pi)*sqrt(1/z)*bessely(-nu - 1/2, z)/2
  764. >>> jn(2, 5.2+0.3j).evalf(20)
  765. 0.099419756723640344491 - 0.054525080242173562897*I
  766. See Also
  767. ========
  768. besselj, bessely, besselk, yn
  769. References
  770. ==========
  771. .. [1] https://dlmf.nist.gov/10.47
  772. """
  773. @classmethod
  774. def eval(cls, nu, z):
  775. if z.is_zero:
  776. if nu.is_zero:
  777. return S.One
  778. elif nu.is_integer:
  779. if nu.is_positive:
  780. return S.Zero
  781. else:
  782. return S.ComplexInfinity
  783. if z in (S.NegativeInfinity, S.Infinity):
  784. return S.Zero
  785. def _eval_rewrite_as_besselj(self, nu, z, **kwargs):
  786. return sqrt(pi/(2*z)) * besselj(nu + S.Half, z)
  787. def _eval_rewrite_as_bessely(self, nu, z, **kwargs):
  788. return S.NegativeOne**nu * sqrt(pi/(2*z)) * bessely(-nu - S.Half, z)
  789. def _eval_rewrite_as_yn(self, nu, z, **kwargs):
  790. return S.NegativeOne**(nu) * yn(-nu - 1, z)
  791. def _expand(self, **hints):
  792. return _jn(self.order, self.argument)
  793. def _eval_evalf(self, prec):
  794. if self.order.is_Integer:
  795. return self.rewrite(besselj)._eval_evalf(prec)
  796. class yn(SphericalBesselBase):
  797. r"""
  798. Spherical Bessel function of the second kind.
  799. Explanation
  800. ===========
  801. This function is another solution to the spherical Bessel equation, and
  802. linearly independent from $j_n$. It can be defined as
  803. .. math ::
  804. y_\nu(z) = \sqrt{\frac{\pi}{2z}} Y_{\nu + \frac{1}{2}}(z),
  805. where $Y_\nu(z)$ is the Bessel function of the second kind.
  806. For integral orders $n$, $y_n$ is calculated using the formula:
  807. .. math:: y_n(z) = (-1)^{n+1} j_{-n-1}(z)
  808. Examples
  809. ========
  810. >>> from sympy import Symbol, yn, sin, cos, expand_func, besselj, bessely
  811. >>> z = Symbol("z")
  812. >>> nu = Symbol("nu", integer=True)
  813. >>> print(expand_func(yn(0, z)))
  814. -cos(z)/z
  815. >>> expand_func(yn(1, z)) == -cos(z)/z**2-sin(z)/z
  816. True
  817. >>> yn(nu, z).rewrite(besselj)
  818. (-1)**(nu + 1)*sqrt(2)*sqrt(pi)*sqrt(1/z)*besselj(-nu - 1/2, z)/2
  819. >>> yn(nu, z).rewrite(bessely)
  820. sqrt(2)*sqrt(pi)*sqrt(1/z)*bessely(nu + 1/2, z)/2
  821. >>> yn(2, 5.2+0.3j).evalf(20)
  822. 0.18525034196069722536 + 0.014895573969924817587*I
  823. See Also
  824. ========
  825. besselj, bessely, besselk, jn
  826. References
  827. ==========
  828. .. [1] https://dlmf.nist.gov/10.47
  829. """
  830. @assume_integer_order
  831. def _eval_rewrite_as_besselj(self, nu, z, **kwargs):
  832. return S.NegativeOne**(nu+1) * sqrt(pi/(2*z)) * besselj(-nu - S.Half, z)
  833. @assume_integer_order
  834. def _eval_rewrite_as_bessely(self, nu, z, **kwargs):
  835. return sqrt(pi/(2*z)) * bessely(nu + S.Half, z)
  836. def _eval_rewrite_as_jn(self, nu, z, **kwargs):
  837. return S.NegativeOne**(nu + 1) * jn(-nu - 1, z)
  838. def _expand(self, **hints):
  839. return _yn(self.order, self.argument)
  840. def _eval_evalf(self, prec):
  841. if self.order.is_Integer:
  842. return self.rewrite(bessely)._eval_evalf(prec)
  843. class SphericalHankelBase(SphericalBesselBase):
  844. @assume_integer_order
  845. def _eval_rewrite_as_besselj(self, nu, z, **kwargs):
  846. # jn +- I*yn
  847. # jn as beeselj: sqrt(pi/(2*z)) * besselj(nu + S.Half, z)
  848. # yn as besselj: (-1)**(nu+1) * sqrt(pi/(2*z)) * besselj(-nu - S.Half, z)
  849. hks = self._hankel_kind_sign
  850. return sqrt(pi/(2*z))*(besselj(nu + S.Half, z) +
  851. hks*I*S.NegativeOne**(nu+1)*besselj(-nu - S.Half, z))
  852. @assume_integer_order
  853. def _eval_rewrite_as_bessely(self, nu, z, **kwargs):
  854. # jn +- I*yn
  855. # jn as bessely: (-1)**nu * sqrt(pi/(2*z)) * bessely(-nu - S.Half, z)
  856. # yn as bessely: sqrt(pi/(2*z)) * bessely(nu + S.Half, z)
  857. hks = self._hankel_kind_sign
  858. return sqrt(pi/(2*z))*(S.NegativeOne**nu*bessely(-nu - S.Half, z) +
  859. hks*I*bessely(nu + S.Half, z))
  860. def _eval_rewrite_as_yn(self, nu, z, **kwargs):
  861. hks = self._hankel_kind_sign
  862. return jn(nu, z).rewrite(yn) + hks*I*yn(nu, z)
  863. def _eval_rewrite_as_jn(self, nu, z, **kwargs):
  864. hks = self._hankel_kind_sign
  865. return jn(nu, z) + hks*I*yn(nu, z).rewrite(jn)
  866. def _eval_expand_func(self, **hints):
  867. if self.order.is_Integer:
  868. return self._expand(**hints)
  869. else:
  870. nu = self.order
  871. z = self.argument
  872. hks = self._hankel_kind_sign
  873. return jn(nu, z) + hks*I*yn(nu, z)
  874. def _expand(self, **hints):
  875. n = self.order
  876. z = self.argument
  877. hks = self._hankel_kind_sign
  878. # fully expanded version
  879. # return ((fn(n, z) * sin(z) +
  880. # (-1)**(n + 1) * fn(-n - 1, z) * cos(z)) + # jn
  881. # (hks * I * (-1)**(n + 1) *
  882. # (fn(-n - 1, z) * hk * I * sin(z) +
  883. # (-1)**(-n) * fn(n, z) * I * cos(z))) # +-I*yn
  884. # )
  885. return (_jn(n, z) + hks*I*_yn(n, z)).expand()
  886. def _eval_evalf(self, prec):
  887. if self.order.is_Integer:
  888. return self.rewrite(besselj)._eval_evalf(prec)
  889. class hn1(SphericalHankelBase):
  890. r"""
  891. Spherical Hankel function of the first kind.
  892. Explanation
  893. ===========
  894. This function is defined as
  895. .. math:: h_\nu^(1)(z) = j_\nu(z) + i y_\nu(z),
  896. where $j_\nu(z)$ and $y_\nu(z)$ are the spherical
  897. Bessel function of the first and second kinds.
  898. For integral orders $n$, $h_n^(1)$ is calculated using the formula:
  899. .. math:: h_n^(1)(z) = j_{n}(z) + i (-1)^{n+1} j_{-n-1}(z)
  900. Examples
  901. ========
  902. >>> from sympy import Symbol, hn1, hankel1, expand_func, yn, jn
  903. >>> z = Symbol("z")
  904. >>> nu = Symbol("nu", integer=True)
  905. >>> print(expand_func(hn1(nu, z)))
  906. jn(nu, z) + I*yn(nu, z)
  907. >>> print(expand_func(hn1(0, z)))
  908. sin(z)/z - I*cos(z)/z
  909. >>> print(expand_func(hn1(1, z)))
  910. -I*sin(z)/z - cos(z)/z + sin(z)/z**2 - I*cos(z)/z**2
  911. >>> hn1(nu, z).rewrite(jn)
  912. (-1)**(nu + 1)*I*jn(-nu - 1, z) + jn(nu, z)
  913. >>> hn1(nu, z).rewrite(yn)
  914. (-1)**nu*yn(-nu - 1, z) + I*yn(nu, z)
  915. >>> hn1(nu, z).rewrite(hankel1)
  916. sqrt(2)*sqrt(pi)*sqrt(1/z)*hankel1(nu, z)/2
  917. See Also
  918. ========
  919. hn2, jn, yn, hankel1, hankel2
  920. References
  921. ==========
  922. .. [1] https://dlmf.nist.gov/10.47
  923. """
  924. _hankel_kind_sign = S.One
  925. @assume_integer_order
  926. def _eval_rewrite_as_hankel1(self, nu, z, **kwargs):
  927. return sqrt(pi/(2*z))*hankel1(nu, z)
  928. class hn2(SphericalHankelBase):
  929. r"""
  930. Spherical Hankel function of the second kind.
  931. Explanation
  932. ===========
  933. This function is defined as
  934. .. math:: h_\nu^(2)(z) = j_\nu(z) - i y_\nu(z),
  935. where $j_\nu(z)$ and $y_\nu(z)$ are the spherical
  936. Bessel function of the first and second kinds.
  937. For integral orders $n$, $h_n^(2)$ is calculated using the formula:
  938. .. math:: h_n^(2)(z) = j_{n} - i (-1)^{n+1} j_{-n-1}(z)
  939. Examples
  940. ========
  941. >>> from sympy import Symbol, hn2, hankel2, expand_func, jn, yn
  942. >>> z = Symbol("z")
  943. >>> nu = Symbol("nu", integer=True)
  944. >>> print(expand_func(hn2(nu, z)))
  945. jn(nu, z) - I*yn(nu, z)
  946. >>> print(expand_func(hn2(0, z)))
  947. sin(z)/z + I*cos(z)/z
  948. >>> print(expand_func(hn2(1, z)))
  949. I*sin(z)/z - cos(z)/z + sin(z)/z**2 + I*cos(z)/z**2
  950. >>> hn2(nu, z).rewrite(hankel2)
  951. sqrt(2)*sqrt(pi)*sqrt(1/z)*hankel2(nu, z)/2
  952. >>> hn2(nu, z).rewrite(jn)
  953. -(-1)**(nu + 1)*I*jn(-nu - 1, z) + jn(nu, z)
  954. >>> hn2(nu, z).rewrite(yn)
  955. (-1)**nu*yn(-nu - 1, z) - I*yn(nu, z)
  956. See Also
  957. ========
  958. hn1, jn, yn, hankel1, hankel2
  959. References
  960. ==========
  961. .. [1] https://dlmf.nist.gov/10.47
  962. """
  963. _hankel_kind_sign = -S.One
  964. @assume_integer_order
  965. def _eval_rewrite_as_hankel2(self, nu, z, **kwargs):
  966. return sqrt(pi/(2*z))*hankel2(nu, z)
  967. def jn_zeros(n, k, method="sympy", dps=15):
  968. """
  969. Zeros of the spherical Bessel function of the first kind.
  970. Explanation
  971. ===========
  972. This returns an array of zeros of $jn$ up to the $k$-th zero.
  973. * method = "sympy": uses `mpmath.besseljzero
  974. <https://mpmath.org/doc/current/functions/bessel.html#mpmath.besseljzero>`_
  975. * method = "scipy": uses the
  976. `SciPy's sph_jn <https://docs.scipy.org/doc/scipy/reference/generated/scipy.special.jn_zeros.html>`_
  977. and
  978. `newton <https://docs.scipy.org/doc/scipy/reference/generated/scipy.optimize.newton.html>`_
  979. to find all
  980. roots, which is faster than computing the zeros using a general
  981. numerical solver, but it requires SciPy and only works with low
  982. precision floating point numbers. (The function used with
  983. method="sympy" is a recent addition to mpmath; before that a general
  984. solver was used.)
  985. Examples
  986. ========
  987. >>> from sympy import jn_zeros
  988. >>> jn_zeros(2, 4, dps=5)
  989. [5.7635, 9.095, 12.323, 15.515]
  990. See Also
  991. ========
  992. jn, yn, besselj, besselk, bessely
  993. Parameters
  994. ==========
  995. n : integer
  996. order of Bessel function
  997. k : integer
  998. number of zeros to return
  999. """
  1000. from math import pi as math_pi
  1001. if method == "sympy":
  1002. from mpmath import besseljzero
  1003. from mpmath.libmp.libmpf import dps_to_prec
  1004. prec = dps_to_prec(dps)
  1005. return [Expr._from_mpmath(besseljzero(S(n + 0.5)._to_mpmath(prec),
  1006. int(l)), prec)
  1007. for l in range(1, k + 1)]
  1008. elif method == "scipy":
  1009. from scipy.optimize import newton
  1010. try:
  1011. from scipy.special import spherical_jn
  1012. f = lambda x: spherical_jn(n, x)
  1013. except ImportError:
  1014. from scipy.special import sph_jn
  1015. f = lambda x: sph_jn(n, x)[0][-1]
  1016. else:
  1017. raise NotImplementedError("Unknown method.")
  1018. def solver(f, x):
  1019. if method == "scipy":
  1020. root = newton(f, x)
  1021. else:
  1022. raise NotImplementedError("Unknown method.")
  1023. return root
  1024. # we need to approximate the position of the first root:
  1025. root = n + math_pi
  1026. # determine the first root exactly:
  1027. root = solver(f, root)
  1028. roots = [root]
  1029. for i in range(k - 1):
  1030. # estimate the position of the next root using the last root + pi:
  1031. root = solver(f, root + math_pi)
  1032. roots.append(root)
  1033. return roots
  1034. class AiryBase(DefinedFunction):
  1035. """
  1036. Abstract base class for Airy functions.
  1037. This class is meant to reduce code duplication.
  1038. """
  1039. def _eval_conjugate(self):
  1040. return self.func(self.args[0].conjugate())
  1041. def _eval_is_extended_real(self):
  1042. return self.args[0].is_extended_real
  1043. def as_real_imag(self, deep=True, **hints):
  1044. z = self.args[0]
  1045. zc = z.conjugate()
  1046. f = self.func
  1047. u = (f(z)+f(zc))/2
  1048. v = I*(f(zc)-f(z))/2
  1049. return u, v
  1050. def _eval_expand_complex(self, deep=True, **hints):
  1051. re_part, im_part = self.as_real_imag(deep=deep, **hints)
  1052. return re_part + im_part*I
  1053. class airyai(AiryBase):
  1054. r"""
  1055. The Airy function $\operatorname{Ai}$ of the first kind.
  1056. Explanation
  1057. ===========
  1058. The Airy function $\operatorname{Ai}(z)$ is defined to be the function
  1059. satisfying Airy's differential equation
  1060. .. math::
  1061. \frac{\mathrm{d}^2 w(z)}{\mathrm{d}z^2} - z w(z) = 0.
  1062. Equivalently, for real $z$
  1063. .. math::
  1064. \operatorname{Ai}(z) := \frac{1}{\pi}
  1065. \int_0^\infty \cos\left(\frac{t^3}{3} + z t\right) \mathrm{d}t.
  1066. Examples
  1067. ========
  1068. Create an Airy function object:
  1069. >>> from sympy import airyai
  1070. >>> from sympy.abc import z
  1071. >>> airyai(z)
  1072. airyai(z)
  1073. Several special values are known:
  1074. >>> airyai(0)
  1075. 3**(1/3)/(3*gamma(2/3))
  1076. >>> from sympy import oo
  1077. >>> airyai(oo)
  1078. 0
  1079. >>> airyai(-oo)
  1080. 0
  1081. The Airy function obeys the mirror symmetry:
  1082. >>> from sympy import conjugate
  1083. >>> conjugate(airyai(z))
  1084. airyai(conjugate(z))
  1085. Differentiation with respect to $z$ is supported:
  1086. >>> from sympy import diff
  1087. >>> diff(airyai(z), z)
  1088. airyaiprime(z)
  1089. >>> diff(airyai(z), z, 2)
  1090. z*airyai(z)
  1091. Series expansion is also supported:
  1092. >>> from sympy import series
  1093. >>> series(airyai(z), z, 0, 3)
  1094. 3**(5/6)*gamma(1/3)/(6*pi) - 3**(1/6)*z*gamma(2/3)/(2*pi) + O(z**3)
  1095. We can numerically evaluate the Airy function to arbitrary precision
  1096. on the whole complex plane:
  1097. >>> airyai(-2).evalf(50)
  1098. 0.22740742820168557599192443603787379946077222541710
  1099. Rewrite $\operatorname{Ai}(z)$ in terms of hypergeometric functions:
  1100. >>> from sympy import hyper
  1101. >>> airyai(z).rewrite(hyper)
  1102. -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))
  1103. See Also
  1104. ========
  1105. airybi: Airy function of the second kind.
  1106. airyaiprime: Derivative of the Airy function of the first kind.
  1107. airybiprime: Derivative of the Airy function of the second kind.
  1108. References
  1109. ==========
  1110. .. [1] https://en.wikipedia.org/wiki/Airy_function
  1111. .. [2] https://dlmf.nist.gov/9
  1112. .. [3] https://encyclopediaofmath.org/wiki/Airy_functions
  1113. .. [4] https://mathworld.wolfram.com/AiryFunctions.html
  1114. """
  1115. nargs = 1
  1116. unbranched = True
  1117. @classmethod
  1118. def eval(cls, arg):
  1119. if arg.is_Number:
  1120. if arg is S.NaN:
  1121. return S.NaN
  1122. elif arg is S.Infinity:
  1123. return S.Zero
  1124. elif arg is S.NegativeInfinity:
  1125. return S.Zero
  1126. elif arg.is_zero:
  1127. return S.One / (3**Rational(2, 3) * gamma(Rational(2, 3)))
  1128. if arg.is_zero:
  1129. return S.One / (3**Rational(2, 3) * gamma(Rational(2, 3)))
  1130. def fdiff(self, argindex=1):
  1131. if argindex == 1:
  1132. return airyaiprime(self.args[0])
  1133. else:
  1134. raise ArgumentIndexError(self, argindex)
  1135. @staticmethod
  1136. @cacheit
  1137. def taylor_term(n, x, *previous_terms):
  1138. if n < 0:
  1139. return S.Zero
  1140. else:
  1141. x = sympify(x)
  1142. if len(previous_terms) > 1:
  1143. p = previous_terms[-1]
  1144. return ((cbrt(3)*x)**(-n)*(cbrt(3)*x)**(n + 1)*sin(pi*(n*Rational(2, 3) + Rational(4, 3)))*factorial(n) *
  1145. 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)
  1146. else:
  1147. return (S.One/(3**Rational(2, 3)*pi) * gamma((n+S.One)/S(3)) * sin(Rational(2, 3)*pi*(n+S.One)) /
  1148. factorial(n) * (cbrt(3)*x)**n)
  1149. def _eval_rewrite_as_besselj(self, z, **kwargs):
  1150. ot = Rational(1, 3)
  1151. tt = Rational(2, 3)
  1152. a = Pow(-z, Rational(3, 2))
  1153. if re(z).is_negative:
  1154. return ot*sqrt(-z) * (besselj(-ot, tt*a) + besselj(ot, tt*a))
  1155. def _eval_rewrite_as_besseli(self, z, **kwargs):
  1156. ot = Rational(1, 3)
  1157. tt = Rational(2, 3)
  1158. a = Pow(z, Rational(3, 2))
  1159. if re(z).is_positive:
  1160. return ot*sqrt(z) * (besseli(-ot, tt*a) - besseli(ot, tt*a))
  1161. else:
  1162. return ot*(Pow(a, ot)*besseli(-ot, tt*a) - z*Pow(a, -ot)*besseli(ot, tt*a))
  1163. def _eval_rewrite_as_hyper(self, z, **kwargs):
  1164. pf1 = S.One / (3**Rational(2, 3)*gamma(Rational(2, 3)))
  1165. pf2 = z / (root(3, 3)*gamma(Rational(1, 3)))
  1166. return pf1 * hyper([], [Rational(2, 3)], z**3/9) - pf2 * hyper([], [Rational(4, 3)], z**3/9)
  1167. def _eval_expand_func(self, **hints):
  1168. arg = self.args[0]
  1169. symbs = arg.free_symbols
  1170. if len(symbs) == 1:
  1171. z = symbs.pop()
  1172. c = Wild("c", exclude=[z])
  1173. d = Wild("d", exclude=[z])
  1174. m = Wild("m", exclude=[z])
  1175. n = Wild("n", exclude=[z])
  1176. M = arg.match(c*(d*z**n)**m)
  1177. if M is not None:
  1178. m = M[m]
  1179. # The transformation is given by 03.05.16.0001.01
  1180. # https://functions.wolfram.com/Bessel-TypeFunctions/AiryAi/16/01/01/0001/
  1181. if (3*m).is_integer:
  1182. c = M[c]
  1183. d = M[d]
  1184. n = M[n]
  1185. pf = (d * z**n)**m / (d**m * z**(m*n))
  1186. newarg = c * d**m * z**(m*n)
  1187. return S.Half * ((pf + S.One)*airyai(newarg) - (pf - S.One)/sqrt(3)*airybi(newarg))
  1188. class airybi(AiryBase):
  1189. r"""
  1190. The Airy function $\operatorname{Bi}$ of the second kind.
  1191. Explanation
  1192. ===========
  1193. The Airy function $\operatorname{Bi}(z)$ is defined to be the function
  1194. satisfying Airy's differential equation
  1195. .. math::
  1196. \frac{\mathrm{d}^2 w(z)}{\mathrm{d}z^2} - z w(z) = 0.
  1197. Equivalently, for real $z$
  1198. .. math::
  1199. \operatorname{Bi}(z) := \frac{1}{\pi}
  1200. \int_0^\infty
  1201. \exp\left(-\frac{t^3}{3} + z t\right)
  1202. + \sin\left(\frac{t^3}{3} + z t\right) \mathrm{d}t.
  1203. Examples
  1204. ========
  1205. Create an Airy function object:
  1206. >>> from sympy import airybi
  1207. >>> from sympy.abc import z
  1208. >>> airybi(z)
  1209. airybi(z)
  1210. Several special values are known:
  1211. >>> airybi(0)
  1212. 3**(5/6)/(3*gamma(2/3))
  1213. >>> from sympy import oo
  1214. >>> airybi(oo)
  1215. oo
  1216. >>> airybi(-oo)
  1217. 0
  1218. The Airy function obeys the mirror symmetry:
  1219. >>> from sympy import conjugate
  1220. >>> conjugate(airybi(z))
  1221. airybi(conjugate(z))
  1222. Differentiation with respect to $z$ is supported:
  1223. >>> from sympy import diff
  1224. >>> diff(airybi(z), z)
  1225. airybiprime(z)
  1226. >>> diff(airybi(z), z, 2)
  1227. z*airybi(z)
  1228. Series expansion is also supported:
  1229. >>> from sympy import series
  1230. >>> series(airybi(z), z, 0, 3)
  1231. 3**(1/3)*gamma(1/3)/(2*pi) + 3**(2/3)*z*gamma(2/3)/(2*pi) + O(z**3)
  1232. We can numerically evaluate the Airy function to arbitrary precision
  1233. on the whole complex plane:
  1234. >>> airybi(-2).evalf(50)
  1235. -0.41230258795639848808323405461146104203453483447240
  1236. Rewrite $\operatorname{Bi}(z)$ in terms of hypergeometric functions:
  1237. >>> from sympy import hyper
  1238. >>> airybi(z).rewrite(hyper)
  1239. 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))
  1240. See Also
  1241. ========
  1242. airyai: Airy function of the first kind.
  1243. airyaiprime: Derivative of the Airy function of the first kind.
  1244. airybiprime: Derivative of the Airy function of the second kind.
  1245. References
  1246. ==========
  1247. .. [1] https://en.wikipedia.org/wiki/Airy_function
  1248. .. [2] https://dlmf.nist.gov/9
  1249. .. [3] https://encyclopediaofmath.org/wiki/Airy_functions
  1250. .. [4] https://mathworld.wolfram.com/AiryFunctions.html
  1251. """
  1252. nargs = 1
  1253. unbranched = True
  1254. @classmethod
  1255. def eval(cls, arg):
  1256. if arg.is_Number:
  1257. if arg is S.NaN:
  1258. return S.NaN
  1259. elif arg is S.Infinity:
  1260. return S.Infinity
  1261. elif arg is S.NegativeInfinity:
  1262. return S.Zero
  1263. elif arg.is_zero:
  1264. return S.One / (3**Rational(1, 6) * gamma(Rational(2, 3)))
  1265. if arg.is_zero:
  1266. return S.One / (3**Rational(1, 6) * gamma(Rational(2, 3)))
  1267. def fdiff(self, argindex=1):
  1268. if argindex == 1:
  1269. return airybiprime(self.args[0])
  1270. else:
  1271. raise ArgumentIndexError(self, argindex)
  1272. @staticmethod
  1273. @cacheit
  1274. def taylor_term(n, x, *previous_terms):
  1275. if n < 0:
  1276. return S.Zero
  1277. else:
  1278. x = sympify(x)
  1279. if len(previous_terms) > 1:
  1280. p = previous_terms[-1]
  1281. return (cbrt(3)*x * Abs(sin(Rational(2, 3)*pi*(n + S.One))) * factorial((n - S.One)/S(3)) /
  1282. ((n + S.One) * Abs(cos(Rational(2, 3)*pi*(n + S.Half))) * factorial((n - 2)/S(3))) * p)
  1283. else:
  1284. return (S.One/(root(3, 6)*pi) * gamma((n + S.One)/S(3)) * Abs(sin(Rational(2, 3)*pi*(n + S.One))) /
  1285. factorial(n) * (cbrt(3)*x)**n)
  1286. def _eval_rewrite_as_besselj(self, z, **kwargs):
  1287. ot = Rational(1, 3)
  1288. tt = Rational(2, 3)
  1289. a = Pow(-z, Rational(3, 2))
  1290. if re(z).is_negative:
  1291. return sqrt(-z/3) * (besselj(-ot, tt*a) - besselj(ot, tt*a))
  1292. def _eval_rewrite_as_besseli(self, z, **kwargs):
  1293. ot = Rational(1, 3)
  1294. tt = Rational(2, 3)
  1295. a = Pow(z, Rational(3, 2))
  1296. if re(z).is_positive:
  1297. return sqrt(z)/sqrt(3) * (besseli(-ot, tt*a) + besseli(ot, tt*a))
  1298. else:
  1299. b = Pow(a, ot)
  1300. c = Pow(a, -ot)
  1301. return sqrt(ot)*(b*besseli(-ot, tt*a) + z*c*besseli(ot, tt*a))
  1302. def _eval_rewrite_as_hyper(self, z, **kwargs):
  1303. pf1 = S.One / (root(3, 6)*gamma(Rational(2, 3)))
  1304. pf2 = z*root(3, 6) / gamma(Rational(1, 3))
  1305. return pf1 * hyper([], [Rational(2, 3)], z**3/9) + pf2 * hyper([], [Rational(4, 3)], z**3/9)
  1306. def _eval_expand_func(self, **hints):
  1307. arg = self.args[0]
  1308. symbs = arg.free_symbols
  1309. if len(symbs) == 1:
  1310. z = symbs.pop()
  1311. c = Wild("c", exclude=[z])
  1312. d = Wild("d", exclude=[z])
  1313. m = Wild("m", exclude=[z])
  1314. n = Wild("n", exclude=[z])
  1315. M = arg.match(c*(d*z**n)**m)
  1316. if M is not None:
  1317. m = M[m]
  1318. # The transformation is given by 03.06.16.0001.01
  1319. # https://functions.wolfram.com/Bessel-TypeFunctions/AiryBi/16/01/01/0001/
  1320. if (3*m).is_integer:
  1321. c = M[c]
  1322. d = M[d]
  1323. n = M[n]
  1324. pf = (d * z**n)**m / (d**m * z**(m*n))
  1325. newarg = c * d**m * z**(m*n)
  1326. return S.Half * (sqrt(3)*(S.One - pf)*airyai(newarg) + (S.One + pf)*airybi(newarg))
  1327. class airyaiprime(AiryBase):
  1328. r"""
  1329. The derivative $\operatorname{Ai}^\prime$ of the Airy function of the first
  1330. kind.
  1331. Explanation
  1332. ===========
  1333. The Airy function $\operatorname{Ai}^\prime(z)$ is defined to be the
  1334. function
  1335. .. math::
  1336. \operatorname{Ai}^\prime(z) := \frac{\mathrm{d} \operatorname{Ai}(z)}{\mathrm{d} z}.
  1337. Examples
  1338. ========
  1339. Create an Airy function object:
  1340. >>> from sympy import airyaiprime
  1341. >>> from sympy.abc import z
  1342. >>> airyaiprime(z)
  1343. airyaiprime(z)
  1344. Several special values are known:
  1345. >>> airyaiprime(0)
  1346. -3**(2/3)/(3*gamma(1/3))
  1347. >>> from sympy import oo
  1348. >>> airyaiprime(oo)
  1349. 0
  1350. The Airy function obeys the mirror symmetry:
  1351. >>> from sympy import conjugate
  1352. >>> conjugate(airyaiprime(z))
  1353. airyaiprime(conjugate(z))
  1354. Differentiation with respect to $z$ is supported:
  1355. >>> from sympy import diff
  1356. >>> diff(airyaiprime(z), z)
  1357. z*airyai(z)
  1358. >>> diff(airyaiprime(z), z, 2)
  1359. z*airyaiprime(z) + airyai(z)
  1360. Series expansion is also supported:
  1361. >>> from sympy import series
  1362. >>> series(airyaiprime(z), z, 0, 3)
  1363. -3**(2/3)/(3*gamma(1/3)) + 3**(1/3)*z**2/(6*gamma(2/3)) + O(z**3)
  1364. We can numerically evaluate the Airy function to arbitrary precision
  1365. on the whole complex plane:
  1366. >>> airyaiprime(-2).evalf(50)
  1367. 0.61825902074169104140626429133247528291577794512415
  1368. Rewrite $\operatorname{Ai}^\prime(z)$ in terms of hypergeometric functions:
  1369. >>> from sympy import hyper
  1370. >>> airyaiprime(z).rewrite(hyper)
  1371. 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))
  1372. See Also
  1373. ========
  1374. airyai: Airy function of the first kind.
  1375. airybi: Airy function of the second kind.
  1376. airybiprime: Derivative of the Airy function of the second kind.
  1377. References
  1378. ==========
  1379. .. [1] https://en.wikipedia.org/wiki/Airy_function
  1380. .. [2] https://dlmf.nist.gov/9
  1381. .. [3] https://encyclopediaofmath.org/wiki/Airy_functions
  1382. .. [4] https://mathworld.wolfram.com/AiryFunctions.html
  1383. """
  1384. nargs = 1
  1385. unbranched = True
  1386. @classmethod
  1387. def eval(cls, arg):
  1388. if arg.is_Number:
  1389. if arg is S.NaN:
  1390. return S.NaN
  1391. elif arg is S.Infinity:
  1392. return S.Zero
  1393. if arg.is_zero:
  1394. return S.NegativeOne / (3**Rational(1, 3) * gamma(Rational(1, 3)))
  1395. def fdiff(self, argindex=1):
  1396. if argindex == 1:
  1397. return self.args[0]*airyai(self.args[0])
  1398. else:
  1399. raise ArgumentIndexError(self, argindex)
  1400. def _eval_evalf(self, prec):
  1401. z = self.args[0]._to_mpmath(prec)
  1402. with workprec(prec):
  1403. res = mp.airyai(z, derivative=1)
  1404. return Expr._from_mpmath(res, prec)
  1405. def _eval_rewrite_as_besselj(self, z, **kwargs):
  1406. tt = Rational(2, 3)
  1407. a = Pow(-z, Rational(3, 2))
  1408. if re(z).is_negative:
  1409. return z/3 * (besselj(-tt, tt*a) - besselj(tt, tt*a))
  1410. def _eval_rewrite_as_besseli(self, z, **kwargs):
  1411. ot = Rational(1, 3)
  1412. tt = Rational(2, 3)
  1413. a = tt * Pow(z, Rational(3, 2))
  1414. if re(z).is_positive:
  1415. return z/3 * (besseli(tt, a) - besseli(-tt, a))
  1416. else:
  1417. a = Pow(z, Rational(3, 2))
  1418. b = Pow(a, tt)
  1419. c = Pow(a, -tt)
  1420. return ot * (z**2*c*besseli(tt, tt*a) - b*besseli(-ot, tt*a))
  1421. def _eval_rewrite_as_hyper(self, z, **kwargs):
  1422. pf1 = z**2 / (2*3**Rational(2, 3)*gamma(Rational(2, 3)))
  1423. pf2 = 1 / (root(3, 3)*gamma(Rational(1, 3)))
  1424. return pf1 * hyper([], [Rational(5, 3)], z**3/9) - pf2 * hyper([], [Rational(1, 3)], z**3/9)
  1425. def _eval_expand_func(self, **hints):
  1426. arg = self.args[0]
  1427. symbs = arg.free_symbols
  1428. if len(symbs) == 1:
  1429. z = symbs.pop()
  1430. c = Wild("c", exclude=[z])
  1431. d = Wild("d", exclude=[z])
  1432. m = Wild("m", exclude=[z])
  1433. n = Wild("n", exclude=[z])
  1434. M = arg.match(c*(d*z**n)**m)
  1435. if M is not None:
  1436. m = M[m]
  1437. # The transformation is in principle
  1438. # given by 03.07.16.0001.01 but note
  1439. # that there is an error in this formula.
  1440. # https://functions.wolfram.com/Bessel-TypeFunctions/AiryAiPrime/16/01/01/0001/
  1441. if (3*m).is_integer:
  1442. c = M[c]
  1443. d = M[d]
  1444. n = M[n]
  1445. pf = (d**m * z**(n*m)) / (d * z**n)**m
  1446. newarg = c * d**m * z**(n*m)
  1447. return S.Half * ((pf + S.One)*airyaiprime(newarg) + (pf - S.One)/sqrt(3)*airybiprime(newarg))
  1448. class airybiprime(AiryBase):
  1449. r"""
  1450. The derivative $\operatorname{Bi}^\prime$ of the Airy function of the first
  1451. kind.
  1452. Explanation
  1453. ===========
  1454. The Airy function $\operatorname{Bi}^\prime(z)$ is defined to be the
  1455. function
  1456. .. math::
  1457. \operatorname{Bi}^\prime(z) := \frac{\mathrm{d} \operatorname{Bi}(z)}{\mathrm{d} z}.
  1458. Examples
  1459. ========
  1460. Create an Airy function object:
  1461. >>> from sympy import airybiprime
  1462. >>> from sympy.abc import z
  1463. >>> airybiprime(z)
  1464. airybiprime(z)
  1465. Several special values are known:
  1466. >>> airybiprime(0)
  1467. 3**(1/6)/gamma(1/3)
  1468. >>> from sympy import oo
  1469. >>> airybiprime(oo)
  1470. oo
  1471. >>> airybiprime(-oo)
  1472. 0
  1473. The Airy function obeys the mirror symmetry:
  1474. >>> from sympy import conjugate
  1475. >>> conjugate(airybiprime(z))
  1476. airybiprime(conjugate(z))
  1477. Differentiation with respect to $z$ is supported:
  1478. >>> from sympy import diff
  1479. >>> diff(airybiprime(z), z)
  1480. z*airybi(z)
  1481. >>> diff(airybiprime(z), z, 2)
  1482. z*airybiprime(z) + airybi(z)
  1483. Series expansion is also supported:
  1484. >>> from sympy import series
  1485. >>> series(airybiprime(z), z, 0, 3)
  1486. 3**(1/6)/gamma(1/3) + 3**(5/6)*z**2/(6*gamma(2/3)) + O(z**3)
  1487. We can numerically evaluate the Airy function to arbitrary precision
  1488. on the whole complex plane:
  1489. >>> airybiprime(-2).evalf(50)
  1490. 0.27879516692116952268509756941098324140300059345163
  1491. Rewrite $\operatorname{Bi}^\prime(z)$ in terms of hypergeometric functions:
  1492. >>> from sympy import hyper
  1493. >>> airybiprime(z).rewrite(hyper)
  1494. 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)
  1495. See Also
  1496. ========
  1497. airyai: Airy function of the first kind.
  1498. airybi: Airy function of the second kind.
  1499. airyaiprime: Derivative of the Airy function of the first kind.
  1500. References
  1501. ==========
  1502. .. [1] https://en.wikipedia.org/wiki/Airy_function
  1503. .. [2] https://dlmf.nist.gov/9
  1504. .. [3] https://encyclopediaofmath.org/wiki/Airy_functions
  1505. .. [4] https://mathworld.wolfram.com/AiryFunctions.html
  1506. """
  1507. nargs = 1
  1508. unbranched = True
  1509. @classmethod
  1510. def eval(cls, arg):
  1511. if arg.is_Number:
  1512. if arg is S.NaN:
  1513. return S.NaN
  1514. elif arg is S.Infinity:
  1515. return S.Infinity
  1516. elif arg is S.NegativeInfinity:
  1517. return S.Zero
  1518. elif arg.is_zero:
  1519. return 3**Rational(1, 6) / gamma(Rational(1, 3))
  1520. if arg.is_zero:
  1521. return 3**Rational(1, 6) / gamma(Rational(1, 3))
  1522. def fdiff(self, argindex=1):
  1523. if argindex == 1:
  1524. return self.args[0]*airybi(self.args[0])
  1525. else:
  1526. raise ArgumentIndexError(self, argindex)
  1527. def _eval_evalf(self, prec):
  1528. z = self.args[0]._to_mpmath(prec)
  1529. with workprec(prec):
  1530. res = mp.airybi(z, derivative=1)
  1531. return Expr._from_mpmath(res, prec)
  1532. def _eval_rewrite_as_besselj(self, z, **kwargs):
  1533. tt = Rational(2, 3)
  1534. a = tt * Pow(-z, Rational(3, 2))
  1535. if re(z).is_negative:
  1536. return -z/sqrt(3) * (besselj(-tt, a) + besselj(tt, a))
  1537. def _eval_rewrite_as_besseli(self, z, **kwargs):
  1538. ot = Rational(1, 3)
  1539. tt = Rational(2, 3)
  1540. a = tt * Pow(z, Rational(3, 2))
  1541. if re(z).is_positive:
  1542. return z/sqrt(3) * (besseli(-tt, a) + besseli(tt, a))
  1543. else:
  1544. a = Pow(z, Rational(3, 2))
  1545. b = Pow(a, tt)
  1546. c = Pow(a, -tt)
  1547. return sqrt(ot) * (b*besseli(-tt, tt*a) + z**2*c*besseli(tt, tt*a))
  1548. def _eval_rewrite_as_hyper(self, z, **kwargs):
  1549. pf1 = z**2 / (2*root(3, 6)*gamma(Rational(2, 3)))
  1550. pf2 = root(3, 6) / gamma(Rational(1, 3))
  1551. return pf1 * hyper([], [Rational(5, 3)], z**3/9) + pf2 * hyper([], [Rational(1, 3)], z**3/9)
  1552. def _eval_expand_func(self, **hints):
  1553. arg = self.args[0]
  1554. symbs = arg.free_symbols
  1555. if len(symbs) == 1:
  1556. z = symbs.pop()
  1557. c = Wild("c", exclude=[z])
  1558. d = Wild("d", exclude=[z])
  1559. m = Wild("m", exclude=[z])
  1560. n = Wild("n", exclude=[z])
  1561. M = arg.match(c*(d*z**n)**m)
  1562. if M is not None:
  1563. m = M[m]
  1564. # The transformation is in principle
  1565. # given by 03.08.16.0001.01 but note
  1566. # that there is an error in this formula.
  1567. # https://functions.wolfram.com/Bessel-TypeFunctions/AiryBiPrime/16/01/01/0001/
  1568. if (3*m).is_integer:
  1569. c = M[c]
  1570. d = M[d]
  1571. n = M[n]
  1572. pf = (d**m * z**(n*m)) / (d * z**n)**m
  1573. newarg = c * d**m * z**(n*m)
  1574. return S.Half * (sqrt(3)*(pf - S.One)*airyaiprime(newarg) + (pf + S.One)*airybiprime(newarg))
  1575. class marcumq(DefinedFunction):
  1576. r"""
  1577. The Marcum Q-function.
  1578. Explanation
  1579. ===========
  1580. The Marcum Q-function is defined by the meromorphic continuation of
  1581. .. math::
  1582. 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
  1583. Examples
  1584. ========
  1585. >>> from sympy import marcumq
  1586. >>> from sympy.abc import m, a, b
  1587. >>> marcumq(m, a, b)
  1588. marcumq(m, a, b)
  1589. Special values:
  1590. >>> marcumq(m, 0, b)
  1591. uppergamma(m, b**2/2)/gamma(m)
  1592. >>> marcumq(0, 0, 0)
  1593. 0
  1594. >>> marcumq(0, a, 0)
  1595. 1 - exp(-a**2/2)
  1596. >>> marcumq(1, a, a)
  1597. 1/2 + exp(-a**2)*besseli(0, a**2)/2
  1598. >>> marcumq(2, a, a)
  1599. 1/2 + exp(-a**2)*besseli(0, a**2)/2 + exp(-a**2)*besseli(1, a**2)
  1600. Differentiation with respect to $a$ and $b$ is supported:
  1601. >>> from sympy import diff
  1602. >>> diff(marcumq(m, a, b), a)
  1603. a*(-marcumq(m, a, b) + marcumq(m + 1, a, b))
  1604. >>> diff(marcumq(m, a, b), b)
  1605. -a**(1 - m)*b**m*exp(-a**2/2 - b**2/2)*besseli(m - 1, a*b)
  1606. References
  1607. ==========
  1608. .. [1] https://en.wikipedia.org/wiki/Marcum_Q-function
  1609. .. [2] https://mathworld.wolfram.com/MarcumQ-Function.html
  1610. """
  1611. @classmethod
  1612. def eval(cls, m, a, b):
  1613. if a is S.Zero:
  1614. if m is S.Zero and b is S.Zero:
  1615. return S.Zero
  1616. return uppergamma(m, b**2 * S.Half) / gamma(m)
  1617. if m is S.Zero and b is S.Zero:
  1618. return 1 - 1 / exp(a**2 * S.Half)
  1619. if a == b:
  1620. if m is S.One:
  1621. return (1 + exp(-a**2) * besseli(0, a**2))*S.Half
  1622. if m == 2:
  1623. return S.Half + S.Half * exp(-a**2) * besseli(0, a**2) + exp(-a**2) * besseli(1, a**2)
  1624. if a.is_zero:
  1625. if m.is_zero and b.is_zero:
  1626. return S.Zero
  1627. return uppergamma(m, b**2*S.Half) / gamma(m)
  1628. if m.is_zero and b.is_zero:
  1629. return 1 - 1 / exp(a**2*S.Half)
  1630. def fdiff(self, argindex=2):
  1631. m, a, b = self.args
  1632. if argindex == 2:
  1633. return a * (-marcumq(m, a, b) + marcumq(1+m, a, b))
  1634. elif argindex == 3:
  1635. return (-b**m / a**(m-1)) * exp(-(a**2 + b**2)/2) * besseli(m-1, a*b)
  1636. else:
  1637. raise ArgumentIndexError(self, argindex)
  1638. def _eval_rewrite_as_Integral(self, m, a, b, **kwargs):
  1639. from sympy.integrals.integrals import Integral
  1640. x = kwargs.get('x', Dummy(uniquely_named_symbol('x').name))
  1641. return a ** (1 - m) * \
  1642. Integral(x**m * exp(-(x**2 + a**2)/2) * besseli(m-1, a*x), [x, b, S.Infinity])
  1643. def _eval_rewrite_as_Sum(self, m, a, b, **kwargs):
  1644. from sympy.concrete.summations import Sum
  1645. k = kwargs.get('k', Dummy('k'))
  1646. return exp(-(a**2 + b**2) / 2) * Sum((a/b)**k * besseli(k, a*b), [k, 1-m, S.Infinity])
  1647. def _eval_rewrite_as_besseli(self, m, a, b, **kwargs):
  1648. if a == b:
  1649. if m == 1:
  1650. return (1 + exp(-a**2) * besseli(0, a**2)) / 2
  1651. if m.is_Integer and m >= 2:
  1652. s = sum(besseli(i, a**2) for i in range(1, m))
  1653. return S.Half + exp(-a**2) * besseli(0, a**2) / 2 + exp(-a**2) * s
  1654. def _eval_is_zero(self):
  1655. if all(arg.is_zero for arg in self.args):
  1656. return True
  1657. class _besseli(DefinedFunction):
  1658. """
  1659. Helper function to make the $\\mathrm{besseli}(nu, z)$
  1660. function tractable for the Gruntz algorithm.
  1661. """
  1662. def _eval_aseries(self, n, args0, x, logx):
  1663. from sympy.functions.combinatorial.factorials import RisingFactorial
  1664. from sympy.series.order import Order
  1665. point = args0[1]
  1666. if point in [S.Infinity, S.NegativeInfinity]:
  1667. nu, z = self.args
  1668. l = [((RisingFactorial(Rational(2*nu - 1, 2), k)*RisingFactorial(
  1669. Rational(2*nu + 1, 2), k))/((2)**(k)*z**(Rational(2*k + 1, 2))*factorial(k))) for k in range(n)]
  1670. return sqrt(pi/(2))*(Add(*l)) + Order(1/z**(Rational(2*n + 1, 2)), x)
  1671. return super()._eval_aseries(n, args0, x, logx)
  1672. def _eval_rewrite_as_intractable(self, nu, z, **kwargs):
  1673. return exp(-z)*besseli(nu, z)
  1674. def _eval_nseries(self, x, n, logx, cdir=0):
  1675. x0 = self.args[0].limit(x, 0)
  1676. if x0.is_zero:
  1677. f = self._eval_rewrite_as_intractable(*self.args)
  1678. return f._eval_nseries(x, n, logx)
  1679. return super()._eval_nseries(x, n, logx)
  1680. class _besselk(DefinedFunction):
  1681. """
  1682. Helper function to make the $\\mathrm{besselk}(nu, z)$
  1683. function tractable for the Gruntz algorithm.
  1684. """
  1685. def _eval_aseries(self, n, args0, x, logx):
  1686. from sympy.functions.combinatorial.factorials import RisingFactorial
  1687. from sympy.series.order import Order
  1688. point = args0[1]
  1689. if point in [S.Infinity, S.NegativeInfinity]:
  1690. nu, z = self.args
  1691. l = [((RisingFactorial(Rational(2*nu - 1, 2), k)*RisingFactorial(
  1692. Rational(2*nu + 1, 2), k))/((-2)**(k)*z**(Rational(2*k + 1, 2))*factorial(k))) for k in range(n)]
  1693. return sqrt(pi/(2))*(Add(*l)) + Order(1/z**(Rational(2*n + 1, 2)), x)
  1694. return super()._eval_aseries(n, args0, x, logx)
  1695. def _eval_rewrite_as_intractable(self,nu, z, **kwargs):
  1696. return exp(z)*besselk(nu, z)
  1697. def _eval_nseries(self, x, n, logx, cdir=0):
  1698. x0 = self.args[0].limit(x, 0)
  1699. if x0.is_zero:
  1700. f = self._eval_rewrite_as_intractable(*self.args)
  1701. return f._eval_nseries(x, n, logx)
  1702. return super()._eval_nseries(x, n, logx)