bessel.py 56 KB

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