error_functions.py 78 KB

12345678910111213141516171819202122232425262728293031323334353637383940414243444546474849505152535455565758596061626364656667686970717273747576777879808182838485868788899091929394959697989910010110210310410510610710810911011111211311411511611711811912012112212312412512612712812913013113213313413513613713813914014114214314414514614714814915015115215315415515615715815916016116216316416516616716816917017117217317417517617717817918018118218318418518618718818919019119219319419519619719819920020120220320420520620720820921021121221321421521621721821922022122222322422522622722822923023123223323423523623723823924024124224324424524624724824925025125225325425525625725825926026126226326426526626726826927027127227327427527627727827928028128228328428528628728828929029129229329429529629729829930030130230330430530630730830931031131231331431531631731831932032132232332432532632732832933033133233333433533633733833934034134234334434534634734834935035135235335435535635735835936036136236336436536636736836937037137237337437537637737837938038138238338438538638738838939039139239339439539639739839940040140240340440540640740840941041141241341441541641741841942042142242342442542642742842943043143243343443543643743843944044144244344444544644744844945045145245345445545645745845946046146246346446546646746846947047147247347447547647747847948048148248348448548648748848949049149249349449549649749849950050150250350450550650750850951051151251351451551651751851952052152252352452552652752852953053153253353453553653753853954054154254354454554654754854955055155255355455555655755855956056156256356456556656756856957057157257357457557657757857958058158258358458558658758858959059159259359459559659759859960060160260360460560660760860961061161261361461561661761861962062162262362462562662762862963063163263363463563663763863964064164264364464564664764864965065165265365465565665765865966066166266366466566666766866967067167267367467567667767867968068168268368468568668768868969069169269369469569669769869970070170270370470570670770870971071171271371471571671771871972072172272372472572672772872973073173273373473573673773873974074174274374474574674774874975075175275375475575675775875976076176276376476576676776876977077177277377477577677777877978078178278378478578678778878979079179279379479579679779879980080180280380480580680780880981081181281381481581681781881982082182282382482582682782882983083183283383483583683783883984084184284384484584684784884985085185285385485585685785885986086186286386486586686786886987087187287387487587687787887988088188288388488588688788888989089189289389489589689789889990090190290390490590690790890991091191291391491591691791891992092192292392492592692792892993093193293393493593693793893994094194294394494594694794894995095195295395495595695795895996096196296396496596696796896997097197297397497597697797897998098198298398498598698798898999099199299399499599699799899910001001100210031004100510061007100810091010101110121013101410151016101710181019102010211022102310241025102610271028102910301031103210331034103510361037103810391040104110421043104410451046104710481049105010511052105310541055105610571058105910601061106210631064106510661067106810691070107110721073107410751076107710781079108010811082108310841085108610871088108910901091109210931094109510961097109810991100110111021103110411051106110711081109111011111112111311141115111611171118111911201121112211231124112511261127112811291130113111321133113411351136113711381139114011411142114311441145114611471148114911501151115211531154115511561157115811591160116111621163116411651166116711681169117011711172117311741175117611771178117911801181118211831184118511861187118811891190119111921193119411951196119711981199120012011202120312041205120612071208120912101211121212131214121512161217121812191220122112221223122412251226122712281229123012311232123312341235123612371238123912401241124212431244124512461247124812491250125112521253125412551256125712581259126012611262126312641265126612671268126912701271127212731274127512761277127812791280128112821283128412851286128712881289129012911292129312941295129612971298129913001301130213031304130513061307130813091310131113121313131413151316131713181319132013211322132313241325132613271328132913301331133213331334133513361337133813391340134113421343134413451346134713481349135013511352135313541355135613571358135913601361136213631364136513661367136813691370137113721373137413751376137713781379138013811382138313841385138613871388138913901391139213931394139513961397139813991400140114021403140414051406140714081409141014111412141314141415141614171418141914201421142214231424142514261427142814291430143114321433143414351436143714381439144014411442144314441445144614471448144914501451145214531454145514561457145814591460146114621463146414651466146714681469147014711472147314741475147614771478147914801481148214831484148514861487148814891490149114921493149414951496149714981499150015011502150315041505150615071508150915101511151215131514151515161517151815191520152115221523152415251526152715281529153015311532153315341535153615371538153915401541154215431544154515461547154815491550155115521553155415551556155715581559156015611562156315641565156615671568156915701571157215731574157515761577157815791580158115821583158415851586158715881589159015911592159315941595159615971598159916001601160216031604160516061607160816091610161116121613161416151616161716181619162016211622162316241625162616271628162916301631163216331634163516361637163816391640164116421643164416451646164716481649165016511652165316541655165616571658165916601661166216631664166516661667166816691670167116721673167416751676167716781679168016811682168316841685168616871688168916901691169216931694169516961697169816991700170117021703170417051706170717081709171017111712171317141715171617171718171917201721172217231724172517261727172817291730173117321733173417351736173717381739174017411742174317441745174617471748174917501751175217531754175517561757175817591760176117621763176417651766176717681769177017711772177317741775177617771778177917801781178217831784178517861787178817891790179117921793179417951796179717981799180018011802180318041805180618071808180918101811181218131814181518161817181818191820182118221823182418251826182718281829183018311832183318341835183618371838183918401841184218431844184518461847184818491850185118521853185418551856185718581859186018611862186318641865186618671868186918701871187218731874187518761877187818791880188118821883188418851886188718881889189018911892189318941895189618971898189919001901190219031904190519061907190819091910191119121913191419151916191719181919192019211922192319241925192619271928192919301931193219331934193519361937193819391940194119421943194419451946194719481949195019511952195319541955195619571958195919601961196219631964196519661967196819691970197119721973197419751976197719781979198019811982198319841985198619871988198919901991199219931994199519961997199819992000200120022003200420052006200720082009201020112012201320142015201620172018201920202021202220232024202520262027202820292030203120322033203420352036203720382039204020412042204320442045204620472048204920502051205220532054205520562057205820592060206120622063206420652066206720682069207020712072207320742075207620772078207920802081208220832084208520862087208820892090209120922093209420952096209720982099210021012102210321042105210621072108210921102111211221132114211521162117211821192120212121222123212421252126212721282129213021312132213321342135213621372138213921402141214221432144214521462147214821492150215121522153215421552156215721582159216021612162216321642165216621672168216921702171217221732174217521762177217821792180218121822183218421852186218721882189219021912192219321942195219621972198219922002201220222032204220522062207220822092210221122122213221422152216221722182219222022212222222322242225222622272228222922302231223222332234223522362237223822392240224122422243224422452246224722482249225022512252225322542255225622572258225922602261226222632264226522662267226822692270227122722273227422752276227722782279228022812282228322842285228622872288228922902291229222932294229522962297229822992300230123022303230423052306230723082309231023112312231323142315231623172318231923202321232223232324232523262327232823292330233123322333233423352336233723382339234023412342234323442345234623472348234923502351235223532354235523562357235823592360236123622363236423652366236723682369237023712372237323742375237623772378237923802381238223832384238523862387238823892390239123922393239423952396239723982399240024012402240324042405240624072408240924102411241224132414241524162417241824192420242124222423242424252426242724282429243024312432243324342435243624372438243924402441244224432444244524462447244824492450245124522453245424552456245724582459246024612462246324642465246624672468246924702471247224732474247524762477247824792480248124822483248424852486248724882489249024912492249324942495249624972498249925002501250225032504250525062507250825092510251125122513251425152516251725182519252025212522252325242525252625272528252925302531253225332534253525362537253825392540254125422543254425452546254725482549255025512552255325542555255625572558255925602561256225632564256525662567256825692570257125722573257425752576257725782579258025812582258325842585258625872588258925902591259225932594259525962597259825992600260126022603260426052606260726082609261026112612261326142615261626172618261926202621262226232624262526262627262826292630263126322633263426352636263726382639264026412642264326442645264626472648264926502651265226532654265526562657265826592660266126622663266426652666266726682669267026712672267326742675267626772678267926802681268226832684268526862687268826892690269126922693269426952696269726982699270027012702270327042705270627072708270927102711271227132714271527162717271827192720272127222723272427252726272727282729273027312732273327342735273627372738273927402741274227432744274527462747274827492750275127522753275427552756275727582759276027612762276327642765276627672768276927702771277227732774277527762777277827792780278127822783278427852786278727882789279027912792279327942795279627972798279928002801
  1. """ This module contains various functions that are special cases
  2. of incomplete gamma functions. It should probably be renamed. """
  3. from sympy.core import EulerGamma # Must be imported from core, not core.numbers
  4. from sympy.core.add import Add
  5. from sympy.core.cache import cacheit
  6. from sympy.core.function import DefinedFunction, ArgumentIndexError, expand_mul
  7. from sympy.core.logic import fuzzy_or
  8. from sympy.core.numbers import I, pi, Rational, Integer
  9. from sympy.core.relational import is_eq
  10. from sympy.core.power import Pow
  11. from sympy.core.singleton import S
  12. from sympy.core.symbol import Dummy, uniquely_named_symbol
  13. from sympy.core.sympify import sympify
  14. from sympy.functions.combinatorial.factorials import factorial, factorial2, RisingFactorial
  15. from sympy.functions.elementary.complexes import polar_lift, re, unpolarify
  16. from sympy.functions.elementary.integers import ceiling, floor
  17. from sympy.functions.elementary.miscellaneous import sqrt, root
  18. from sympy.functions.elementary.exponential import exp, log, exp_polar
  19. from sympy.functions.elementary.hyperbolic import cosh, sinh
  20. from sympy.functions.elementary.trigonometric import cos, sin, sinc
  21. from sympy.functions.special.hyper import hyper, meijerg
  22. # TODO series expansions
  23. # TODO see the "Note:" in Ei
  24. # Helper function
  25. def real_to_real_as_real_imag(self, deep=True, **hints):
  26. if self.args[0].is_extended_real:
  27. if deep:
  28. hints['complex'] = False
  29. return (self.expand(deep, **hints), S.Zero)
  30. else:
  31. return (self, S.Zero)
  32. if deep:
  33. x, y = self.args[0].expand(deep, **hints).as_real_imag()
  34. else:
  35. x, y = self.args[0].as_real_imag()
  36. re = (self.func(x + I*y) + self.func(x - I*y))/2
  37. im = (self.func(x + I*y) - self.func(x - I*y))/(2*I)
  38. return (re, im)
  39. ###############################################################################
  40. ################################ ERROR FUNCTION ###############################
  41. ###############################################################################
  42. class erf(DefinedFunction):
  43. r"""
  44. The Gauss error function.
  45. Explanation
  46. ===========
  47. This function is defined as:
  48. .. math ::
  49. \mathrm{erf}(x) = \frac{2}{\sqrt{\pi}} \int_0^x e^{-t^2} \mathrm{d}t.
  50. Examples
  51. ========
  52. >>> from sympy import I, oo, erf
  53. >>> from sympy.abc import z
  54. Several special values are known:
  55. >>> erf(0)
  56. 0
  57. >>> erf(oo)
  58. 1
  59. >>> erf(-oo)
  60. -1
  61. >>> erf(I*oo)
  62. oo*I
  63. >>> erf(-I*oo)
  64. -oo*I
  65. In general one can pull out factors of -1 and $I$ from the argument:
  66. >>> erf(-z)
  67. -erf(z)
  68. The error function obeys the mirror symmetry:
  69. >>> from sympy import conjugate
  70. >>> conjugate(erf(z))
  71. erf(conjugate(z))
  72. Differentiation with respect to $z$ is supported:
  73. >>> from sympy import diff
  74. >>> diff(erf(z), z)
  75. 2*exp(-z**2)/sqrt(pi)
  76. We can numerically evaluate the error function to arbitrary precision
  77. on the whole complex plane:
  78. >>> erf(4).evalf(30)
  79. 0.999999984582742099719981147840
  80. >>> erf(-4*I).evalf(30)
  81. -1296959.73071763923152794095062*I
  82. See Also
  83. ========
  84. erfc: Complementary error function.
  85. erfi: Imaginary error function.
  86. erf2: Two-argument error function.
  87. erfinv: Inverse error function.
  88. erfcinv: Inverse Complementary error function.
  89. erf2inv: Inverse two-argument error function.
  90. References
  91. ==========
  92. .. [1] https://en.wikipedia.org/wiki/Error_function
  93. .. [2] https://dlmf.nist.gov/7
  94. .. [3] https://mathworld.wolfram.com/Erf.html
  95. .. [4] https://functions.wolfram.com/GammaBetaErf/Erf
  96. """
  97. unbranched = True
  98. def fdiff(self, argindex=1):
  99. if argindex == 1:
  100. return 2*exp(-self.args[0]**2)/sqrt(pi)
  101. else:
  102. raise ArgumentIndexError(self, argindex)
  103. def inverse(self, argindex=1):
  104. """
  105. Returns the inverse of this function.
  106. """
  107. return erfinv
  108. @classmethod
  109. def eval(cls, arg):
  110. if arg.is_Number:
  111. if arg is S.NaN:
  112. return S.NaN
  113. elif arg is S.Infinity:
  114. return S.One
  115. elif arg is S.NegativeInfinity:
  116. return S.NegativeOne
  117. elif arg.is_zero:
  118. return S.Zero
  119. if isinstance(arg, erfinv):
  120. return arg.args[0]
  121. if isinstance(arg, erfcinv):
  122. return S.One - arg.args[0]
  123. if arg.is_zero:
  124. return S.Zero
  125. # Only happens with unevaluated erf2inv
  126. if isinstance(arg, erf2inv) and arg.args[0].is_zero:
  127. return arg.args[1]
  128. # Try to pull out factors of I
  129. t = arg.extract_multiplicatively(I)
  130. if t in (S.Infinity, S.NegativeInfinity):
  131. return arg
  132. # Try to pull out factors of -1
  133. if arg.could_extract_minus_sign():
  134. return -cls(-arg)
  135. @staticmethod
  136. @cacheit
  137. def taylor_term(n, x, *previous_terms):
  138. if n < 0 or n % 2 == 0:
  139. return S.Zero
  140. else:
  141. x = sympify(x)
  142. k = floor((n - 1)/S(2))
  143. if len(previous_terms) > 2:
  144. return -previous_terms[-2] * x**2 * (n - 2)/(n*k)
  145. else:
  146. return 2*S.NegativeOne**k * x**n/(n*factorial(k)*sqrt(pi))
  147. def _eval_conjugate(self):
  148. return self.func(self.args[0].conjugate())
  149. def _eval_is_real(self):
  150. if self.args[0].is_extended_real is True:
  151. return True
  152. # There are cases where erf(z) becomes a real number
  153. # even if z is a complex number
  154. def _eval_is_imaginary(self):
  155. if self.args[0].is_imaginary is True:
  156. return True
  157. def _eval_is_finite(self):
  158. z = self.args[0]
  159. return fuzzy_or([z.is_finite, z.is_extended_real])
  160. def _eval_is_zero(self):
  161. if self.args[0].is_extended_real is True:
  162. return self.args[0].is_zero
  163. def _eval_is_positive(self):
  164. if self.args[0].is_extended_real is True:
  165. return self.args[0].is_extended_positive
  166. def _eval_is_negative(self):
  167. if self.args[0].is_extended_real is True:
  168. return self.args[0].is_extended_negative
  169. def _eval_rewrite_as_uppergamma(self, z, **kwargs):
  170. from sympy.functions.special.gamma_functions import uppergamma
  171. return sqrt(z**2)/z*(S.One - uppergamma(S.Half, z**2)/sqrt(pi))
  172. def _eval_rewrite_as_fresnels(self, z, **kwargs):
  173. arg = (S.One - I)*z/sqrt(pi)
  174. return (S.One + I)*(fresnelc(arg) - I*fresnels(arg))
  175. def _eval_rewrite_as_fresnelc(self, z, **kwargs):
  176. arg = (S.One - I)*z/sqrt(pi)
  177. return (S.One + I)*(fresnelc(arg) - I*fresnels(arg))
  178. def _eval_rewrite_as_meijerg(self, z, **kwargs):
  179. return z/sqrt(pi)*meijerg([S.Half], [], [0], [Rational(-1, 2)], z**2)
  180. def _eval_rewrite_as_hyper(self, z, **kwargs):
  181. return 2*z/sqrt(pi)*hyper([S.Half], [3*S.Half], -z**2)
  182. def _eval_rewrite_as_expint(self, z, **kwargs):
  183. return sqrt(z**2)/z - z*expint(S.Half, z**2)/sqrt(pi)
  184. def _eval_rewrite_as_tractable(self, z, limitvar=None, **kwargs):
  185. from sympy.series.limits import limit
  186. if limitvar:
  187. lim = limit(z, limitvar, S.Infinity)
  188. if lim is S.NegativeInfinity:
  189. return S.NegativeOne + _erfs(-z)*exp(-z**2)
  190. return S.One - _erfs(z)*exp(-z**2)
  191. def _eval_rewrite_as_erfc(self, z, **kwargs):
  192. return S.One - erfc(z)
  193. def _eval_rewrite_as_erfi(self, z, **kwargs):
  194. return -I*erfi(I*z)
  195. def _eval_as_leading_term(self, x, logx, cdir):
  196. arg = self.args[0].as_leading_term(x, logx=logx, cdir=cdir)
  197. arg0 = arg.subs(x, 0)
  198. if arg0 is S.ComplexInfinity:
  199. arg0 = arg.limit(x, 0, dir='-' if cdir == -1 else '+')
  200. if x in arg.free_symbols and arg0.is_zero:
  201. return 2*arg/sqrt(pi)
  202. else:
  203. return self.func(arg0)
  204. def _eval_aseries(self, n, args0, x, logx):
  205. from sympy.series.order import Order
  206. point = args0[0]
  207. if point in [S.Infinity, S.NegativeInfinity]:
  208. z = self.args[0]
  209. try:
  210. _, ex = z.leadterm(x)
  211. except (ValueError, NotImplementedError):
  212. return self
  213. ex = -ex # as x->1/x for aseries
  214. if ex.is_positive:
  215. newn = ceiling(n/ex)
  216. s = [S.NegativeOne**k * factorial2(2*k - 1) / (z**(2*k + 1) * 2**k)
  217. for k in range(newn)] + [Order(1/z**newn, x)]
  218. return S.One - (exp(-z**2)/sqrt(pi)) * Add(*s)
  219. return super(erf, self)._eval_aseries(n, args0, x, logx)
  220. as_real_imag = real_to_real_as_real_imag
  221. class erfc(DefinedFunction):
  222. r"""
  223. Complementary Error Function.
  224. Explanation
  225. ===========
  226. The function is defined as:
  227. .. math ::
  228. \mathrm{erfc}(x) = \frac{2}{\sqrt{\pi}} \int_x^\infty e^{-t^2} \mathrm{d}t
  229. Examples
  230. ========
  231. >>> from sympy import I, oo, erfc
  232. >>> from sympy.abc import z
  233. Several special values are known:
  234. >>> erfc(0)
  235. 1
  236. >>> erfc(oo)
  237. 0
  238. >>> erfc(-oo)
  239. 2
  240. >>> erfc(I*oo)
  241. -oo*I
  242. >>> erfc(-I*oo)
  243. oo*I
  244. The error function obeys the mirror symmetry:
  245. >>> from sympy import conjugate
  246. >>> conjugate(erfc(z))
  247. erfc(conjugate(z))
  248. Differentiation with respect to $z$ is supported:
  249. >>> from sympy import diff
  250. >>> diff(erfc(z), z)
  251. -2*exp(-z**2)/sqrt(pi)
  252. It also follows
  253. >>> erfc(-z)
  254. 2 - erfc(z)
  255. We can numerically evaluate the complementary error function to arbitrary
  256. precision on the whole complex plane:
  257. >>> erfc(4).evalf(30)
  258. 0.0000000154172579002800188521596734869
  259. >>> erfc(4*I).evalf(30)
  260. 1.0 - 1296959.73071763923152794095062*I
  261. See Also
  262. ========
  263. erf: Gaussian error function.
  264. erfi: Imaginary error function.
  265. erf2: Two-argument error function.
  266. erfinv: Inverse error function.
  267. erfcinv: Inverse Complementary error function.
  268. erf2inv: Inverse two-argument error function.
  269. References
  270. ==========
  271. .. [1] https://en.wikipedia.org/wiki/Error_function
  272. .. [2] https://dlmf.nist.gov/7
  273. .. [3] https://mathworld.wolfram.com/Erfc.html
  274. .. [4] https://functions.wolfram.com/GammaBetaErf/Erfc
  275. """
  276. unbranched = True
  277. def fdiff(self, argindex=1):
  278. if argindex == 1:
  279. return -2*exp(-self.args[0]**2)/sqrt(pi)
  280. else:
  281. raise ArgumentIndexError(self, argindex)
  282. def inverse(self, argindex=1):
  283. """
  284. Returns the inverse of this function.
  285. """
  286. return erfcinv
  287. @classmethod
  288. def eval(cls, arg):
  289. if arg.is_Number:
  290. if arg is S.NaN:
  291. return S.NaN
  292. elif arg is S.Infinity:
  293. return S.Zero
  294. elif arg.is_zero:
  295. return S.One
  296. if isinstance(arg, erfinv):
  297. return S.One - arg.args[0]
  298. if isinstance(arg, erfcinv):
  299. return arg.args[0]
  300. if arg.is_zero:
  301. return S.One
  302. # Try to pull out factors of I
  303. t = arg.extract_multiplicatively(I)
  304. if t in (S.Infinity, S.NegativeInfinity):
  305. return -arg
  306. # Try to pull out factors of -1
  307. if arg.could_extract_minus_sign():
  308. return 2 - cls(-arg)
  309. @staticmethod
  310. @cacheit
  311. def taylor_term(n, x, *previous_terms):
  312. if n == 0:
  313. return S.One
  314. elif n < 0 or n % 2 == 0:
  315. return S.Zero
  316. else:
  317. x = sympify(x)
  318. k = floor((n - 1)/S(2))
  319. if len(previous_terms) > 2:
  320. return -previous_terms[-2] * x**2 * (n - 2)/(n*k)
  321. else:
  322. return -2*S.NegativeOne**k * x**n/(n*factorial(k)*sqrt(pi))
  323. def _eval_conjugate(self):
  324. return self.func(self.args[0].conjugate())
  325. def _eval_is_real(self):
  326. if self.args[0].is_extended_real is True:
  327. return True
  328. if self.args[0].is_imaginary is True:
  329. return False
  330. def _eval_rewrite_as_tractable(self, z, limitvar=None, **kwargs):
  331. return self.rewrite(erf).rewrite("tractable", deep=True, limitvar=limitvar)
  332. def _eval_rewrite_as_erf(self, z, **kwargs):
  333. return S.One - erf(z)
  334. def _eval_rewrite_as_erfi(self, z, **kwargs):
  335. return S.One + I*erfi(I*z)
  336. def _eval_rewrite_as_fresnels(self, z, **kwargs):
  337. arg = (S.One - I)*z/sqrt(pi)
  338. return S.One - (S.One + I)*(fresnelc(arg) - I*fresnels(arg))
  339. def _eval_rewrite_as_fresnelc(self, z, **kwargs):
  340. arg = (S.One-I)*z/sqrt(pi)
  341. return S.One - (S.One + I)*(fresnelc(arg) - I*fresnels(arg))
  342. def _eval_rewrite_as_meijerg(self, z, **kwargs):
  343. return S.One - z/sqrt(pi)*meijerg([S.Half], [], [0], [Rational(-1, 2)], z**2)
  344. def _eval_rewrite_as_hyper(self, z, **kwargs):
  345. return S.One - 2*z/sqrt(pi)*hyper([S.Half], [3*S.Half], -z**2)
  346. def _eval_rewrite_as_uppergamma(self, z, **kwargs):
  347. from sympy.functions.special.gamma_functions import uppergamma
  348. return S.One - sqrt(z**2)/z*(S.One - uppergamma(S.Half, z**2)/sqrt(pi))
  349. def _eval_rewrite_as_expint(self, z, **kwargs):
  350. return S.One - sqrt(z**2)/z + z*expint(S.Half, z**2)/sqrt(pi)
  351. def _eval_expand_func(self, **hints):
  352. return self.rewrite(erf)
  353. def _eval_as_leading_term(self, x, logx, cdir):
  354. arg = self.args[0].as_leading_term(x, logx=logx, cdir=cdir)
  355. arg0 = arg.subs(x, 0)
  356. if arg0 is S.ComplexInfinity:
  357. arg0 = arg.limit(x, 0, dir='-' if cdir == -1 else '+')
  358. if arg0.is_zero:
  359. return S.One
  360. else:
  361. return self.func(arg0)
  362. as_real_imag = real_to_real_as_real_imag
  363. def _eval_aseries(self, n, args0, x, logx):
  364. return S.One - erf(*self.args)._eval_aseries(n, args0, x, logx)
  365. class erfi(DefinedFunction):
  366. r"""
  367. Imaginary error function.
  368. Explanation
  369. ===========
  370. The function erfi is defined as:
  371. .. math ::
  372. \mathrm{erfi}(x) = \frac{2}{\sqrt{\pi}} \int_0^x e^{t^2} \mathrm{d}t
  373. Examples
  374. ========
  375. >>> from sympy import I, oo, erfi
  376. >>> from sympy.abc import z
  377. Several special values are known:
  378. >>> erfi(0)
  379. 0
  380. >>> erfi(oo)
  381. oo
  382. >>> erfi(-oo)
  383. -oo
  384. >>> erfi(I*oo)
  385. I
  386. >>> erfi(-I*oo)
  387. -I
  388. In general one can pull out factors of -1 and $I$ from the argument:
  389. >>> erfi(-z)
  390. -erfi(z)
  391. >>> from sympy import conjugate
  392. >>> conjugate(erfi(z))
  393. erfi(conjugate(z))
  394. Differentiation with respect to $z$ is supported:
  395. >>> from sympy import diff
  396. >>> diff(erfi(z), z)
  397. 2*exp(z**2)/sqrt(pi)
  398. We can numerically evaluate the imaginary error function to arbitrary
  399. precision on the whole complex plane:
  400. >>> erfi(2).evalf(30)
  401. 18.5648024145755525987042919132
  402. >>> erfi(-2*I).evalf(30)
  403. -0.995322265018952734162069256367*I
  404. See Also
  405. ========
  406. erf: Gaussian error function.
  407. erfc: Complementary error function.
  408. erf2: Two-argument error function.
  409. erfinv: Inverse error function.
  410. erfcinv: Inverse Complementary error function.
  411. erf2inv: Inverse two-argument error function.
  412. References
  413. ==========
  414. .. [1] https://en.wikipedia.org/wiki/Error_function
  415. .. [2] https://mathworld.wolfram.com/Erfi.html
  416. .. [3] https://functions.wolfram.com/GammaBetaErf/Erfi
  417. """
  418. unbranched = True
  419. def fdiff(self, argindex=1):
  420. if argindex == 1:
  421. return 2*exp(self.args[0]**2)/sqrt(pi)
  422. else:
  423. raise ArgumentIndexError(self, argindex)
  424. @classmethod
  425. def eval(cls, z):
  426. if z.is_Number:
  427. if z is S.NaN:
  428. return S.NaN
  429. elif z.is_zero:
  430. return S.Zero
  431. elif z is S.Infinity:
  432. return S.Infinity
  433. if z.is_zero:
  434. return S.Zero
  435. # Try to pull out factors of -1
  436. if z.could_extract_minus_sign():
  437. return -cls(-z)
  438. # Try to pull out factors of I
  439. nz = z.extract_multiplicatively(I)
  440. if nz is not None:
  441. if nz is S.Infinity:
  442. return I
  443. if isinstance(nz, erfinv):
  444. return I*nz.args[0]
  445. if isinstance(nz, erfcinv):
  446. return I*(S.One - nz.args[0])
  447. # Only happens with unevaluated erf2inv
  448. if isinstance(nz, erf2inv) and nz.args[0].is_zero:
  449. return I*nz.args[1]
  450. @staticmethod
  451. @cacheit
  452. def taylor_term(n, x, *previous_terms):
  453. if n < 0 or n % 2 == 0:
  454. return S.Zero
  455. else:
  456. x = sympify(x)
  457. k = floor((n - 1)/S(2))
  458. if len(previous_terms) > 2:
  459. return previous_terms[-2] * x**2 * (n - 2)/(n*k)
  460. else:
  461. return 2 * x**n/(n*factorial(k)*sqrt(pi))
  462. def _eval_conjugate(self):
  463. return self.func(self.args[0].conjugate())
  464. def _eval_is_extended_real(self):
  465. return self.args[0].is_extended_real
  466. def _eval_is_zero(self):
  467. return self.args[0].is_zero
  468. def _eval_rewrite_as_tractable(self, z, limitvar=None, **kwargs):
  469. return self.rewrite(erf).rewrite("tractable", deep=True, limitvar=limitvar)
  470. def _eval_rewrite_as_erf(self, z, **kwargs):
  471. return -I*erf(I*z)
  472. def _eval_rewrite_as_erfc(self, z, **kwargs):
  473. return I*erfc(I*z) - I
  474. def _eval_rewrite_as_fresnels(self, z, **kwargs):
  475. arg = (S.One + I)*z/sqrt(pi)
  476. return (S.One - I)*(fresnelc(arg) - I*fresnels(arg))
  477. def _eval_rewrite_as_fresnelc(self, z, **kwargs):
  478. arg = (S.One + I)*z/sqrt(pi)
  479. return (S.One - I)*(fresnelc(arg) - I*fresnels(arg))
  480. def _eval_rewrite_as_meijerg(self, z, **kwargs):
  481. return z/sqrt(pi)*meijerg([S.Half], [], [0], [Rational(-1, 2)], -z**2)
  482. def _eval_rewrite_as_hyper(self, z, **kwargs):
  483. return 2*z/sqrt(pi)*hyper([S.Half], [3*S.Half], z**2)
  484. def _eval_rewrite_as_uppergamma(self, z, **kwargs):
  485. from sympy.functions.special.gamma_functions import uppergamma
  486. return sqrt(-z**2)/z*(uppergamma(S.Half, -z**2)/sqrt(pi) - S.One)
  487. def _eval_rewrite_as_expint(self, z, **kwargs):
  488. return sqrt(-z**2)/z - z*expint(S.Half, -z**2)/sqrt(pi)
  489. def _eval_expand_func(self, **hints):
  490. return self.rewrite(erf)
  491. as_real_imag = real_to_real_as_real_imag
  492. def _eval_as_leading_term(self, x, logx, cdir):
  493. arg = self.args[0].as_leading_term(x, logx=logx, cdir=cdir)
  494. arg0 = arg.subs(x, 0)
  495. if x in arg.free_symbols and arg0.is_zero:
  496. return 2*arg/sqrt(pi)
  497. elif arg0.is_finite:
  498. return self.func(arg0)
  499. return self.func(arg)
  500. def _eval_aseries(self, n, args0, x, logx):
  501. from sympy.series.order import Order
  502. point = args0[0]
  503. if point is S.Infinity:
  504. z = self.args[0]
  505. s = [factorial2(2*k - 1) / (2**k * z**(2*k + 1))
  506. for k in range(n)] + [Order(1/z**n, x)]
  507. return -I + (exp(z**2)/sqrt(pi)) * Add(*s)
  508. return super(erfi, self)._eval_aseries(n, args0, x, logx)
  509. class erf2(DefinedFunction):
  510. r"""
  511. Two-argument error function.
  512. Explanation
  513. ===========
  514. This function is defined as:
  515. .. math ::
  516. \mathrm{erf2}(x, y) = \frac{2}{\sqrt{\pi}} \int_x^y e^{-t^2} \mathrm{d}t
  517. Examples
  518. ========
  519. >>> from sympy import oo, erf2
  520. >>> from sympy.abc import x, y
  521. Several special values are known:
  522. >>> erf2(0, 0)
  523. 0
  524. >>> erf2(x, x)
  525. 0
  526. >>> erf2(x, oo)
  527. 1 - erf(x)
  528. >>> erf2(x, -oo)
  529. -erf(x) - 1
  530. >>> erf2(oo, y)
  531. erf(y) - 1
  532. >>> erf2(-oo, y)
  533. erf(y) + 1
  534. In general one can pull out factors of -1:
  535. >>> erf2(-x, -y)
  536. -erf2(x, y)
  537. The error function obeys the mirror symmetry:
  538. >>> from sympy import conjugate
  539. >>> conjugate(erf2(x, y))
  540. erf2(conjugate(x), conjugate(y))
  541. Differentiation with respect to $x$, $y$ is supported:
  542. >>> from sympy import diff
  543. >>> diff(erf2(x, y), x)
  544. -2*exp(-x**2)/sqrt(pi)
  545. >>> diff(erf2(x, y), y)
  546. 2*exp(-y**2)/sqrt(pi)
  547. See Also
  548. ========
  549. erf: Gaussian error function.
  550. erfc: Complementary error function.
  551. erfi: Imaginary error function.
  552. erfinv: Inverse error function.
  553. erfcinv: Inverse Complementary error function.
  554. erf2inv: Inverse two-argument error function.
  555. References
  556. ==========
  557. .. [1] https://functions.wolfram.com/GammaBetaErf/Erf2/
  558. """
  559. def fdiff(self, argindex):
  560. x, y = self.args
  561. if argindex == 1:
  562. return -2*exp(-x**2)/sqrt(pi)
  563. elif argindex == 2:
  564. return 2*exp(-y**2)/sqrt(pi)
  565. else:
  566. raise ArgumentIndexError(self, argindex)
  567. @classmethod
  568. def eval(cls, x, y):
  569. chk = (S.Infinity, S.NegativeInfinity, S.Zero)
  570. if x is S.NaN or y is S.NaN:
  571. return S.NaN
  572. elif x == y:
  573. return S.Zero
  574. elif x in chk or y in chk:
  575. return erf(y) - erf(x)
  576. if isinstance(y, erf2inv) and y.args[0] == x:
  577. return y.args[1]
  578. if x.is_zero or y.is_zero or x.is_extended_real and x.is_infinite or \
  579. y.is_extended_real and y.is_infinite:
  580. return erf(y) - erf(x)
  581. #Try to pull out -1 factor
  582. sign_x = x.could_extract_minus_sign()
  583. sign_y = y.could_extract_minus_sign()
  584. if (sign_x and sign_y):
  585. return -cls(-x, -y)
  586. elif (sign_x or sign_y):
  587. return erf(y)-erf(x)
  588. def _eval_conjugate(self):
  589. return self.func(self.args[0].conjugate(), self.args[1].conjugate())
  590. def _eval_is_extended_real(self):
  591. return self.args[0].is_extended_real and self.args[1].is_extended_real
  592. def _eval_rewrite_as_erf(self, x, y, **kwargs):
  593. return erf(y) - erf(x)
  594. def _eval_rewrite_as_erfc(self, x, y, **kwargs):
  595. return erfc(x) - erfc(y)
  596. def _eval_rewrite_as_erfi(self, x, y, **kwargs):
  597. return I*(erfi(I*x)-erfi(I*y))
  598. def _eval_rewrite_as_fresnels(self, x, y, **kwargs):
  599. return erf(y).rewrite(fresnels) - erf(x).rewrite(fresnels)
  600. def _eval_rewrite_as_fresnelc(self, x, y, **kwargs):
  601. return erf(y).rewrite(fresnelc) - erf(x).rewrite(fresnelc)
  602. def _eval_rewrite_as_meijerg(self, x, y, **kwargs):
  603. return erf(y).rewrite(meijerg) - erf(x).rewrite(meijerg)
  604. def _eval_rewrite_as_hyper(self, x, y, **kwargs):
  605. return erf(y).rewrite(hyper) - erf(x).rewrite(hyper)
  606. def _eval_rewrite_as_uppergamma(self, x, y, **kwargs):
  607. from sympy.functions.special.gamma_functions import uppergamma
  608. return (sqrt(y**2)/y*(S.One - uppergamma(S.Half, y**2)/sqrt(pi)) -
  609. sqrt(x**2)/x*(S.One - uppergamma(S.Half, x**2)/sqrt(pi)))
  610. def _eval_rewrite_as_expint(self, x, y, **kwargs):
  611. return erf(y).rewrite(expint) - erf(x).rewrite(expint)
  612. def _eval_expand_func(self, **hints):
  613. return self.rewrite(erf)
  614. def _eval_is_zero(self):
  615. return is_eq(*self.args)
  616. class erfinv(DefinedFunction):
  617. r"""
  618. Inverse Error Function. The erfinv function is defined as:
  619. .. math ::
  620. \mathrm{erf}(x) = y \quad \Rightarrow \quad \mathrm{erfinv}(y) = x
  621. Examples
  622. ========
  623. >>> from sympy import erfinv
  624. >>> from sympy.abc import x
  625. Several special values are known:
  626. >>> erfinv(0)
  627. 0
  628. >>> erfinv(1)
  629. oo
  630. Differentiation with respect to $x$ is supported:
  631. >>> from sympy import diff
  632. >>> diff(erfinv(x), x)
  633. sqrt(pi)*exp(erfinv(x)**2)/2
  634. We can numerically evaluate the inverse error function to arbitrary
  635. precision on [-1, 1]:
  636. >>> erfinv(0.2).evalf(30)
  637. 0.179143454621291692285822705344
  638. See Also
  639. ========
  640. erf: Gaussian error function.
  641. erfc: Complementary error function.
  642. erfi: Imaginary error function.
  643. erf2: Two-argument error function.
  644. erfcinv: Inverse Complementary error function.
  645. erf2inv: Inverse two-argument error function.
  646. References
  647. ==========
  648. .. [1] https://en.wikipedia.org/wiki/Error_function#Inverse_functions
  649. .. [2] https://functions.wolfram.com/GammaBetaErf/InverseErf/
  650. """
  651. def fdiff(self, argindex =1):
  652. if argindex == 1:
  653. return sqrt(pi)*exp(self.func(self.args[0])**2)*S.Half
  654. else :
  655. raise ArgumentIndexError(self, argindex)
  656. def inverse(self, argindex=1):
  657. """
  658. Returns the inverse of this function.
  659. """
  660. return erf
  661. @classmethod
  662. def eval(cls, z):
  663. if z is S.NaN:
  664. return S.NaN
  665. elif z is S.NegativeOne:
  666. return S.NegativeInfinity
  667. elif z.is_zero:
  668. return S.Zero
  669. elif z is S.One:
  670. return S.Infinity
  671. if isinstance(z, erf) and z.args[0].is_extended_real:
  672. return z.args[0]
  673. if z.is_zero:
  674. return S.Zero
  675. # Try to pull out factors of -1
  676. nz = z.extract_multiplicatively(-1)
  677. if nz is not None and (isinstance(nz, erf) and (nz.args[0]).is_extended_real):
  678. return -nz.args[0]
  679. def _eval_rewrite_as_erfcinv(self, z, **kwargs):
  680. return erfcinv(1-z)
  681. def _eval_is_zero(self):
  682. return self.args[0].is_zero
  683. class erfcinv (DefinedFunction):
  684. r"""
  685. Inverse Complementary Error Function. The erfcinv function is defined as:
  686. .. math ::
  687. \mathrm{erfc}(x) = y \quad \Rightarrow \quad \mathrm{erfcinv}(y) = x
  688. Examples
  689. ========
  690. >>> from sympy import erfcinv
  691. >>> from sympy.abc import x
  692. Several special values are known:
  693. >>> erfcinv(1)
  694. 0
  695. >>> erfcinv(0)
  696. oo
  697. Differentiation with respect to $x$ is supported:
  698. >>> from sympy import diff
  699. >>> diff(erfcinv(x), x)
  700. -sqrt(pi)*exp(erfcinv(x)**2)/2
  701. See Also
  702. ========
  703. erf: Gaussian error function.
  704. erfc: Complementary error function.
  705. erfi: Imaginary error function.
  706. erf2: Two-argument error function.
  707. erfinv: Inverse error function.
  708. erf2inv: Inverse two-argument error function.
  709. References
  710. ==========
  711. .. [1] https://en.wikipedia.org/wiki/Error_function#Inverse_functions
  712. .. [2] https://functions.wolfram.com/GammaBetaErf/InverseErfc/
  713. """
  714. def fdiff(self, argindex =1):
  715. if argindex == 1:
  716. return -sqrt(pi)*exp(self.func(self.args[0])**2)*S.Half
  717. else:
  718. raise ArgumentIndexError(self, argindex)
  719. def inverse(self, argindex=1):
  720. """
  721. Returns the inverse of this function.
  722. """
  723. return erfc
  724. @classmethod
  725. def eval(cls, z):
  726. if z is S.NaN:
  727. return S.NaN
  728. elif z.is_zero:
  729. return S.Infinity
  730. elif z is S.One:
  731. return S.Zero
  732. elif z == 2:
  733. return S.NegativeInfinity
  734. if z.is_zero:
  735. return S.Infinity
  736. def _eval_rewrite_as_erfinv(self, z, **kwargs):
  737. return erfinv(1-z)
  738. def _eval_is_zero(self):
  739. return (self.args[0] - 1).is_zero
  740. def _eval_is_infinite(self):
  741. z = self.args[0]
  742. return fuzzy_or([z.is_zero, is_eq(z, Integer(2))])
  743. class erf2inv(DefinedFunction):
  744. r"""
  745. Two-argument Inverse error function. The erf2inv function is defined as:
  746. .. math ::
  747. \mathrm{erf2}(x, w) = y \quad \Rightarrow \quad \mathrm{erf2inv}(x, y) = w
  748. Examples
  749. ========
  750. >>> from sympy import erf2inv, oo
  751. >>> from sympy.abc import x, y
  752. Several special values are known:
  753. >>> erf2inv(0, 0)
  754. 0
  755. >>> erf2inv(1, 0)
  756. 1
  757. >>> erf2inv(0, 1)
  758. oo
  759. >>> erf2inv(0, y)
  760. erfinv(y)
  761. >>> erf2inv(oo, y)
  762. erfcinv(-y)
  763. Differentiation with respect to $x$ and $y$ is supported:
  764. >>> from sympy import diff
  765. >>> diff(erf2inv(x, y), x)
  766. exp(-x**2 + erf2inv(x, y)**2)
  767. >>> diff(erf2inv(x, y), y)
  768. sqrt(pi)*exp(erf2inv(x, y)**2)/2
  769. See Also
  770. ========
  771. erf: Gaussian error function.
  772. erfc: Complementary error function.
  773. erfi: Imaginary error function.
  774. erf2: Two-argument error function.
  775. erfinv: Inverse error function.
  776. erfcinv: Inverse complementary error function.
  777. References
  778. ==========
  779. .. [1] https://functions.wolfram.com/GammaBetaErf/InverseErf2/
  780. """
  781. def fdiff(self, argindex):
  782. x, y = self.args
  783. if argindex == 1:
  784. return exp(self.func(x,y)**2-x**2)
  785. elif argindex == 2:
  786. return sqrt(pi)*S.Half*exp(self.func(x,y)**2)
  787. else:
  788. raise ArgumentIndexError(self, argindex)
  789. @classmethod
  790. def eval(cls, x, y):
  791. if x is S.NaN or y is S.NaN:
  792. return S.NaN
  793. elif x.is_zero and y.is_zero:
  794. return S.Zero
  795. elif x.is_zero and y is S.One:
  796. return S.Infinity
  797. elif x is S.One and y.is_zero:
  798. return S.One
  799. elif x.is_zero:
  800. return erfinv(y)
  801. elif x is S.Infinity:
  802. return erfcinv(-y)
  803. elif y.is_zero:
  804. return x
  805. elif y is S.Infinity:
  806. return erfinv(x)
  807. if x.is_zero:
  808. if y.is_zero:
  809. return S.Zero
  810. else:
  811. return erfinv(y)
  812. if y.is_zero:
  813. return x
  814. def _eval_is_zero(self):
  815. x, y = self.args
  816. if x.is_zero and y.is_zero:
  817. return True
  818. ###############################################################################
  819. #################### EXPONENTIAL INTEGRALS ####################################
  820. ###############################################################################
  821. class Ei(DefinedFunction):
  822. r"""
  823. The classical exponential integral.
  824. Explanation
  825. ===========
  826. For use in SymPy, this function is defined as
  827. .. math:: \operatorname{Ei}(x) = \sum_{n=1}^\infty \frac{x^n}{n\, n!}
  828. + \log(x) + \gamma,
  829. where $\gamma$ is the Euler-Mascheroni constant.
  830. If $x$ is a polar number, this defines an analytic function on the
  831. Riemann surface of the logarithm. Otherwise this defines an analytic
  832. function in the cut plane $\mathbb{C} \setminus (-\infty, 0]$.
  833. **Background**
  834. The name exponential integral comes from the following statement:
  835. .. math:: \operatorname{Ei}(x) = \int_{-\infty}^x \frac{e^t}{t} \mathrm{d}t
  836. If the integral is interpreted as a Cauchy principal value, this statement
  837. holds for $x > 0$ and $\operatorname{Ei}(x)$ as defined above.
  838. Examples
  839. ========
  840. >>> from sympy import Ei, polar_lift, exp_polar, I, pi
  841. >>> from sympy.abc import x
  842. >>> Ei(-1)
  843. Ei(-1)
  844. This yields a real value:
  845. >>> Ei(-1).n(chop=True)
  846. -0.219383934395520
  847. On the other hand the analytic continuation is not real:
  848. >>> Ei(polar_lift(-1)).n(chop=True)
  849. -0.21938393439552 + 3.14159265358979*I
  850. The exponential integral has a logarithmic branch point at the origin:
  851. >>> Ei(x*exp_polar(2*I*pi))
  852. Ei(x) + 2*I*pi
  853. Differentiation is supported:
  854. >>> Ei(x).diff(x)
  855. exp(x)/x
  856. The exponential integral is related to many other special functions.
  857. For example:
  858. >>> from sympy import expint, Shi
  859. >>> Ei(x).rewrite(expint)
  860. -expint(1, x*exp_polar(I*pi)) - I*pi
  861. >>> Ei(x).rewrite(Shi)
  862. Chi(x) + Shi(x)
  863. See Also
  864. ========
  865. expint: Generalised exponential integral.
  866. E1: Special case of the generalised exponential integral.
  867. li: Logarithmic integral.
  868. Li: Offset logarithmic integral.
  869. Si: Sine integral.
  870. Ci: Cosine integral.
  871. Shi: Hyperbolic sine integral.
  872. Chi: Hyperbolic cosine integral.
  873. uppergamma: Upper incomplete gamma function.
  874. References
  875. ==========
  876. .. [1] https://dlmf.nist.gov/6.6
  877. .. [2] https://en.wikipedia.org/wiki/Exponential_integral
  878. .. [3] Abramowitz & Stegun, section 5: https://web.archive.org/web/20201128173312/http://people.math.sfu.ca/~cbm/aands/page_228.htm
  879. """
  880. @classmethod
  881. def eval(cls, z):
  882. if z.is_zero:
  883. return S.NegativeInfinity
  884. elif z is S.Infinity:
  885. return S.Infinity
  886. elif z is S.NegativeInfinity:
  887. return S.Zero
  888. if z.is_zero:
  889. return S.NegativeInfinity
  890. nz, n = z.extract_branch_factor()
  891. if n:
  892. return Ei(nz) + 2*I*pi*n
  893. def fdiff(self, argindex=1):
  894. arg = unpolarify(self.args[0])
  895. if argindex == 1:
  896. return exp(arg)/arg
  897. else:
  898. raise ArgumentIndexError(self, argindex)
  899. def _eval_evalf(self, prec):
  900. if (self.args[0]/polar_lift(-1)).is_positive:
  901. return super()._eval_evalf(prec) + (I*pi)._eval_evalf(prec)
  902. return super()._eval_evalf(prec)
  903. def _eval_rewrite_as_uppergamma(self, z, **kwargs):
  904. from sympy.functions.special.gamma_functions import uppergamma
  905. # XXX this does not currently work usefully because uppergamma
  906. # immediately turns into expint
  907. return -uppergamma(0, polar_lift(-1)*z) - I*pi
  908. def _eval_rewrite_as_expint(self, z, **kwargs):
  909. return -expint(1, polar_lift(-1)*z) - I*pi
  910. def _eval_rewrite_as_li(self, z, **kwargs):
  911. if isinstance(z, log):
  912. return li(z.args[0])
  913. # TODO:
  914. # Actually it only holds that:
  915. # Ei(z) = li(exp(z))
  916. # for -pi < imag(z) <= pi
  917. return li(exp(z))
  918. def _eval_rewrite_as_Si(self, z, **kwargs):
  919. if z.is_negative:
  920. return Shi(z) + Chi(z) - I*pi
  921. else:
  922. return Shi(z) + Chi(z)
  923. _eval_rewrite_as_Ci = _eval_rewrite_as_Si
  924. _eval_rewrite_as_Chi = _eval_rewrite_as_Si
  925. _eval_rewrite_as_Shi = _eval_rewrite_as_Si
  926. def _eval_rewrite_as_tractable(self, z, limitvar=None, **kwargs):
  927. return exp(z) * _eis(z)
  928. def _eval_rewrite_as_Integral(self, z, **kwargs):
  929. from sympy.integrals.integrals import Integral
  930. t = Dummy(uniquely_named_symbol('t', [z]).name)
  931. return Integral(S.Exp1**t/t, (t, S.NegativeInfinity, z))
  932. def _eval_as_leading_term(self, x, logx, cdir):
  933. from sympy import re
  934. x0 = self.args[0].limit(x, 0)
  935. arg = self.args[0].as_leading_term(x, cdir=cdir)
  936. cdir = arg.dir(x, cdir)
  937. if x0.is_zero:
  938. c, e = arg.as_coeff_exponent(x)
  939. logx = log(x) if logx is None else logx
  940. return log(c) + e*logx + EulerGamma - (
  941. I*pi if re(cdir).is_negative else S.Zero)
  942. return super()._eval_as_leading_term(x, logx=logx, cdir=cdir)
  943. def _eval_nseries(self, x, n, logx, cdir=0):
  944. x0 = self.args[0].limit(x, 0)
  945. if x0.is_zero:
  946. f = self._eval_rewrite_as_Si(*self.args)
  947. return f._eval_nseries(x, n, logx)
  948. return super()._eval_nseries(x, n, logx)
  949. def _eval_aseries(self, n, args0, x, logx):
  950. from sympy.series.order import Order
  951. point = args0[0]
  952. if point in (S.Infinity, S.NegativeInfinity):
  953. z = self.args[0]
  954. s = [factorial(k) / (z)**k for k in range(n)] + \
  955. [Order(1/z**n, x)]
  956. return (exp(z)/z) * Add(*s)
  957. return super(Ei, self)._eval_aseries(n, args0, x, logx)
  958. class expint(DefinedFunction):
  959. r"""
  960. Generalized exponential integral.
  961. Explanation
  962. ===========
  963. This function is defined as
  964. .. math:: \operatorname{E}_\nu(z) = z^{\nu - 1} \Gamma(1 - \nu, z),
  965. where $\Gamma(1 - \nu, z)$ is the upper incomplete gamma function
  966. (``uppergamma``).
  967. Hence for $z$ with positive real part we have
  968. .. math:: \operatorname{E}_\nu(z)
  969. = \int_1^\infty \frac{e^{-zt}}{t^\nu} \mathrm{d}t,
  970. which explains the name.
  971. The representation as an incomplete gamma function provides an analytic
  972. continuation for $\operatorname{E}_\nu(z)$. If $\nu$ is a
  973. non-positive integer, the exponential integral is thus an unbranched
  974. function of $z$, otherwise there is a branch point at the origin.
  975. Refer to the incomplete gamma function documentation for details of the
  976. branching behavior.
  977. Examples
  978. ========
  979. >>> from sympy import expint, S
  980. >>> from sympy.abc import nu, z
  981. Differentiation is supported. Differentiation with respect to $z$ further
  982. explains the name: for integral orders, the exponential integral is an
  983. iterated integral of the exponential function.
  984. >>> expint(nu, z).diff(z)
  985. -expint(nu - 1, z)
  986. Differentiation with respect to $\nu$ has no classical expression:
  987. >>> expint(nu, z).diff(nu)
  988. -z**(nu - 1)*meijerg(((), (1, 1)), ((0, 0, 1 - nu), ()), z)
  989. At non-postive integer orders, the exponential integral reduces to the
  990. exponential function:
  991. >>> expint(0, z)
  992. exp(-z)/z
  993. >>> expint(-1, z)
  994. exp(-z)/z + exp(-z)/z**2
  995. At half-integers it reduces to error functions:
  996. >>> expint(S(1)/2, z)
  997. sqrt(pi)*erfc(sqrt(z))/sqrt(z)
  998. At positive integer orders it can be rewritten in terms of exponentials
  999. and ``expint(1, z)``. Use ``expand_func()`` to do this:
  1000. >>> from sympy import expand_func
  1001. >>> expand_func(expint(5, z))
  1002. z**4*expint(1, z)/24 + (-z**3 + z**2 - 2*z + 6)*exp(-z)/24
  1003. The generalised exponential integral is essentially equivalent to the
  1004. incomplete gamma function:
  1005. >>> from sympy import uppergamma
  1006. >>> expint(nu, z).rewrite(uppergamma)
  1007. z**(nu - 1)*uppergamma(1 - nu, z)
  1008. As such it is branched at the origin:
  1009. >>> from sympy import exp_polar, pi, I
  1010. >>> expint(4, z*exp_polar(2*pi*I))
  1011. I*pi*z**3/3 + expint(4, z)
  1012. >>> expint(nu, z*exp_polar(2*pi*I))
  1013. z**(nu - 1)*(exp(2*I*pi*nu) - 1)*gamma(1 - nu) + expint(nu, z)
  1014. See Also
  1015. ========
  1016. Ei: Another related function called exponential integral.
  1017. E1: The classical case, returns expint(1, z).
  1018. li: Logarithmic integral.
  1019. Li: Offset logarithmic integral.
  1020. Si: Sine integral.
  1021. Ci: Cosine integral.
  1022. Shi: Hyperbolic sine integral.
  1023. Chi: Hyperbolic cosine integral.
  1024. uppergamma
  1025. References
  1026. ==========
  1027. .. [1] https://dlmf.nist.gov/8.19
  1028. .. [2] https://functions.wolfram.com/GammaBetaErf/ExpIntegralE/
  1029. .. [3] https://en.wikipedia.org/wiki/Exponential_integral
  1030. """
  1031. @classmethod
  1032. def eval(cls, nu, z):
  1033. from sympy.functions.special.gamma_functions import (gamma, uppergamma)
  1034. nu2 = unpolarify(nu)
  1035. if nu != nu2:
  1036. return expint(nu2, z)
  1037. if nu.is_Integer and nu <= 0 or (not nu.is_Integer and (2*nu).is_Integer):
  1038. return unpolarify(expand_mul(z**(nu - 1)*uppergamma(1 - nu, z)))
  1039. # Extract branching information. This can be deduced from what is
  1040. # explained in lowergamma.eval().
  1041. z, n = z.extract_branch_factor()
  1042. if n is S.Zero:
  1043. return
  1044. if nu.is_integer:
  1045. if not nu > 0:
  1046. return
  1047. return expint(nu, z) \
  1048. - 2*pi*I*n*S.NegativeOne**(nu - 1)/factorial(nu - 1)*unpolarify(z)**(nu - 1)
  1049. else:
  1050. return (exp(2*I*pi*nu*n) - 1)*z**(nu - 1)*gamma(1 - nu) + expint(nu, z)
  1051. def fdiff(self, argindex):
  1052. nu, z = self.args
  1053. if argindex == 1:
  1054. return -z**(nu - 1)*meijerg([], [1, 1], [0, 0, 1 - nu], [], z)
  1055. elif argindex == 2:
  1056. return -expint(nu - 1, z)
  1057. else:
  1058. raise ArgumentIndexError(self, argindex)
  1059. def _eval_rewrite_as_uppergamma(self, nu, z, **kwargs):
  1060. from sympy.functions.special.gamma_functions import uppergamma
  1061. return z**(nu - 1)*uppergamma(1 - nu, z)
  1062. def _eval_rewrite_as_Ei(self, nu, z, **kwargs):
  1063. if nu == 1:
  1064. return -Ei(z*exp_polar(-I*pi)) - I*pi
  1065. elif nu.is_Integer and nu > 1:
  1066. # DLMF, 8.19.7
  1067. x = -unpolarify(z)
  1068. return x**(nu - 1)/factorial(nu - 1)*E1(z).rewrite(Ei) + \
  1069. exp(x)/factorial(nu - 1) * \
  1070. Add(*[factorial(nu - k - 2)*x**k for k in range(nu - 1)])
  1071. else:
  1072. return self
  1073. def _eval_expand_func(self, **hints):
  1074. return self.rewrite(Ei).rewrite(expint, **hints)
  1075. def _eval_rewrite_as_Si(self, nu, z, **kwargs):
  1076. if nu != 1:
  1077. return self
  1078. return Shi(z) - Chi(z)
  1079. _eval_rewrite_as_Ci = _eval_rewrite_as_Si
  1080. _eval_rewrite_as_Chi = _eval_rewrite_as_Si
  1081. _eval_rewrite_as_Shi = _eval_rewrite_as_Si
  1082. def _eval_nseries(self, x, n, logx, cdir=0):
  1083. if not self.args[0].has(x):
  1084. nu = self.args[0]
  1085. if nu == 1:
  1086. f = self._eval_rewrite_as_Si(*self.args)
  1087. return f._eval_nseries(x, n, logx)
  1088. elif nu.is_Integer and nu > 1:
  1089. f = self._eval_rewrite_as_Ei(*self.args)
  1090. return f._eval_nseries(x, n, logx)
  1091. return super()._eval_nseries(x, n, logx)
  1092. def _eval_aseries(self, n, args0, x, logx):
  1093. from sympy.series.order import Order
  1094. point = args0[1]
  1095. nu = self.args[0]
  1096. if point is S.Infinity:
  1097. z = self.args[1]
  1098. s = [S.NegativeOne**k * RisingFactorial(nu, k) / z**k for k in range(n)] + [Order(1/z**n, x)]
  1099. return (exp(-z)/z) * Add(*s)
  1100. return super(expint, self)._eval_aseries(n, args0, x, logx)
  1101. def _eval_rewrite_as_Integral(self, *args, **kwargs):
  1102. from sympy.integrals.integrals import Integral
  1103. n, x = self.args
  1104. t = Dummy(uniquely_named_symbol('t', args).name)
  1105. return Integral(t**-n * exp(-t*x), (t, 1, S.Infinity))
  1106. def E1(z):
  1107. """
  1108. Classical case of the generalized exponential integral.
  1109. Explanation
  1110. ===========
  1111. This is equivalent to ``expint(1, z)``.
  1112. Examples
  1113. ========
  1114. >>> from sympy import E1
  1115. >>> E1(0)
  1116. expint(1, 0)
  1117. >>> E1(5)
  1118. expint(1, 5)
  1119. See Also
  1120. ========
  1121. Ei: Exponential integral.
  1122. expint: Generalised exponential integral.
  1123. li: Logarithmic integral.
  1124. Li: Offset logarithmic integral.
  1125. Si: Sine integral.
  1126. Ci: Cosine integral.
  1127. Shi: Hyperbolic sine integral.
  1128. Chi: Hyperbolic cosine integral.
  1129. """
  1130. return expint(1, z)
  1131. class li(DefinedFunction):
  1132. r"""
  1133. The classical logarithmic integral.
  1134. Explanation
  1135. ===========
  1136. For use in SymPy, this function is defined as
  1137. .. math:: \operatorname{li}(x) = \int_0^x \frac{1}{\log(t)} \mathrm{d}t \,.
  1138. Examples
  1139. ========
  1140. >>> from sympy import I, oo, li
  1141. >>> from sympy.abc import z
  1142. Several special values are known:
  1143. >>> li(0)
  1144. 0
  1145. >>> li(1)
  1146. -oo
  1147. >>> li(oo)
  1148. oo
  1149. Differentiation with respect to $z$ is supported:
  1150. >>> from sympy import diff
  1151. >>> diff(li(z), z)
  1152. 1/log(z)
  1153. Defining the ``li`` function via an integral:
  1154. >>> from sympy import integrate
  1155. >>> integrate(li(z))
  1156. z*li(z) - Ei(2*log(z))
  1157. >>> integrate(li(z),z)
  1158. z*li(z) - Ei(2*log(z))
  1159. The logarithmic integral can also be defined in terms of ``Ei``:
  1160. >>> from sympy import Ei
  1161. >>> li(z).rewrite(Ei)
  1162. Ei(log(z))
  1163. >>> diff(li(z).rewrite(Ei), z)
  1164. 1/log(z)
  1165. We can numerically evaluate the logarithmic integral to arbitrary precision
  1166. on the whole complex plane (except the singular points):
  1167. >>> li(2).evalf(30)
  1168. 1.04516378011749278484458888919
  1169. >>> li(2*I).evalf(30)
  1170. 1.0652795784357498247001125598 + 3.08346052231061726610939702133*I
  1171. We can even compute Soldner's constant by the help of mpmath:
  1172. >>> from mpmath import findroot
  1173. >>> findroot(li, 2)
  1174. 1.45136923488338
  1175. Further transformations include rewriting ``li`` in terms of
  1176. the trigonometric integrals ``Si``, ``Ci``, ``Shi`` and ``Chi``:
  1177. >>> from sympy import Si, Ci, Shi, Chi
  1178. >>> li(z).rewrite(Si)
  1179. -log(I*log(z)) - log(1/log(z))/2 + log(log(z))/2 + Ci(I*log(z)) + Shi(log(z))
  1180. >>> li(z).rewrite(Ci)
  1181. -log(I*log(z)) - log(1/log(z))/2 + log(log(z))/2 + Ci(I*log(z)) + Shi(log(z))
  1182. >>> li(z).rewrite(Shi)
  1183. -log(1/log(z))/2 + log(log(z))/2 + Chi(log(z)) - Shi(log(z))
  1184. >>> li(z).rewrite(Chi)
  1185. -log(1/log(z))/2 + log(log(z))/2 + Chi(log(z)) - Shi(log(z))
  1186. See Also
  1187. ========
  1188. Li: Offset logarithmic integral.
  1189. Ei: Exponential integral.
  1190. expint: Generalised exponential integral.
  1191. E1: Special case of the generalised exponential integral.
  1192. Si: Sine integral.
  1193. Ci: Cosine integral.
  1194. Shi: Hyperbolic sine integral.
  1195. Chi: Hyperbolic cosine integral.
  1196. References
  1197. ==========
  1198. .. [1] https://en.wikipedia.org/wiki/Logarithmic_integral
  1199. .. [2] https://mathworld.wolfram.com/LogarithmicIntegral.html
  1200. .. [3] https://dlmf.nist.gov/6
  1201. .. [4] https://mathworld.wolfram.com/SoldnersConstant.html
  1202. """
  1203. @classmethod
  1204. def eval(cls, z):
  1205. if z.is_zero:
  1206. return S.Zero
  1207. elif z is S.One:
  1208. return S.NegativeInfinity
  1209. elif z is S.Infinity:
  1210. return S.Infinity
  1211. if z.is_zero:
  1212. return S.Zero
  1213. def fdiff(self, argindex=1):
  1214. arg = self.args[0]
  1215. if argindex == 1:
  1216. return S.One / log(arg)
  1217. else:
  1218. raise ArgumentIndexError(self, argindex)
  1219. def _eval_conjugate(self):
  1220. z = self.args[0]
  1221. # Exclude values on the branch cut (-oo, 0)
  1222. if not z.is_extended_negative:
  1223. return self.func(z.conjugate())
  1224. def _eval_rewrite_as_Li(self, z, **kwargs):
  1225. return Li(z) + li(2)
  1226. def _eval_rewrite_as_Ei(self, z, **kwargs):
  1227. return Ei(log(z))
  1228. def _eval_rewrite_as_uppergamma(self, z, **kwargs):
  1229. from sympy.functions.special.gamma_functions import uppergamma
  1230. return (-uppergamma(0, -log(z)) +
  1231. S.Half*(log(log(z)) - log(S.One/log(z))) - log(-log(z)))
  1232. def _eval_rewrite_as_Si(self, z, **kwargs):
  1233. return (Ci(I*log(z)) - I*Si(I*log(z)) -
  1234. S.Half*(log(S.One/log(z)) - log(log(z))) - log(I*log(z)))
  1235. _eval_rewrite_as_Ci = _eval_rewrite_as_Si
  1236. def _eval_rewrite_as_Shi(self, z, **kwargs):
  1237. return (Chi(log(z)) - Shi(log(z)) - S.Half*(log(S.One/log(z)) - log(log(z))))
  1238. _eval_rewrite_as_Chi = _eval_rewrite_as_Shi
  1239. def _eval_rewrite_as_hyper(self, z, **kwargs):
  1240. return (log(z)*hyper((1, 1), (2, 2), log(z)) +
  1241. S.Half*(log(log(z)) - log(S.One/log(z))) + EulerGamma)
  1242. def _eval_rewrite_as_meijerg(self, z, **kwargs):
  1243. return (-log(-log(z)) - S.Half*(log(S.One/log(z)) - log(log(z)))
  1244. - meijerg(((), (1,)), ((0, 0), ()), -log(z)))
  1245. def _eval_rewrite_as_tractable(self, z, limitvar=None, **kwargs):
  1246. return z * _eis(log(z))
  1247. def _eval_nseries(self, x, n, logx, cdir=0):
  1248. z = self.args[0]
  1249. s = [(log(z))**k / (factorial(k) * k) for k in range(1, n)]
  1250. return EulerGamma + log(log(z)) + Add(*s)
  1251. def _eval_is_zero(self):
  1252. z = self.args[0]
  1253. if z.is_zero:
  1254. return True
  1255. class Li(DefinedFunction):
  1256. r"""
  1257. The offset logarithmic integral.
  1258. Explanation
  1259. ===========
  1260. For use in SymPy, this function is defined as
  1261. .. math:: \operatorname{Li}(x) = \operatorname{li}(x) - \operatorname{li}(2)
  1262. Examples
  1263. ========
  1264. >>> from sympy import Li
  1265. >>> from sympy.abc import z
  1266. The following special value is known:
  1267. >>> Li(2)
  1268. 0
  1269. Differentiation with respect to $z$ is supported:
  1270. >>> from sympy import diff
  1271. >>> diff(Li(z), z)
  1272. 1/log(z)
  1273. The shifted logarithmic integral can be written in terms of $li(z)$:
  1274. >>> from sympy import li
  1275. >>> Li(z).rewrite(li)
  1276. li(z) - li(2)
  1277. We can numerically evaluate the logarithmic integral to arbitrary precision
  1278. on the whole complex plane (except the singular points):
  1279. >>> Li(2).evalf(30)
  1280. 0
  1281. >>> Li(4).evalf(30)
  1282. 1.92242131492155809316615998938
  1283. See Also
  1284. ========
  1285. li: Logarithmic integral.
  1286. Ei: Exponential integral.
  1287. expint: Generalised exponential integral.
  1288. E1: Special case of the generalised exponential integral.
  1289. Si: Sine integral.
  1290. Ci: Cosine integral.
  1291. Shi: Hyperbolic sine integral.
  1292. Chi: Hyperbolic cosine integral.
  1293. References
  1294. ==========
  1295. .. [1] https://en.wikipedia.org/wiki/Logarithmic_integral
  1296. .. [2] https://mathworld.wolfram.com/LogarithmicIntegral.html
  1297. .. [3] https://dlmf.nist.gov/6
  1298. """
  1299. @classmethod
  1300. def eval(cls, z):
  1301. if z is S.Infinity:
  1302. return S.Infinity
  1303. elif z == S(2):
  1304. return S.Zero
  1305. def fdiff(self, argindex=1):
  1306. arg = self.args[0]
  1307. if argindex == 1:
  1308. return S.One / log(arg)
  1309. else:
  1310. raise ArgumentIndexError(self, argindex)
  1311. def _eval_evalf(self, prec):
  1312. return self.rewrite(li).evalf(prec)
  1313. def _eval_rewrite_as_li(self, z, **kwargs):
  1314. return li(z) - li(2)
  1315. def _eval_rewrite_as_tractable(self, z, limitvar=None, **kwargs):
  1316. return self.rewrite(li).rewrite("tractable", deep=True)
  1317. def _eval_nseries(self, x, n, logx, cdir=0):
  1318. f = self._eval_rewrite_as_li(*self.args)
  1319. return f._eval_nseries(x, n, logx)
  1320. ###############################################################################
  1321. #################### TRIGONOMETRIC INTEGRALS ##################################
  1322. ###############################################################################
  1323. class TrigonometricIntegral(DefinedFunction):
  1324. """ Base class for trigonometric integrals. """
  1325. @classmethod
  1326. def eval(cls, z):
  1327. if z is S.Zero:
  1328. return cls._atzero
  1329. elif z is S.Infinity:
  1330. return cls._atinf()
  1331. elif z is S.NegativeInfinity:
  1332. return cls._atneginf()
  1333. if z.is_zero:
  1334. return cls._atzero
  1335. nz = z.extract_multiplicatively(polar_lift(I))
  1336. if nz is None and cls._trigfunc(0) == 0:
  1337. nz = z.extract_multiplicatively(I)
  1338. if nz is not None:
  1339. return cls._Ifactor(nz, 1)
  1340. nz = z.extract_multiplicatively(polar_lift(-I))
  1341. if nz is not None:
  1342. return cls._Ifactor(nz, -1)
  1343. nz = z.extract_multiplicatively(polar_lift(-1))
  1344. if nz is None and cls._trigfunc(0) == 0:
  1345. nz = z.extract_multiplicatively(-1)
  1346. if nz is not None:
  1347. return cls._minusfactor(nz)
  1348. nz, n = z.extract_branch_factor()
  1349. if n == 0 and nz == z:
  1350. return
  1351. return 2*pi*I*n*cls._trigfunc(0) + cls(nz)
  1352. def fdiff(self, argindex=1):
  1353. arg = unpolarify(self.args[0])
  1354. if argindex == 1:
  1355. return self._trigfunc(arg)/arg
  1356. else:
  1357. raise ArgumentIndexError(self, argindex)
  1358. def _eval_rewrite_as_Ei(self, z, **kwargs):
  1359. return self._eval_rewrite_as_expint(z).rewrite(Ei)
  1360. def _eval_rewrite_as_uppergamma(self, z, **kwargs):
  1361. from sympy.functions.special.gamma_functions import uppergamma
  1362. return self._eval_rewrite_as_expint(z).rewrite(uppergamma)
  1363. def _eval_nseries(self, x, n, logx, cdir=0):
  1364. # NOTE this is fairly inefficient
  1365. if self.args[0].subs(x, 0) != 0:
  1366. return super()._eval_nseries(x, n, logx)
  1367. baseseries = self._trigfunc(x)._eval_nseries(x, n, logx)
  1368. if self._trigfunc(0) != 0:
  1369. baseseries -= 1
  1370. baseseries = baseseries.replace(Pow, lambda t, n: t**n/n, simultaneous=False)
  1371. if self._trigfunc(0) != 0:
  1372. baseseries += EulerGamma + log(x)
  1373. return baseseries.subs(x, self.args[0])._eval_nseries(x, n, logx)
  1374. class Si(TrigonometricIntegral):
  1375. r"""
  1376. Sine integral.
  1377. Explanation
  1378. ===========
  1379. This function is defined by
  1380. .. math:: \operatorname{Si}(z) = \int_0^z \frac{\sin{t}}{t} \mathrm{d}t.
  1381. It is an entire function.
  1382. Examples
  1383. ========
  1384. >>> from sympy import Si
  1385. >>> from sympy.abc import z
  1386. The sine integral is an antiderivative of $sin(z)/z$:
  1387. >>> Si(z).diff(z)
  1388. sin(z)/z
  1389. It is unbranched:
  1390. >>> from sympy import exp_polar, I, pi
  1391. >>> Si(z*exp_polar(2*I*pi))
  1392. Si(z)
  1393. Sine integral behaves much like ordinary sine under multiplication by ``I``:
  1394. >>> Si(I*z)
  1395. I*Shi(z)
  1396. >>> Si(-z)
  1397. -Si(z)
  1398. It can also be expressed in terms of exponential integrals, but beware
  1399. that the latter is branched:
  1400. >>> from sympy import expint
  1401. >>> Si(z).rewrite(expint)
  1402. -I*(-expint(1, z*exp_polar(-I*pi/2))/2 +
  1403. expint(1, z*exp_polar(I*pi/2))/2) + pi/2
  1404. It can be rewritten in the form of sinc function (by definition):
  1405. >>> from sympy import sinc
  1406. >>> Si(z).rewrite(sinc)
  1407. Integral(sinc(_t), (_t, 0, z))
  1408. See Also
  1409. ========
  1410. Ci: Cosine integral.
  1411. Shi: Hyperbolic sine integral.
  1412. Chi: Hyperbolic cosine integral.
  1413. Ei: Exponential integral.
  1414. expint: Generalised exponential integral.
  1415. sinc: unnormalized sinc function
  1416. E1: Special case of the generalised exponential integral.
  1417. li: Logarithmic integral.
  1418. Li: Offset logarithmic integral.
  1419. References
  1420. ==========
  1421. .. [1] https://en.wikipedia.org/wiki/Trigonometric_integral
  1422. """
  1423. _trigfunc = sin
  1424. _atzero = S.Zero
  1425. @classmethod
  1426. def _atinf(cls):
  1427. return pi*S.Half
  1428. @classmethod
  1429. def _atneginf(cls):
  1430. return -pi*S.Half
  1431. @classmethod
  1432. def _minusfactor(cls, z):
  1433. return -Si(z)
  1434. @classmethod
  1435. def _Ifactor(cls, z, sign):
  1436. return I*Shi(z)*sign
  1437. def _eval_rewrite_as_expint(self, z, **kwargs):
  1438. # XXX should we polarify z?
  1439. return pi/2 + (E1(polar_lift(I)*z) - E1(polar_lift(-I)*z))/2/I
  1440. def _eval_rewrite_as_Integral(self, z, **kwargs):
  1441. from sympy.integrals.integrals import Integral
  1442. t = Dummy(uniquely_named_symbol('t', [z]).name)
  1443. return Integral(sinc(t), (t, 0, z))
  1444. _eval_rewrite_as_sinc = _eval_rewrite_as_Integral
  1445. def _eval_as_leading_term(self, x, logx, cdir):
  1446. arg = self.args[0].as_leading_term(x, logx=logx, cdir=cdir)
  1447. arg0 = arg.subs(x, 0)
  1448. if arg0 is S.NaN:
  1449. arg0 = arg.limit(x, 0, dir='-' if re(cdir).is_negative else '+')
  1450. if arg0.is_zero:
  1451. return arg
  1452. elif not arg0.is_infinite:
  1453. return self.func(arg0)
  1454. else:
  1455. return self
  1456. def _eval_aseries(self, n, args0, x, logx):
  1457. from sympy.series.order import Order
  1458. point = args0[0]
  1459. # Expansion at oo
  1460. if point is S.Infinity:
  1461. z = self.args[0]
  1462. p = [S.NegativeOne**k * factorial(2*k) / z**(2*k + 1)
  1463. for k in range(n//2 + 1)] + [Order(1/z**n, x)]
  1464. q = [S.NegativeOne**k * factorial(2*k + 1) / z**(2*(k + 1))
  1465. for k in range(n//2)] + [Order(1/z**n, x)]
  1466. return pi/2 - cos(z)*Add(*p) - sin(z)*Add(*q)
  1467. # All other points are not handled
  1468. return super(Si, self)._eval_aseries(n, args0, x, logx)
  1469. def _eval_is_zero(self):
  1470. z = self.args[0]
  1471. if z.is_zero:
  1472. return True
  1473. class Ci(TrigonometricIntegral):
  1474. r"""
  1475. Cosine integral.
  1476. Explanation
  1477. ===========
  1478. This function is defined for positive $x$ by
  1479. .. math:: \operatorname{Ci}(x) = \gamma + \log{x}
  1480. + \int_0^x \frac{\cos{t} - 1}{t} \mathrm{d}t
  1481. = -\int_x^\infty \frac{\cos{t}}{t} \mathrm{d}t,
  1482. where $\gamma$ is the Euler-Mascheroni constant.
  1483. We have
  1484. .. math:: \operatorname{Ci}(z) =
  1485. -\frac{\operatorname{E}_1\left(e^{i\pi/2} z\right)
  1486. + \operatorname{E}_1\left(e^{-i \pi/2} z\right)}{2}
  1487. which holds for all polar $z$ and thus provides an analytic
  1488. continuation to the Riemann surface of the logarithm.
  1489. The formula also holds as stated
  1490. for $z \in \mathbb{C}$ with $\Re(z) > 0$.
  1491. By lifting to the principal branch, we obtain an analytic function on the
  1492. cut complex plane.
  1493. Examples
  1494. ========
  1495. >>> from sympy import Ci
  1496. >>> from sympy.abc import z
  1497. The cosine integral is a primitive of $\cos(z)/z$:
  1498. >>> Ci(z).diff(z)
  1499. cos(z)/z
  1500. It has a logarithmic branch point at the origin:
  1501. >>> from sympy import exp_polar, I, pi
  1502. >>> Ci(z*exp_polar(2*I*pi))
  1503. Ci(z) + 2*I*pi
  1504. The cosine integral behaves somewhat like ordinary $\cos$ under
  1505. multiplication by $i$:
  1506. >>> from sympy import polar_lift
  1507. >>> Ci(polar_lift(I)*z)
  1508. Chi(z) + I*pi/2
  1509. >>> Ci(polar_lift(-1)*z)
  1510. Ci(z) + I*pi
  1511. It can also be expressed in terms of exponential integrals:
  1512. >>> from sympy import expint
  1513. >>> Ci(z).rewrite(expint)
  1514. -expint(1, z*exp_polar(-I*pi/2))/2 - expint(1, z*exp_polar(I*pi/2))/2
  1515. See Also
  1516. ========
  1517. Si: Sine integral.
  1518. Shi: Hyperbolic sine integral.
  1519. Chi: Hyperbolic cosine integral.
  1520. Ei: Exponential integral.
  1521. expint: Generalised exponential integral.
  1522. E1: Special case of the generalised exponential integral.
  1523. li: Logarithmic integral.
  1524. Li: Offset logarithmic integral.
  1525. References
  1526. ==========
  1527. .. [1] https://en.wikipedia.org/wiki/Trigonometric_integral
  1528. """
  1529. _trigfunc = cos
  1530. _atzero = S.ComplexInfinity
  1531. @classmethod
  1532. def _atinf(cls):
  1533. return S.Zero
  1534. @classmethod
  1535. def _atneginf(cls):
  1536. return I*pi
  1537. @classmethod
  1538. def _minusfactor(cls, z):
  1539. return Ci(z) + I*pi
  1540. @classmethod
  1541. def _Ifactor(cls, z, sign):
  1542. return Chi(z) + I*pi/2*sign
  1543. def _eval_rewrite_as_expint(self, z, **kwargs):
  1544. return -(E1(polar_lift(I)*z) + E1(polar_lift(-I)*z))/2
  1545. def _eval_rewrite_as_Integral(self, z, **kwargs):
  1546. from sympy.integrals.integrals import Integral
  1547. t = Dummy(uniquely_named_symbol('t', [z]).name)
  1548. return S.EulerGamma + log(z) - Integral((1-cos(t))/t, (t, 0, z))
  1549. def _eval_as_leading_term(self, x, logx, cdir):
  1550. arg = self.args[0].as_leading_term(x, logx=logx, cdir=cdir)
  1551. arg0 = arg.subs(x, 0)
  1552. if arg0 is S.NaN:
  1553. arg0 = arg.limit(x, 0, dir='-' if re(cdir).is_negative else '+')
  1554. if arg0.is_zero:
  1555. c, e = arg.as_coeff_exponent(x)
  1556. logx = log(x) if logx is None else logx
  1557. return log(c) + e*logx + EulerGamma
  1558. elif arg0.is_finite:
  1559. return self.func(arg0)
  1560. else:
  1561. return self
  1562. def _eval_aseries(self, n, args0, x, logx):
  1563. from sympy.series.order import Order
  1564. point = args0[0]
  1565. if point in (S.Infinity, S.NegativeInfinity):
  1566. z = self.args[0]
  1567. p = [S.NegativeOne**k * factorial(2*k) / z**(2*k + 1)
  1568. for k in range(n//2 + 1)] + [Order(1/z**n, x)]
  1569. q = [S.NegativeOne**k * factorial(2*k + 1) / z**(2*(k + 1))
  1570. for k in range(n//2)] + [Order(1/z**n, x)]
  1571. result = sin(z)*(Add(*p)) - cos(z)*(Add(*q))
  1572. if point is S.NegativeInfinity:
  1573. result += I*pi
  1574. return result
  1575. return super(Ci, self)._eval_aseries(n, args0, x, logx)
  1576. class Shi(TrigonometricIntegral):
  1577. r"""
  1578. Sinh integral.
  1579. Explanation
  1580. ===========
  1581. This function is defined by
  1582. .. math:: \operatorname{Shi}(z) = \int_0^z \frac{\sinh{t}}{t} \mathrm{d}t.
  1583. It is an entire function.
  1584. Examples
  1585. ========
  1586. >>> from sympy import Shi
  1587. >>> from sympy.abc import z
  1588. The Sinh integral is a primitive of $\sinh(z)/z$:
  1589. >>> Shi(z).diff(z)
  1590. sinh(z)/z
  1591. It is unbranched:
  1592. >>> from sympy import exp_polar, I, pi
  1593. >>> Shi(z*exp_polar(2*I*pi))
  1594. Shi(z)
  1595. The $\sinh$ integral behaves much like ordinary $\sinh$ under
  1596. multiplication by $i$:
  1597. >>> Shi(I*z)
  1598. I*Si(z)
  1599. >>> Shi(-z)
  1600. -Shi(z)
  1601. It can also be expressed in terms of exponential integrals, but beware
  1602. that the latter is branched:
  1603. >>> from sympy import expint
  1604. >>> Shi(z).rewrite(expint)
  1605. expint(1, z)/2 - expint(1, z*exp_polar(I*pi))/2 - I*pi/2
  1606. See Also
  1607. ========
  1608. Si: Sine integral.
  1609. Ci: Cosine integral.
  1610. Chi: Hyperbolic cosine integral.
  1611. Ei: Exponential integral.
  1612. expint: Generalised exponential integral.
  1613. E1: Special case of the generalised exponential integral.
  1614. li: Logarithmic integral.
  1615. Li: Offset logarithmic integral.
  1616. References
  1617. ==========
  1618. .. [1] https://en.wikipedia.org/wiki/Trigonometric_integral
  1619. """
  1620. _trigfunc = sinh
  1621. _atzero = S.Zero
  1622. @classmethod
  1623. def _atinf(cls):
  1624. return S.Infinity
  1625. @classmethod
  1626. def _atneginf(cls):
  1627. return S.NegativeInfinity
  1628. @classmethod
  1629. def _minusfactor(cls, z):
  1630. return -Shi(z)
  1631. @classmethod
  1632. def _Ifactor(cls, z, sign):
  1633. return I*Si(z)*sign
  1634. def _eval_rewrite_as_expint(self, z, **kwargs):
  1635. # XXX should we polarify z?
  1636. return (E1(z) - E1(exp_polar(I*pi)*z))/2 - I*pi/2
  1637. def _eval_is_zero(self):
  1638. z = self.args[0]
  1639. if z.is_zero:
  1640. return True
  1641. def _eval_as_leading_term(self, x, logx, cdir):
  1642. arg = self.args[0].as_leading_term(x)
  1643. arg0 = arg.subs(x, 0)
  1644. if arg0 is S.NaN:
  1645. arg0 = arg.limit(x, 0, dir='-' if re(cdir).is_negative else '+')
  1646. if arg0.is_zero:
  1647. return arg
  1648. elif not arg0.is_infinite:
  1649. return self.func(arg0)
  1650. else:
  1651. return self
  1652. class Chi(TrigonometricIntegral):
  1653. r"""
  1654. Cosh integral.
  1655. Explanation
  1656. ===========
  1657. This function is defined for positive $x$ by
  1658. .. math:: \operatorname{Chi}(x) = \gamma + \log{x}
  1659. + \int_0^x \frac{\cosh{t} - 1}{t} \mathrm{d}t,
  1660. where $\gamma$ is the Euler-Mascheroni constant.
  1661. We have
  1662. .. math:: \operatorname{Chi}(z) = \operatorname{Ci}\left(e^{i \pi/2}z\right)
  1663. - i\frac{\pi}{2},
  1664. which holds for all polar $z$ and thus provides an analytic
  1665. continuation to the Riemann surface of the logarithm.
  1666. By lifting to the principal branch we obtain an analytic function on the
  1667. cut complex plane.
  1668. Examples
  1669. ========
  1670. >>> from sympy import Chi
  1671. >>> from sympy.abc import z
  1672. The $\cosh$ integral is a primitive of $\cosh(z)/z$:
  1673. >>> Chi(z).diff(z)
  1674. cosh(z)/z
  1675. It has a logarithmic branch point at the origin:
  1676. >>> from sympy import exp_polar, I, pi
  1677. >>> Chi(z*exp_polar(2*I*pi))
  1678. Chi(z) + 2*I*pi
  1679. The $\cosh$ integral behaves somewhat like ordinary $\cosh$ under
  1680. multiplication by $i$:
  1681. >>> from sympy import polar_lift
  1682. >>> Chi(polar_lift(I)*z)
  1683. Ci(z) + I*pi/2
  1684. >>> Chi(polar_lift(-1)*z)
  1685. Chi(z) + I*pi
  1686. It can also be expressed in terms of exponential integrals:
  1687. >>> from sympy import expint
  1688. >>> Chi(z).rewrite(expint)
  1689. -expint(1, z)/2 - expint(1, z*exp_polar(I*pi))/2 - I*pi/2
  1690. See Also
  1691. ========
  1692. Si: Sine integral.
  1693. Ci: Cosine integral.
  1694. Shi: Hyperbolic sine integral.
  1695. Ei: Exponential integral.
  1696. expint: Generalised exponential integral.
  1697. E1: Special case of the generalised exponential integral.
  1698. li: Logarithmic integral.
  1699. Li: Offset logarithmic integral.
  1700. References
  1701. ==========
  1702. .. [1] https://en.wikipedia.org/wiki/Trigonometric_integral
  1703. """
  1704. _trigfunc = cosh
  1705. _atzero = S.ComplexInfinity
  1706. @classmethod
  1707. def _atinf(cls):
  1708. return S.Infinity
  1709. @classmethod
  1710. def _atneginf(cls):
  1711. return S.Infinity
  1712. @classmethod
  1713. def _minusfactor(cls, z):
  1714. return Chi(z) + I*pi
  1715. @classmethod
  1716. def _Ifactor(cls, z, sign):
  1717. return Ci(z) + I*pi/2*sign
  1718. def _eval_rewrite_as_expint(self, z, **kwargs):
  1719. return -I*pi/2 - (E1(z) + E1(exp_polar(I*pi)*z))/2
  1720. def _eval_as_leading_term(self, x, logx, cdir):
  1721. arg = self.args[0].as_leading_term(x, logx=logx, cdir=cdir)
  1722. arg0 = arg.subs(x, 0)
  1723. if arg0 is S.NaN:
  1724. arg0 = arg.limit(x, 0, dir='-' if re(cdir).is_negative else '+')
  1725. if arg0.is_zero:
  1726. c, e = arg.as_coeff_exponent(x)
  1727. logx = log(x) if logx is None else logx
  1728. return log(c) + e*logx + EulerGamma
  1729. elif arg0.is_finite:
  1730. return self.func(arg0)
  1731. else:
  1732. return self
  1733. ###############################################################################
  1734. #################### FRESNEL INTEGRALS ########################################
  1735. ###############################################################################
  1736. class FresnelIntegral(DefinedFunction):
  1737. """ Base class for the Fresnel integrals."""
  1738. unbranched = True
  1739. @classmethod
  1740. def eval(cls, z):
  1741. # Values at positive infinities signs
  1742. # if any were extracted automatically
  1743. if z is S.Infinity:
  1744. return S.Half
  1745. # Value at zero
  1746. if z.is_zero:
  1747. return S.Zero
  1748. # Try to pull out factors of -1 and I
  1749. prefact = S.One
  1750. newarg = z
  1751. changed = False
  1752. nz = newarg.extract_multiplicatively(-1)
  1753. if nz is not None:
  1754. prefact = -prefact
  1755. newarg = nz
  1756. changed = True
  1757. nz = newarg.extract_multiplicatively(I)
  1758. if nz is not None:
  1759. prefact = cls._sign*I*prefact
  1760. newarg = nz
  1761. changed = True
  1762. if changed:
  1763. return prefact*cls(newarg)
  1764. def fdiff(self, argindex=1):
  1765. if argindex == 1:
  1766. return self._trigfunc(S.Half*pi*self.args[0]**2)
  1767. else:
  1768. raise ArgumentIndexError(self, argindex)
  1769. def _eval_is_extended_real(self):
  1770. return self.args[0].is_extended_real
  1771. _eval_is_finite = _eval_is_extended_real
  1772. def _eval_is_zero(self):
  1773. return self.args[0].is_zero
  1774. def _eval_conjugate(self):
  1775. return self.func(self.args[0].conjugate())
  1776. as_real_imag = real_to_real_as_real_imag
  1777. class fresnels(FresnelIntegral):
  1778. r"""
  1779. Fresnel integral S.
  1780. Explanation
  1781. ===========
  1782. This function is defined by
  1783. .. math:: \operatorname{S}(z) = \int_0^z \sin{\frac{\pi}{2} t^2} \mathrm{d}t.
  1784. It is an entire function.
  1785. Examples
  1786. ========
  1787. >>> from sympy import I, oo, fresnels
  1788. >>> from sympy.abc import z
  1789. Several special values are known:
  1790. >>> fresnels(0)
  1791. 0
  1792. >>> fresnels(oo)
  1793. 1/2
  1794. >>> fresnels(-oo)
  1795. -1/2
  1796. >>> fresnels(I*oo)
  1797. -I/2
  1798. >>> fresnels(-I*oo)
  1799. I/2
  1800. In general one can pull out factors of -1 and $i$ from the argument:
  1801. >>> fresnels(-z)
  1802. -fresnels(z)
  1803. >>> fresnels(I*z)
  1804. -I*fresnels(z)
  1805. The Fresnel S integral obeys the mirror symmetry
  1806. $\overline{S(z)} = S(\bar{z})$:
  1807. >>> from sympy import conjugate
  1808. >>> conjugate(fresnels(z))
  1809. fresnels(conjugate(z))
  1810. Differentiation with respect to $z$ is supported:
  1811. >>> from sympy import diff
  1812. >>> diff(fresnels(z), z)
  1813. sin(pi*z**2/2)
  1814. Defining the Fresnel functions via an integral:
  1815. >>> from sympy import integrate, pi, sin, expand_func
  1816. >>> integrate(sin(pi*z**2/2), z)
  1817. 3*fresnels(z)*gamma(3/4)/(4*gamma(7/4))
  1818. >>> expand_func(integrate(sin(pi*z**2/2), z))
  1819. fresnels(z)
  1820. We can numerically evaluate the Fresnel integral to arbitrary precision
  1821. on the whole complex plane:
  1822. >>> fresnels(2).evalf(30)
  1823. 0.343415678363698242195300815958
  1824. >>> fresnels(-2*I).evalf(30)
  1825. 0.343415678363698242195300815958*I
  1826. See Also
  1827. ========
  1828. fresnelc: Fresnel cosine integral.
  1829. References
  1830. ==========
  1831. .. [1] https://en.wikipedia.org/wiki/Fresnel_integral
  1832. .. [2] https://dlmf.nist.gov/7
  1833. .. [3] https://mathworld.wolfram.com/FresnelIntegrals.html
  1834. .. [4] https://functions.wolfram.com/GammaBetaErf/FresnelS
  1835. .. [5] The converging factors for the fresnel integrals
  1836. by John W. Wrench Jr. and Vicki Alley
  1837. """
  1838. _trigfunc = sin
  1839. _sign = -S.One
  1840. @staticmethod
  1841. @cacheit
  1842. def taylor_term(n, x, *previous_terms):
  1843. if n < 0:
  1844. return S.Zero
  1845. else:
  1846. x = sympify(x)
  1847. if len(previous_terms) > 1:
  1848. p = previous_terms[-1]
  1849. return (-pi**2*x**4*(4*n - 1)/(8*n*(2*n + 1)*(4*n + 3))) * p
  1850. else:
  1851. return x**3 * (-x**4)**n * (S(2)**(-2*n - 1)*pi**(2*n + 1)) / ((4*n + 3)*factorial(2*n + 1))
  1852. def _eval_rewrite_as_erf(self, z, **kwargs):
  1853. return (S.One + I)/4 * (erf((S.One + I)/2*sqrt(pi)*z) - I*erf((S.One - I)/2*sqrt(pi)*z))
  1854. def _eval_rewrite_as_hyper(self, z, **kwargs):
  1855. return pi*z**3/6 * hyper([Rational(3, 4)], [Rational(3, 2), Rational(7, 4)], -pi**2*z**4/16)
  1856. def _eval_rewrite_as_meijerg(self, z, **kwargs):
  1857. return (pi*z**Rational(9, 4) / (sqrt(2)*(z**2)**Rational(3, 4)*(-z)**Rational(3, 4))
  1858. * meijerg([], [1], [Rational(3, 4)], [Rational(1, 4), 0], -pi**2*z**4/16))
  1859. def _eval_rewrite_as_Integral(self, z, **kwargs):
  1860. from sympy.integrals.integrals import Integral
  1861. t = Dummy(uniquely_named_symbol('t', [z]).name)
  1862. return Integral(sin(pi*t**2/2), (t, 0, z))
  1863. def _eval_as_leading_term(self, x, logx, cdir):
  1864. from sympy.series.order import Order
  1865. arg = self.args[0].as_leading_term(x, logx=logx, cdir=cdir)
  1866. arg0 = arg.subs(x, 0)
  1867. if arg0 is S.ComplexInfinity:
  1868. arg0 = arg.limit(x, 0, dir='-' if re(cdir).is_negative else '+')
  1869. if arg0.is_zero:
  1870. return pi*arg**3/6
  1871. elif arg0 in [S.Infinity, S.NegativeInfinity]:
  1872. s = 1 if arg0 is S.Infinity else -1
  1873. return s*S.Half + Order(x, x)
  1874. else:
  1875. return self.func(arg0)
  1876. def _eval_aseries(self, n, args0, x, logx):
  1877. from sympy.series.order import Order
  1878. point = args0[0]
  1879. # Expansion at oo and -oo
  1880. if point in [S.Infinity, -S.Infinity]:
  1881. z = self.args[0]
  1882. # expansion of S(x) = S1(x*sqrt(pi/2)), see reference[5] page 1-8
  1883. # as only real infinities are dealt with, sin and cos are O(1)
  1884. p = [S.NegativeOne**k * factorial(4*k + 1) /
  1885. (2**(2*k + 2) * z**(4*k + 3) * 2**(2*k)*factorial(2*k))
  1886. for k in range(0, n) if 4*k + 3 < n]
  1887. q = [1/(2*z)] + [S.NegativeOne**k * factorial(4*k - 1) /
  1888. (2**(2*k + 1) * z**(4*k + 1) * 2**(2*k - 1)*factorial(2*k - 1))
  1889. for k in range(1, n) if 4*k + 1 < n]
  1890. p = [-sqrt(2/pi)*t for t in p]
  1891. q = [-sqrt(2/pi)*t for t in q]
  1892. s = 1 if point is S.Infinity else -1
  1893. # The expansion at oo is 1/2 + some odd powers of z
  1894. # To get the expansion at -oo, replace z by -z and flip the sign
  1895. # The result -1/2 + the same odd powers of z as before.
  1896. return s*S.Half + (sin(z**2)*Add(*p) + cos(z**2)*Add(*q)
  1897. ).subs(x, sqrt(2/pi)*x) + Order(1/z**n, x)
  1898. # All other points are not handled
  1899. return super()._eval_aseries(n, args0, x, logx)
  1900. class fresnelc(FresnelIntegral):
  1901. r"""
  1902. Fresnel integral C.
  1903. Explanation
  1904. ===========
  1905. This function is defined by
  1906. .. math:: \operatorname{C}(z) = \int_0^z \cos{\frac{\pi}{2} t^2} \mathrm{d}t.
  1907. It is an entire function.
  1908. Examples
  1909. ========
  1910. >>> from sympy import I, oo, fresnelc
  1911. >>> from sympy.abc import z
  1912. Several special values are known:
  1913. >>> fresnelc(0)
  1914. 0
  1915. >>> fresnelc(oo)
  1916. 1/2
  1917. >>> fresnelc(-oo)
  1918. -1/2
  1919. >>> fresnelc(I*oo)
  1920. I/2
  1921. >>> fresnelc(-I*oo)
  1922. -I/2
  1923. In general one can pull out factors of -1 and $i$ from the argument:
  1924. >>> fresnelc(-z)
  1925. -fresnelc(z)
  1926. >>> fresnelc(I*z)
  1927. I*fresnelc(z)
  1928. The Fresnel C integral obeys the mirror symmetry
  1929. $\overline{C(z)} = C(\bar{z})$:
  1930. >>> from sympy import conjugate
  1931. >>> conjugate(fresnelc(z))
  1932. fresnelc(conjugate(z))
  1933. Differentiation with respect to $z$ is supported:
  1934. >>> from sympy import diff
  1935. >>> diff(fresnelc(z), z)
  1936. cos(pi*z**2/2)
  1937. Defining the Fresnel functions via an integral:
  1938. >>> from sympy import integrate, pi, cos, expand_func
  1939. >>> integrate(cos(pi*z**2/2), z)
  1940. fresnelc(z)*gamma(1/4)/(4*gamma(5/4))
  1941. >>> expand_func(integrate(cos(pi*z**2/2), z))
  1942. fresnelc(z)
  1943. We can numerically evaluate the Fresnel integral to arbitrary precision
  1944. on the whole complex plane:
  1945. >>> fresnelc(2).evalf(30)
  1946. 0.488253406075340754500223503357
  1947. >>> fresnelc(-2*I).evalf(30)
  1948. -0.488253406075340754500223503357*I
  1949. See Also
  1950. ========
  1951. fresnels: Fresnel sine integral.
  1952. References
  1953. ==========
  1954. .. [1] https://en.wikipedia.org/wiki/Fresnel_integral
  1955. .. [2] https://dlmf.nist.gov/7
  1956. .. [3] https://mathworld.wolfram.com/FresnelIntegrals.html
  1957. .. [4] https://functions.wolfram.com/GammaBetaErf/FresnelC
  1958. .. [5] The converging factors for the fresnel integrals
  1959. by John W. Wrench Jr. and Vicki Alley
  1960. """
  1961. _trigfunc = cos
  1962. _sign = S.One
  1963. @staticmethod
  1964. @cacheit
  1965. def taylor_term(n, x, *previous_terms):
  1966. if n < 0:
  1967. return S.Zero
  1968. else:
  1969. x = sympify(x)
  1970. if len(previous_terms) > 1:
  1971. p = previous_terms[-1]
  1972. return (-pi**2*x**4*(4*n - 3)/(8*n*(2*n - 1)*(4*n + 1))) * p
  1973. else:
  1974. return x * (-x**4)**n * (S(2)**(-2*n)*pi**(2*n)) / ((4*n + 1)*factorial(2*n))
  1975. def _eval_rewrite_as_erf(self, z, **kwargs):
  1976. return (S.One - I)/4 * (erf((S.One + I)/2*sqrt(pi)*z) + I*erf((S.One - I)/2*sqrt(pi)*z))
  1977. def _eval_rewrite_as_hyper(self, z, **kwargs):
  1978. return z * hyper([Rational(1, 4)], [S.Half, Rational(5, 4)], -pi**2*z**4/16)
  1979. def _eval_rewrite_as_meijerg(self, z, **kwargs):
  1980. return (pi*z**Rational(3, 4) / (sqrt(2)*root(z**2, 4)*root(-z, 4))
  1981. * meijerg([], [1], [Rational(1, 4)], [Rational(3, 4), 0], -pi**2*z**4/16))
  1982. def _eval_rewrite_as_Integral(self, z, **kwargs):
  1983. from sympy.integrals.integrals import Integral
  1984. t = Dummy(uniquely_named_symbol('t', [z]).name)
  1985. return Integral(cos(pi*t**2/2), (t, 0, z))
  1986. def _eval_as_leading_term(self, x, logx, cdir):
  1987. from sympy.series.order import Order
  1988. arg = self.args[0].as_leading_term(x, logx=logx, cdir=cdir)
  1989. arg0 = arg.subs(x, 0)
  1990. if arg0 is S.ComplexInfinity:
  1991. arg0 = arg.limit(x, 0, dir='-' if re(cdir).is_negative else '+')
  1992. if arg0.is_zero:
  1993. return arg
  1994. elif arg0 in [S.Infinity, S.NegativeInfinity]:
  1995. s = 1 if arg0 is S.Infinity else -1
  1996. return s*S.Half + Order(x, x)
  1997. else:
  1998. return self.func(arg0)
  1999. def _eval_aseries(self, n, args0, x, logx):
  2000. from sympy.series.order import Order
  2001. point = args0[0]
  2002. # Expansion at oo
  2003. if point in [S.Infinity, -S.Infinity]:
  2004. z = self.args[0]
  2005. # expansion of C(x) = C1(x*sqrt(pi/2)), see reference[5] page 1-8
  2006. # as only real infinities are dealt with, sin and cos are O(1)
  2007. p = [S.NegativeOne**k * factorial(4*k + 1) /
  2008. (2**(2*k + 2) * z**(4*k + 3) * 2**(2*k)*factorial(2*k))
  2009. for k in range(n) if 4*k + 3 < n]
  2010. q = [1/(2*z)] + [S.NegativeOne**k * factorial(4*k - 1) /
  2011. (2**(2*k + 1) * z**(4*k + 1) * 2**(2*k - 1)*factorial(2*k - 1))
  2012. for k in range(1, n) if 4*k + 1 < n]
  2013. p = [-sqrt(2/pi)*t for t in p]
  2014. q = [ sqrt(2/pi)*t for t in q]
  2015. s = 1 if point is S.Infinity else -1
  2016. # The expansion at oo is 1/2 + some odd powers of z
  2017. # To get the expansion at -oo, replace z by -z and flip the sign
  2018. # The result -1/2 + the same odd powers of z as before.
  2019. return s*S.Half + (cos(z**2)*Add(*p) + sin(z**2)*Add(*q)
  2020. ).subs(x, sqrt(2/pi)*x) + Order(1/z**n, x)
  2021. # All other points are not handled
  2022. return super()._eval_aseries(n, args0, x, logx)
  2023. ###############################################################################
  2024. #################### HELPER FUNCTIONS #########################################
  2025. ###############################################################################
  2026. class _erfs(DefinedFunction):
  2027. """
  2028. Helper function to make the $\\mathrm{erf}(z)$ function
  2029. tractable for the Gruntz algorithm.
  2030. """
  2031. @classmethod
  2032. def eval(cls, arg):
  2033. if arg.is_zero:
  2034. return S.One
  2035. def _eval_aseries(self, n, args0, x, logx):
  2036. from sympy.series.order import Order
  2037. point = args0[0]
  2038. # Expansion at oo
  2039. if point is S.Infinity:
  2040. z = self.args[0]
  2041. l = [1/sqrt(pi) * factorial(2*k)*(-S(
  2042. 4))**(-k)/factorial(k) * (1/z)**(2*k + 1) for k in range(n)]
  2043. o = Order(1/z**(2*n + 1), x)
  2044. # It is very inefficient to first add the order and then do the nseries
  2045. return (Add(*l))._eval_nseries(x, n, logx) + o
  2046. # Expansion at I*oo
  2047. t = point.extract_multiplicatively(I)
  2048. if t is S.Infinity:
  2049. z = self.args[0]
  2050. # TODO: is the series really correct?
  2051. l = [1/sqrt(pi) * factorial(2*k)*(-S(
  2052. 4))**(-k)/factorial(k) * (1/z)**(2*k + 1) for k in range(n)]
  2053. o = Order(1/z**(2*n + 1), x)
  2054. # It is very inefficient to first add the order and then do the nseries
  2055. return (Add(*l))._eval_nseries(x, n, logx) + o
  2056. # All other points are not handled
  2057. return super()._eval_aseries(n, args0, x, logx)
  2058. def fdiff(self, argindex=1):
  2059. if argindex == 1:
  2060. z = self.args[0]
  2061. return -2/sqrt(pi) + 2*z*_erfs(z)
  2062. else:
  2063. raise ArgumentIndexError(self, argindex)
  2064. def _eval_rewrite_as_intractable(self, z, **kwargs):
  2065. return (S.One - erf(z))*exp(z**2)
  2066. class _eis(DefinedFunction):
  2067. """
  2068. Helper function to make the $\\mathrm{Ei}(z)$ and $\\mathrm{li}(z)$
  2069. functions tractable for the Gruntz algorithm.
  2070. """
  2071. def _eval_aseries(self, n, args0, x, logx):
  2072. from sympy.series.order import Order
  2073. if args0[0] not in (S.Infinity, S.NegativeInfinity):
  2074. return super()._eval_aseries(n, args0, x, logx)
  2075. z = self.args[0]
  2076. l = [factorial(k) * (1/z)**(k + 1) for k in range(n)]
  2077. o = Order(1/z**(n + 1), x)
  2078. # It is very inefficient to first add the order and then do the nseries
  2079. return (Add(*l))._eval_nseries(x, n, logx) + o
  2080. def fdiff(self, argindex=1):
  2081. if argindex == 1:
  2082. z = self.args[0]
  2083. return S.One / z - _eis(z)
  2084. else:
  2085. raise ArgumentIndexError(self, argindex)
  2086. def _eval_rewrite_as_intractable(self, z, **kwargs):
  2087. return exp(-z)*Ei(z)
  2088. def _eval_as_leading_term(self, x, logx, cdir):
  2089. x0 = self.args[0].limit(x, 0)
  2090. if x0.is_zero:
  2091. f = self._eval_rewrite_as_intractable(*self.args)
  2092. return f._eval_as_leading_term(x, logx=logx, cdir=cdir)
  2093. return super()._eval_as_leading_term(x, logx=logx, cdir=cdir)
  2094. def _eval_nseries(self, x, n, logx, cdir=0):
  2095. x0 = self.args[0].limit(x, 0)
  2096. if x0.is_zero:
  2097. f = self._eval_rewrite_as_intractable(*self.args)
  2098. return f._eval_nseries(x, n, logx)
  2099. return super()._eval_nseries(x, n, logx)