meijerint.py 78 KB

12345678910111213141516171819202122232425262728293031323334353637383940414243444546474849505152535455565758596061626364656667686970717273747576777879808182838485868788899091929394959697989910010110210310410510610710810911011111211311411511611711811912012112212312412512612712812913013113213313413513613713813914014114214314414514614714814915015115215315415515615715815916016116216316416516616716816917017117217317417517617717817918018118218318418518618718818919019119219319419519619719819920020120220320420520620720820921021121221321421521621721821922022122222322422522622722822923023123223323423523623723823924024124224324424524624724824925025125225325425525625725825926026126226326426526626726826927027127227327427527627727827928028128228328428528628728828929029129229329429529629729829930030130230330430530630730830931031131231331431531631731831932032132232332432532632732832933033133233333433533633733833934034134234334434534634734834935035135235335435535635735835936036136236336436536636736836937037137237337437537637737837938038138238338438538638738838939039139239339439539639739839940040140240340440540640740840941041141241341441541641741841942042142242342442542642742842943043143243343443543643743843944044144244344444544644744844945045145245345445545645745845946046146246346446546646746846947047147247347447547647747847948048148248348448548648748848949049149249349449549649749849950050150250350450550650750850951051151251351451551651751851952052152252352452552652752852953053153253353453553653753853954054154254354454554654754854955055155255355455555655755855956056156256356456556656756856957057157257357457557657757857958058158258358458558658758858959059159259359459559659759859960060160260360460560660760860961061161261361461561661761861962062162262362462562662762862963063163263363463563663763863964064164264364464564664764864965065165265365465565665765865966066166266366466566666766866967067167267367467567667767867968068168268368468568668768868969069169269369469569669769869970070170270370470570670770870971071171271371471571671771871972072172272372472572672772872973073173273373473573673773873974074174274374474574674774874975075175275375475575675775875976076176276376476576676776876977077177277377477577677777877978078178278378478578678778878979079179279379479579679779879980080180280380480580680780880981081181281381481581681781881982082182282382482582682782882983083183283383483583683783883984084184284384484584684784884985085185285385485585685785885986086186286386486586686786886987087187287387487587687787887988088188288388488588688788888989089189289389489589689789889990090190290390490590690790890991091191291391491591691791891992092192292392492592692792892993093193293393493593693793893994094194294394494594694794894995095195295395495595695795895996096196296396496596696796896997097197297397497597697797897998098198298398498598698798898999099199299399499599699799899910001001100210031004100510061007100810091010101110121013101410151016101710181019102010211022102310241025102610271028102910301031103210331034103510361037103810391040104110421043104410451046104710481049105010511052105310541055105610571058105910601061106210631064106510661067106810691070107110721073107410751076107710781079108010811082108310841085108610871088108910901091109210931094109510961097109810991100110111021103110411051106110711081109111011111112111311141115111611171118111911201121112211231124112511261127112811291130113111321133113411351136113711381139114011411142114311441145114611471148114911501151115211531154115511561157115811591160116111621163116411651166116711681169117011711172117311741175117611771178117911801181118211831184118511861187118811891190119111921193119411951196119711981199120012011202120312041205120612071208120912101211121212131214121512161217121812191220122112221223122412251226122712281229123012311232123312341235123612371238123912401241124212431244124512461247124812491250125112521253125412551256125712581259126012611262126312641265126612671268126912701271127212731274127512761277127812791280128112821283128412851286128712881289129012911292129312941295129612971298129913001301130213031304130513061307130813091310131113121313131413151316131713181319132013211322132313241325132613271328132913301331133213331334133513361337133813391340134113421343134413451346134713481349135013511352135313541355135613571358135913601361136213631364136513661367136813691370137113721373137413751376137713781379138013811382138313841385138613871388138913901391139213931394139513961397139813991400140114021403140414051406140714081409141014111412141314141415141614171418141914201421142214231424142514261427142814291430143114321433143414351436143714381439144014411442144314441445144614471448144914501451145214531454145514561457145814591460146114621463146414651466146714681469147014711472147314741475147614771478147914801481148214831484148514861487148814891490149114921493149414951496149714981499150015011502150315041505150615071508150915101511151215131514151515161517151815191520152115221523152415251526152715281529153015311532153315341535153615371538153915401541154215431544154515461547154815491550155115521553155415551556155715581559156015611562156315641565156615671568156915701571157215731574157515761577157815791580158115821583158415851586158715881589159015911592159315941595159615971598159916001601160216031604160516061607160816091610161116121613161416151616161716181619162016211622162316241625162616271628162916301631163216331634163516361637163816391640164116421643164416451646164716481649165016511652165316541655165616571658165916601661166216631664166516661667166816691670167116721673167416751676167716781679168016811682168316841685168616871688168916901691169216931694169516961697169816991700170117021703170417051706170717081709171017111712171317141715171617171718171917201721172217231724172517261727172817291730173117321733173417351736173717381739174017411742174317441745174617471748174917501751175217531754175517561757175817591760176117621763176417651766176717681769177017711772177317741775177617771778177917801781178217831784178517861787178817891790179117921793179417951796179717981799180018011802180318041805180618071808180918101811181218131814181518161817181818191820182118221823182418251826182718281829183018311832183318341835183618371838183918401841184218431844184518461847184818491850185118521853185418551856185718581859186018611862186318641865186618671868186918701871187218731874187518761877187818791880188118821883188418851886188718881889189018911892189318941895189618971898189919001901190219031904190519061907190819091910191119121913191419151916191719181919192019211922192319241925192619271928192919301931193219331934193519361937193819391940194119421943194419451946194719481949195019511952195319541955195619571958195919601961196219631964196519661967196819691970197119721973197419751976197719781979198019811982198319841985198619871988198919901991199219931994199519961997199819992000200120022003200420052006200720082009201020112012201320142015201620172018201920202021202220232024202520262027202820292030203120322033203420352036203720382039204020412042204320442045204620472048204920502051205220532054205520562057205820592060206120622063206420652066206720682069207020712072207320742075207620772078207920802081208220832084208520862087208820892090209120922093209420952096209720982099210021012102210321042105210621072108210921102111211221132114211521162117211821192120212121222123212421252126212721282129213021312132213321342135213621372138213921402141214221432144214521462147214821492150215121522153215421552156215721582159216021612162216321642165216621672168216921702171217221732174217521762177217821792180218121822183218421852186218721882189
  1. """
  2. Integrate functions by rewriting them as Meijer G-functions.
  3. There are three user-visible functions that can be used by other parts of the
  4. sympy library to solve various integration problems:
  5. - meijerint_indefinite
  6. - meijerint_definite
  7. - meijerint_inversion
  8. They can be used to compute, respectively, indefinite integrals, definite
  9. integrals over intervals of the real line, and inverse laplace-type integrals
  10. (from c-I*oo to c+I*oo). See the respective docstrings for details.
  11. The main references for this are:
  12. [L] Luke, Y. L. (1969), The Special Functions and Their Approximations,
  13. Volume 1
  14. [R] Kelly B. Roach. Meijer G Function Representations.
  15. In: Proceedings of the 1997 International Symposium on Symbolic and
  16. Algebraic Computation, pages 205-211, New York, 1997. ACM.
  17. [P] A. P. Prudnikov, Yu. A. Brychkov and O. I. Marichev (1990).
  18. Integrals and Series: More Special Functions, Vol. 3,.
  19. Gordon and Breach Science Publisher
  20. """
  21. from typing import Dict as tDict, Tuple as tTuple
  22. from sympy import SYMPY_DEBUG
  23. from sympy.core import S, Expr
  24. from sympy.core.add import Add
  25. from sympy.core.cache import cacheit
  26. from sympy.core.containers import Tuple
  27. from sympy.core.exprtools import factor_terms
  28. from sympy.core.function import (expand, expand_mul, expand_power_base,
  29. expand_trig, Function)
  30. from sympy.core.mul import Mul
  31. from sympy.core.numbers import ilcm, Rational, pi
  32. from sympy.core.relational import Eq, Ne, _canonical_coeff
  33. from sympy.core.sorting import default_sort_key, ordered
  34. from sympy.core.symbol import Dummy, symbols, Wild
  35. from sympy.functions.combinatorial.factorials import factorial
  36. from sympy.functions.elementary.complexes import (re, im, arg, Abs, sign,
  37. unpolarify, polarify, polar_lift, principal_branch, unbranched_argument,
  38. periodic_argument)
  39. from sympy.functions.elementary.exponential import exp, exp_polar, log
  40. from sympy.functions.elementary.integers import ceiling
  41. from sympy.functions.elementary.hyperbolic import (cosh, sinh,
  42. _rewrite_hyperbolics_as_exp, HyperbolicFunction)
  43. from sympy.functions.elementary.miscellaneous import sqrt
  44. from sympy.functions.elementary.piecewise import Piecewise, piecewise_fold
  45. from sympy.functions.elementary.trigonometric import (cos, sin, sinc,
  46. TrigonometricFunction)
  47. from sympy.functions.special.bessel import besselj, bessely, besseli, besselk
  48. from sympy.functions.special.delta_functions import DiracDelta, Heaviside
  49. from sympy.functions.special.elliptic_integrals import elliptic_k, elliptic_e
  50. from sympy.functions.special.error_functions import (erf, erfc, erfi, Ei,
  51. expint, Si, Ci, Shi, Chi, fresnels, fresnelc)
  52. from sympy.functions.special.gamma_functions import gamma
  53. from sympy.functions.special.hyper import hyper, meijerg
  54. from sympy.functions.special.singularity_functions import SingularityFunction
  55. from .integrals import Integral
  56. from sympy.logic.boolalg import And, Or, BooleanAtom, Not, BooleanFunction
  57. from sympy.polys import cancel, factor
  58. from sympy.simplify.fu import sincos_to_sum
  59. from sympy.simplify import (collect, gammasimp, hyperexpand, powdenest,
  60. powsimp, simplify)
  61. from sympy.utilities.iterables import multiset_partitions
  62. from sympy.utilities.misc import debug as _debug
  63. # keep this at top for easy reference
  64. z = Dummy('z')
  65. def _has(res, *f):
  66. # return True if res has f; in the case of Piecewise
  67. # only return True if *all* pieces have f
  68. res = piecewise_fold(res)
  69. if getattr(res, 'is_Piecewise', False):
  70. return all(_has(i, *f) for i in res.args)
  71. return res.has(*f)
  72. def _create_lookup_table(table):
  73. """ Add formulae for the function -> meijerg lookup table. """
  74. def wild(n):
  75. return Wild(n, exclude=[z])
  76. p, q, a, b, c = list(map(wild, 'pqabc'))
  77. n = Wild('n', properties=[lambda x: x.is_Integer and x > 0])
  78. t = p*z**q
  79. def add(formula, an, ap, bm, bq, arg=t, fac=S.One, cond=True, hint=True):
  80. table.setdefault(_mytype(formula, z), []).append((formula,
  81. [(fac, meijerg(an, ap, bm, bq, arg))], cond, hint))
  82. def addi(formula, inst, cond, hint=True):
  83. table.setdefault(
  84. _mytype(formula, z), []).append((formula, inst, cond, hint))
  85. def constant(a):
  86. return [(a, meijerg([1], [], [], [0], z)),
  87. (a, meijerg([], [1], [0], [], z))]
  88. table[()] = [(a, constant(a), True, True)]
  89. # [P], Section 8.
  90. class IsNonPositiveInteger(Function):
  91. @classmethod
  92. def eval(cls, arg):
  93. arg = unpolarify(arg)
  94. if arg.is_Integer is True:
  95. return arg <= 0
  96. # Section 8.4.2
  97. # TODO this needs more polar_lift (c/f entry for exp)
  98. add(Heaviside(t - b)*(t - b)**(a - 1), [a], [], [], [0], t/b,
  99. gamma(a)*b**(a - 1), And(b > 0))
  100. add(Heaviside(b - t)*(b - t)**(a - 1), [], [a], [0], [], t/b,
  101. gamma(a)*b**(a - 1), And(b > 0))
  102. add(Heaviside(z - (b/p)**(1/q))*(t - b)**(a - 1), [a], [], [], [0], t/b,
  103. gamma(a)*b**(a - 1), And(b > 0))
  104. add(Heaviside((b/p)**(1/q) - z)*(b - t)**(a - 1), [], [a], [0], [], t/b,
  105. gamma(a)*b**(a - 1), And(b > 0))
  106. add((b + t)**(-a), [1 - a], [], [0], [], t/b, b**(-a)/gamma(a),
  107. hint=Not(IsNonPositiveInteger(a)))
  108. add(Abs(b - t)**(-a), [1 - a], [(1 - a)/2], [0], [(1 - a)/2], t/b,
  109. 2*sin(pi*a/2)*gamma(1 - a)*Abs(b)**(-a), re(a) < 1)
  110. add((t**a - b**a)/(t - b), [0, a], [], [0, a], [], t/b,
  111. b**(a - 1)*sin(a*pi)/pi)
  112. # 12
  113. def A1(r, sign, nu):
  114. return pi**Rational(-1, 2)*(-sign*nu/2)**(1 - 2*r)
  115. def tmpadd(r, sgn):
  116. # XXX the a**2 is bad for matching
  117. add((sqrt(a**2 + t) + sgn*a)**b/(a**2 + t)**r,
  118. [(1 + b)/2, 1 - 2*r + b/2], [],
  119. [(b - sgn*b)/2], [(b + sgn*b)/2], t/a**2,
  120. a**(b - 2*r)*A1(r, sgn, b))
  121. tmpadd(0, 1)
  122. tmpadd(0, -1)
  123. tmpadd(S.Half, 1)
  124. tmpadd(S.Half, -1)
  125. # 13
  126. def tmpadd(r, sgn):
  127. add((sqrt(a + p*z**q) + sgn*sqrt(p)*z**(q/2))**b/(a + p*z**q)**r,
  128. [1 - r + sgn*b/2], [1 - r - sgn*b/2], [0, S.Half], [],
  129. p*z**q/a, a**(b/2 - r)*A1(r, sgn, b))
  130. tmpadd(0, 1)
  131. tmpadd(0, -1)
  132. tmpadd(S.Half, 1)
  133. tmpadd(S.Half, -1)
  134. # (those after look obscure)
  135. # Section 8.4.3
  136. add(exp(polar_lift(-1)*t), [], [], [0], [])
  137. # TODO can do sin^n, sinh^n by expansion ... where?
  138. # 8.4.4 (hyperbolic functions)
  139. add(sinh(t), [], [1], [S.Half], [1, 0], t**2/4, pi**Rational(3, 2))
  140. add(cosh(t), [], [S.Half], [0], [S.Half, S.Half], t**2/4, pi**Rational(3, 2))
  141. # Section 8.4.5
  142. # TODO can do t + a. but can also do by expansion... (XXX not really)
  143. add(sin(t), [], [], [S.Half], [0], t**2/4, sqrt(pi))
  144. add(cos(t), [], [], [0], [S.Half], t**2/4, sqrt(pi))
  145. # Section 8.4.6 (sinc function)
  146. add(sinc(t), [], [], [0], [Rational(-1, 2)], t**2/4, sqrt(pi)/2)
  147. # Section 8.5.5
  148. def make_log1(subs):
  149. N = subs[n]
  150. return [(S.NegativeOne**N*factorial(N),
  151. meijerg([], [1]*(N + 1), [0]*(N + 1), [], t))]
  152. def make_log2(subs):
  153. N = subs[n]
  154. return [(factorial(N),
  155. meijerg([1]*(N + 1), [], [], [0]*(N + 1), t))]
  156. # TODO these only hold for positive p, and can be made more general
  157. # but who uses log(x)*Heaviside(a-x) anyway ...
  158. # TODO also it would be nice to derive them recursively ...
  159. addi(log(t)**n*Heaviside(1 - t), make_log1, True)
  160. addi(log(t)**n*Heaviside(t - 1), make_log2, True)
  161. def make_log3(subs):
  162. return make_log1(subs) + make_log2(subs)
  163. addi(log(t)**n, make_log3, True)
  164. addi(log(t + a),
  165. constant(log(a)) + [(S.One, meijerg([1, 1], [], [1], [0], t/a))],
  166. True)
  167. addi(log(Abs(t - a)), constant(log(Abs(a))) +
  168. [(pi, meijerg([1, 1], [S.Half], [1], [0, S.Half], t/a))],
  169. True)
  170. # TODO log(x)/(x+a) and log(x)/(x-1) can also be done. should they
  171. # be derivable?
  172. # TODO further formulae in this section seem obscure
  173. # Sections 8.4.9-10
  174. # TODO
  175. # Section 8.4.11
  176. addi(Ei(t),
  177. constant(-S.ImaginaryUnit*pi) + [(S.NegativeOne, meijerg([], [1], [0, 0], [],
  178. t*polar_lift(-1)))],
  179. True)
  180. # Section 8.4.12
  181. add(Si(t), [1], [], [S.Half], [0, 0], t**2/4, sqrt(pi)/2)
  182. add(Ci(t), [], [1], [0, 0], [S.Half], t**2/4, -sqrt(pi)/2)
  183. # Section 8.4.13
  184. add(Shi(t), [S.Half], [], [0], [Rational(-1, 2), Rational(-1, 2)], polar_lift(-1)*t**2/4,
  185. t*sqrt(pi)/4)
  186. add(Chi(t), [], [S.Half, 1], [0, 0], [S.Half, S.Half], t**2/4, -
  187. pi**S('3/2')/2)
  188. # generalized exponential integral
  189. add(expint(a, t), [], [a], [a - 1, 0], [], t)
  190. # Section 8.4.14
  191. add(erf(t), [1], [], [S.Half], [0], t**2, 1/sqrt(pi))
  192. # TODO exp(-x)*erf(I*x) does not work
  193. add(erfc(t), [], [1], [0, S.Half], [], t**2, 1/sqrt(pi))
  194. # This formula for erfi(z) yields a wrong(?) minus sign
  195. #add(erfi(t), [1], [], [S.Half], [0], -t**2, I/sqrt(pi))
  196. add(erfi(t), [S.Half], [], [0], [Rational(-1, 2)], -t**2, t/sqrt(pi))
  197. # Fresnel Integrals
  198. add(fresnels(t), [1], [], [Rational(3, 4)], [0, Rational(1, 4)], pi**2*t**4/16, S.Half)
  199. add(fresnelc(t), [1], [], [Rational(1, 4)], [0, Rational(3, 4)], pi**2*t**4/16, S.Half)
  200. ##### bessel-type functions #####
  201. # Section 8.4.19
  202. add(besselj(a, t), [], [], [a/2], [-a/2], t**2/4)
  203. # all of the following are derivable
  204. #add(sin(t)*besselj(a, t), [Rational(1, 4), Rational(3, 4)], [], [(1+a)/2],
  205. # [-a/2, a/2, (1-a)/2], t**2, 1/sqrt(2))
  206. #add(cos(t)*besselj(a, t), [Rational(1, 4), Rational(3, 4)], [], [a/2],
  207. # [-a/2, (1+a)/2, (1-a)/2], t**2, 1/sqrt(2))
  208. #add(besselj(a, t)**2, [S.Half], [], [a], [-a, 0], t**2, 1/sqrt(pi))
  209. #add(besselj(a, t)*besselj(b, t), [0, S.Half], [], [(a + b)/2],
  210. # [-(a+b)/2, (a - b)/2, (b - a)/2], t**2, 1/sqrt(pi))
  211. # Section 8.4.20
  212. add(bessely(a, t), [], [-(a + 1)/2], [a/2, -a/2], [-(a + 1)/2], t**2/4)
  213. # TODO all of the following should be derivable
  214. #add(sin(t)*bessely(a, t), [Rational(1, 4), Rational(3, 4)], [(1 - a - 1)/2],
  215. # [(1 + a)/2, (1 - a)/2], [(1 - a - 1)/2, (1 - 1 - a)/2, (1 - 1 + a)/2],
  216. # t**2, 1/sqrt(2))
  217. #add(cos(t)*bessely(a, t), [Rational(1, 4), Rational(3, 4)], [(0 - a - 1)/2],
  218. # [(0 + a)/2, (0 - a)/2], [(0 - a - 1)/2, (1 - 0 - a)/2, (1 - 0 + a)/2],
  219. # t**2, 1/sqrt(2))
  220. #add(besselj(a, t)*bessely(b, t), [0, S.Half], [(a - b - 1)/2],
  221. # [(a + b)/2, (a - b)/2], [(a - b - 1)/2, -(a + b)/2, (b - a)/2],
  222. # t**2, 1/sqrt(pi))
  223. #addi(bessely(a, t)**2,
  224. # [(2/sqrt(pi), meijerg([], [S.Half, S.Half - a], [0, a, -a],
  225. # [S.Half - a], t**2)),
  226. # (1/sqrt(pi), meijerg([S.Half], [], [a], [-a, 0], t**2))],
  227. # True)
  228. #addi(bessely(a, t)*bessely(b, t),
  229. # [(2/sqrt(pi), meijerg([], [0, S.Half, (1 - a - b)/2],
  230. # [(a + b)/2, (a - b)/2, (b - a)/2, -(a + b)/2],
  231. # [(1 - a - b)/2], t**2)),
  232. # (1/sqrt(pi), meijerg([0, S.Half], [], [(a + b)/2],
  233. # [-(a + b)/2, (a - b)/2, (b - a)/2], t**2))],
  234. # True)
  235. # Section 8.4.21 ?
  236. # Section 8.4.22
  237. add(besseli(a, t), [], [(1 + a)/2], [a/2], [-a/2, (1 + a)/2], t**2/4, pi)
  238. # TODO many more formulas. should all be derivable
  239. # Section 8.4.23
  240. add(besselk(a, t), [], [], [a/2, -a/2], [], t**2/4, S.Half)
  241. # TODO many more formulas. should all be derivable
  242. # Complete elliptic integrals K(z) and E(z)
  243. add(elliptic_k(t), [S.Half, S.Half], [], [0], [0], -t, S.Half)
  244. add(elliptic_e(t), [S.Half, 3*S.Half], [], [0], [0], -t, Rational(-1, 2)/2)
  245. ####################################################################
  246. # First some helper functions.
  247. ####################################################################
  248. from sympy.utilities.timeutils import timethis
  249. timeit = timethis('meijerg')
  250. def _mytype(f, x):
  251. """ Create a hashable entity describing the type of f. """
  252. if x not in f.free_symbols:
  253. return ()
  254. elif f.is_Function:
  255. return (type(f),)
  256. else:
  257. types = [_mytype(a, x) for a in f.args]
  258. res = []
  259. for t in types:
  260. res += list(t)
  261. res.sort()
  262. return tuple(res)
  263. class _CoeffExpValueError(ValueError):
  264. """
  265. Exception raised by _get_coeff_exp, for internal use only.
  266. """
  267. pass
  268. def _get_coeff_exp(expr, x):
  269. """
  270. When expr is known to be of the form c*x**b, with c and/or b possibly 1,
  271. return c, b.
  272. Examples
  273. ========
  274. >>> from sympy.abc import x, a, b
  275. >>> from sympy.integrals.meijerint import _get_coeff_exp
  276. >>> _get_coeff_exp(a*x**b, x)
  277. (a, b)
  278. >>> _get_coeff_exp(x, x)
  279. (1, 1)
  280. >>> _get_coeff_exp(2*x, x)
  281. (2, 1)
  282. >>> _get_coeff_exp(x**3, x)
  283. (1, 3)
  284. """
  285. (c, m) = expand_power_base(powsimp(expr)).as_coeff_mul(x)
  286. if not m:
  287. return c, S.Zero
  288. [m] = m
  289. if m.is_Pow:
  290. if m.base != x:
  291. raise _CoeffExpValueError('expr not of form a*x**b')
  292. return c, m.exp
  293. elif m == x:
  294. return c, S.One
  295. else:
  296. raise _CoeffExpValueError('expr not of form a*x**b: %s' % expr)
  297. def _exponents(expr, x):
  298. """
  299. Find the exponents of ``x`` (not including zero) in ``expr``.
  300. Examples
  301. ========
  302. >>> from sympy.integrals.meijerint import _exponents
  303. >>> from sympy.abc import x, y
  304. >>> from sympy import sin
  305. >>> _exponents(x, x)
  306. {1}
  307. >>> _exponents(x**2, x)
  308. {2}
  309. >>> _exponents(x**2 + x, x)
  310. {1, 2}
  311. >>> _exponents(x**3*sin(x + x**y) + 1/x, x)
  312. {-1, 1, 3, y}
  313. """
  314. def _exponents_(expr, x, res):
  315. if expr == x:
  316. res.update([1])
  317. return
  318. if expr.is_Pow and expr.base == x:
  319. res.update([expr.exp])
  320. return
  321. for argument in expr.args:
  322. _exponents_(argument, x, res)
  323. res = set()
  324. _exponents_(expr, x, res)
  325. return res
  326. def _functions(expr, x):
  327. """ Find the types of functions in expr, to estimate the complexity. """
  328. return {e.func for e in expr.atoms(Function) if x in e.free_symbols}
  329. def _find_splitting_points(expr, x):
  330. """
  331. Find numbers a such that a linear substitution x -> x + a would
  332. (hopefully) simplify expr.
  333. Examples
  334. ========
  335. >>> from sympy.integrals.meijerint import _find_splitting_points as fsp
  336. >>> from sympy import sin
  337. >>> from sympy.abc import x
  338. >>> fsp(x, x)
  339. {0}
  340. >>> fsp((x-1)**3, x)
  341. {1}
  342. >>> fsp(sin(x+3)*x, x)
  343. {-3, 0}
  344. """
  345. p, q = [Wild(n, exclude=[x]) for n in 'pq']
  346. def compute_innermost(expr, res):
  347. if not isinstance(expr, Expr):
  348. return
  349. m = expr.match(p*x + q)
  350. if m and m[p] != 0:
  351. res.add(-m[q]/m[p])
  352. return
  353. if expr.is_Atom:
  354. return
  355. for argument in expr.args:
  356. compute_innermost(argument, res)
  357. innermost = set()
  358. compute_innermost(expr, innermost)
  359. return innermost
  360. def _split_mul(f, x):
  361. """
  362. Split expression ``f`` into fac, po, g, where fac is a constant factor,
  363. po = x**s for some s independent of s, and g is "the rest".
  364. Examples
  365. ========
  366. >>> from sympy.integrals.meijerint import _split_mul
  367. >>> from sympy import sin
  368. >>> from sympy.abc import s, x
  369. >>> _split_mul((3*x)**s*sin(x**2)*x, x)
  370. (3**s, x*x**s, sin(x**2))
  371. """
  372. fac = S.One
  373. po = S.One
  374. g = S.One
  375. f = expand_power_base(f)
  376. args = Mul.make_args(f)
  377. for a in args:
  378. if a == x:
  379. po *= x
  380. elif x not in a.free_symbols:
  381. fac *= a
  382. else:
  383. if a.is_Pow and x not in a.exp.free_symbols:
  384. c, t = a.base.as_coeff_mul(x)
  385. if t != (x,):
  386. c, t = expand_mul(a.base).as_coeff_mul(x)
  387. if t == (x,):
  388. po *= x**a.exp
  389. fac *= unpolarify(polarify(c**a.exp, subs=False))
  390. continue
  391. g *= a
  392. return fac, po, g
  393. def _mul_args(f):
  394. """
  395. Return a list ``L`` such that ``Mul(*L) == f``.
  396. If ``f`` is not a ``Mul`` or ``Pow``, ``L=[f]``.
  397. If ``f=g**n`` for an integer ``n``, ``L=[g]*n``.
  398. If ``f`` is a ``Mul``, ``L`` comes from applying ``_mul_args`` to all factors of ``f``.
  399. """
  400. args = Mul.make_args(f)
  401. gs = []
  402. for g in args:
  403. if g.is_Pow and g.exp.is_Integer:
  404. n = g.exp
  405. base = g.base
  406. if n < 0:
  407. n = -n
  408. base = 1/base
  409. gs += [base]*n
  410. else:
  411. gs.append(g)
  412. return gs
  413. def _mul_as_two_parts(f):
  414. """
  415. Find all the ways to split ``f`` into a product of two terms.
  416. Return None on failure.
  417. Explanation
  418. ===========
  419. Although the order is canonical from multiset_partitions, this is
  420. not necessarily the best order to process the terms. For example,
  421. if the case of len(gs) == 2 is removed and multiset is allowed to
  422. sort the terms, some tests fail.
  423. Examples
  424. ========
  425. >>> from sympy.integrals.meijerint import _mul_as_two_parts
  426. >>> from sympy import sin, exp, ordered
  427. >>> from sympy.abc import x
  428. >>> list(ordered(_mul_as_two_parts(x*sin(x)*exp(x))))
  429. [(x, exp(x)*sin(x)), (x*exp(x), sin(x)), (x*sin(x), exp(x))]
  430. """
  431. gs = _mul_args(f)
  432. if len(gs) < 2:
  433. return None
  434. if len(gs) == 2:
  435. return [tuple(gs)]
  436. return [(Mul(*x), Mul(*y)) for (x, y) in multiset_partitions(gs, 2)]
  437. def _inflate_g(g, n):
  438. """ Return C, h such that h is a G function of argument z**n and
  439. g = C*h. """
  440. # TODO should this be a method of meijerg?
  441. # See: [L, page 150, equation (5)]
  442. def inflate(params, n):
  443. """ (a1, .., ak) -> (a1/n, (a1+1)/n, ..., (ak + n-1)/n) """
  444. res = []
  445. for a in params:
  446. for i in range(n):
  447. res.append((a + i)/n)
  448. return res
  449. v = S(len(g.ap) - len(g.bq))
  450. C = n**(1 + g.nu + v/2)
  451. C /= (2*pi)**((n - 1)*g.delta)
  452. return C, meijerg(inflate(g.an, n), inflate(g.aother, n),
  453. inflate(g.bm, n), inflate(g.bother, n),
  454. g.argument**n * n**(n*v))
  455. def _flip_g(g):
  456. """ Turn the G function into one of inverse argument
  457. (i.e. G(1/x) -> G'(x)) """
  458. # See [L], section 5.2
  459. def tr(l):
  460. return [1 - a for a in l]
  461. return meijerg(tr(g.bm), tr(g.bother), tr(g.an), tr(g.aother), 1/g.argument)
  462. def _inflate_fox_h(g, a):
  463. r"""
  464. Let d denote the integrand in the definition of the G function ``g``.
  465. Consider the function H which is defined in the same way, but with
  466. integrand d/Gamma(a*s) (contour conventions as usual).
  467. If ``a`` is rational, the function H can be written as C*G, for a constant C
  468. and a G-function G.
  469. This function returns C, G.
  470. """
  471. if a < 0:
  472. return _inflate_fox_h(_flip_g(g), -a)
  473. p = S(a.p)
  474. q = S(a.q)
  475. # We use the substitution s->qs, i.e. inflate g by q. We are left with an
  476. # extra factor of Gamma(p*s), for which we use Gauss' multiplication
  477. # theorem.
  478. D, g = _inflate_g(g, q)
  479. z = g.argument
  480. D /= (2*pi)**((1 - p)/2)*p**Rational(-1, 2)
  481. z /= p**p
  482. bs = [(n + 1)/p for n in range(p)]
  483. return D, meijerg(g.an, g.aother, g.bm, list(g.bother) + bs, z)
  484. _dummies = {} # type: tDict[tTuple[str, str], Dummy]
  485. def _dummy(name, token, expr, **kwargs):
  486. """
  487. Return a dummy. This will return the same dummy if the same token+name is
  488. requested more than once, and it is not already in expr.
  489. This is for being cache-friendly.
  490. """
  491. d = _dummy_(name, token, **kwargs)
  492. if d in expr.free_symbols:
  493. return Dummy(name, **kwargs)
  494. return d
  495. def _dummy_(name, token, **kwargs):
  496. """
  497. Return a dummy associated to name and token. Same effect as declaring
  498. it globally.
  499. """
  500. global _dummies
  501. if not (name, token) in _dummies:
  502. _dummies[(name, token)] = Dummy(name, **kwargs)
  503. return _dummies[(name, token)]
  504. def _is_analytic(f, x):
  505. """ Check if f(x), when expressed using G functions on the positive reals,
  506. will in fact agree with the G functions almost everywhere """
  507. return not any(x in expr.free_symbols for expr in f.atoms(Heaviside, Abs))
  508. def _condsimp(cond, first=True):
  509. """
  510. Do naive simplifications on ``cond``.
  511. Explanation
  512. ===========
  513. Note that this routine is completely ad-hoc, simplification rules being
  514. added as need arises rather than following any logical pattern.
  515. Examples
  516. ========
  517. >>> from sympy.integrals.meijerint import _condsimp as simp
  518. >>> from sympy import Or, Eq
  519. >>> from sympy.abc import x, y
  520. >>> simp(Or(x < y, Eq(x, y)))
  521. x <= y
  522. """
  523. if first:
  524. cond = cond.replace(lambda _: _.is_Relational, _canonical_coeff)
  525. first = False
  526. if not isinstance(cond, BooleanFunction):
  527. return cond
  528. p, q, r = symbols('p q r', cls=Wild)
  529. # transforms tests use 0, 4, 5 and 11-14
  530. # meijer tests use 0, 2, 11, 14
  531. # joint_rv uses 6, 7
  532. rules = [
  533. (Or(p < q, Eq(p, q)), p <= q), # 0
  534. # The next two obviously are instances of a general pattern, but it is
  535. # easier to spell out the few cases we care about.
  536. (And(Abs(arg(p)) <= pi, Abs(arg(p) - 2*pi) <= pi),
  537. Eq(arg(p) - pi, 0)), # 1
  538. (And(Abs(2*arg(p) + pi) <= pi, Abs(2*arg(p) - pi) <= pi),
  539. Eq(arg(p), 0)), # 2
  540. (And(Abs(2*arg(p) + pi) < pi, Abs(2*arg(p) - pi) <= pi),
  541. S.false), # 3
  542. (And(Abs(arg(p) - pi/2) <= pi/2, Abs(arg(p) + pi/2) <= pi/2),
  543. Eq(arg(p), 0)), # 4
  544. (And(Abs(arg(p) - pi/2) <= pi/2, Abs(arg(p) + pi/2) < pi/2),
  545. S.false), # 5
  546. (And(Abs(arg(p**2/2 + 1)) < pi, Ne(Abs(arg(p**2/2 + 1)), pi)),
  547. S.true), # 6
  548. (Or(Abs(arg(p**2/2 + 1)) < pi, Ne(1/(p**2/2 + 1), 0)),
  549. S.true), # 7
  550. (And(Abs(unbranched_argument(p)) <= pi,
  551. Abs(unbranched_argument(exp_polar(-2*pi*S.ImaginaryUnit)*p)) <= pi),
  552. Eq(unbranched_argument(exp_polar(-S.ImaginaryUnit*pi)*p), 0)), # 8
  553. (And(Abs(unbranched_argument(p)) <= pi/2,
  554. Abs(unbranched_argument(exp_polar(-pi*S.ImaginaryUnit)*p)) <= pi/2),
  555. Eq(unbranched_argument(exp_polar(-S.ImaginaryUnit*pi/2)*p), 0)), # 9
  556. (Or(p <= q, And(p < q, r)), p <= q), # 10
  557. (Ne(p**2, 1) & (p**2 > 1), p**2 > 1), # 11
  558. (Ne(1/p, 1) & (cos(Abs(arg(p)))*Abs(p) > 1), Abs(p) > 1), # 12
  559. (Ne(p, 2) & (cos(Abs(arg(p)))*Abs(p) > 2), Abs(p) > 2), # 13
  560. ((Abs(arg(p)) < pi/2) & (cos(Abs(arg(p)))*sqrt(Abs(p**2)) > 1), p**2 > 1), # 14
  561. ]
  562. cond = cond.func(*list(map(lambda _: _condsimp(_, first), cond.args)))
  563. change = True
  564. while change:
  565. change = False
  566. for irule, (fro, to) in enumerate(rules):
  567. if fro.func != cond.func:
  568. continue
  569. for n, arg1 in enumerate(cond.args):
  570. if r in fro.args[0].free_symbols:
  571. m = arg1.match(fro.args[1])
  572. num = 1
  573. else:
  574. num = 0
  575. m = arg1.match(fro.args[0])
  576. if not m:
  577. continue
  578. otherargs = [x.subs(m) for x in fro.args[:num] + fro.args[num + 1:]]
  579. otherlist = [n]
  580. for arg2 in otherargs:
  581. for k, arg3 in enumerate(cond.args):
  582. if k in otherlist:
  583. continue
  584. if arg2 == arg3:
  585. otherlist += [k]
  586. break
  587. if isinstance(arg3, And) and arg2.args[1] == r and \
  588. isinstance(arg2, And) and arg2.args[0] in arg3.args:
  589. otherlist += [k]
  590. break
  591. if isinstance(arg3, And) and arg2.args[0] == r and \
  592. isinstance(arg2, And) and arg2.args[1] in arg3.args:
  593. otherlist += [k]
  594. break
  595. if len(otherlist) != len(otherargs) + 1:
  596. continue
  597. newargs = [arg_ for (k, arg_) in enumerate(cond.args)
  598. if k not in otherlist] + [to.subs(m)]
  599. if SYMPY_DEBUG:
  600. if irule not in (0, 2, 4, 5, 6, 7, 11, 12, 13, 14):
  601. print('used new rule:', irule)
  602. cond = cond.func(*newargs)
  603. change = True
  604. break
  605. # final tweak
  606. def rel_touchup(rel):
  607. if rel.rel_op != '==' or rel.rhs != 0:
  608. return rel
  609. # handle Eq(*, 0)
  610. LHS = rel.lhs
  611. m = LHS.match(arg(p)**q)
  612. if not m:
  613. m = LHS.match(unbranched_argument(polar_lift(p)**q))
  614. if not m:
  615. if isinstance(LHS, periodic_argument) and not LHS.args[0].is_polar \
  616. and LHS.args[1] is S.Infinity:
  617. return (LHS.args[0] > 0)
  618. return rel
  619. return (m[p] > 0)
  620. cond = cond.replace(lambda _: _.is_Relational, rel_touchup)
  621. if SYMPY_DEBUG:
  622. print('_condsimp: ', cond)
  623. return cond
  624. def _eval_cond(cond):
  625. """ Re-evaluate the conditions. """
  626. if isinstance(cond, bool):
  627. return cond
  628. return _condsimp(cond.doit())
  629. ####################################################################
  630. # Now the "backbone" functions to do actual integration.
  631. ####################################################################
  632. def _my_principal_branch(expr, period, full_pb=False):
  633. """ Bring expr nearer to its principal branch by removing superfluous
  634. factors.
  635. This function does *not* guarantee to yield the principal branch,
  636. to avoid introducing opaque principal_branch() objects,
  637. unless full_pb=True. """
  638. res = principal_branch(expr, period)
  639. if not full_pb:
  640. res = res.replace(principal_branch, lambda x, y: x)
  641. return res
  642. def _rewrite_saxena_1(fac, po, g, x):
  643. """
  644. Rewrite the integral fac*po*g dx, from zero to infinity, as
  645. integral fac*G, where G has argument a*x. Note po=x**s.
  646. Return fac, G.
  647. """
  648. _, s = _get_coeff_exp(po, x)
  649. a, b = _get_coeff_exp(g.argument, x)
  650. period = g.get_period()
  651. a = _my_principal_branch(a, period)
  652. # We substitute t = x**b.
  653. C = fac/(Abs(b)*a**((s + 1)/b - 1))
  654. # Absorb a factor of (at)**((1 + s)/b - 1).
  655. def tr(l):
  656. return [a + (1 + s)/b - 1 for a in l]
  657. return C, meijerg(tr(g.an), tr(g.aother), tr(g.bm), tr(g.bother),
  658. a*x)
  659. def _check_antecedents_1(g, x, helper=False):
  660. r"""
  661. Return a condition under which the mellin transform of g exists.
  662. Any power of x has already been absorbed into the G function,
  663. so this is just $\int_0^\infty g\, dx$.
  664. See [L, section 5.6.1]. (Note that s=1.)
  665. If ``helper`` is True, only check if the MT exists at infinity, i.e. if
  666. $\int_1^\infty g\, dx$ exists.
  667. """
  668. # NOTE if you update these conditions, please update the documentation as well
  669. delta = g.delta
  670. eta, _ = _get_coeff_exp(g.argument, x)
  671. m, n, p, q = S([len(g.bm), len(g.an), len(g.ap), len(g.bq)])
  672. if p > q:
  673. def tr(l):
  674. return [1 - x for x in l]
  675. return _check_antecedents_1(meijerg(tr(g.bm), tr(g.bother),
  676. tr(g.an), tr(g.aother), x/eta),
  677. x)
  678. tmp = []
  679. for b in g.bm:
  680. tmp += [-re(b) < 1]
  681. for a in g.an:
  682. tmp += [1 < 1 - re(a)]
  683. cond_3 = And(*tmp)
  684. for b in g.bother:
  685. tmp += [-re(b) < 1]
  686. for a in g.aother:
  687. tmp += [1 < 1 - re(a)]
  688. cond_3_star = And(*tmp)
  689. cond_4 = (-re(g.nu) + (q + 1 - p)/2 > q - p)
  690. def debug(*msg):
  691. _debug(*msg)
  692. debug('Checking antecedents for 1 function:')
  693. debug(' delta=%s, eta=%s, m=%s, n=%s, p=%s, q=%s'
  694. % (delta, eta, m, n, p, q))
  695. debug(' ap = %s, %s' % (list(g.an), list(g.aother)))
  696. debug(' bq = %s, %s' % (list(g.bm), list(g.bother)))
  697. debug(' cond_3=%s, cond_3*=%s, cond_4=%s' % (cond_3, cond_3_star, cond_4))
  698. conds = []
  699. # case 1
  700. case1 = []
  701. tmp1 = [1 <= n, p < q, 1 <= m]
  702. tmp2 = [1 <= p, 1 <= m, Eq(q, p + 1), Not(And(Eq(n, 0), Eq(m, p + 1)))]
  703. tmp3 = [1 <= p, Eq(q, p)]
  704. for k in range(ceiling(delta/2) + 1):
  705. tmp3 += [Ne(Abs(unbranched_argument(eta)), (delta - 2*k)*pi)]
  706. tmp = [delta > 0, Abs(unbranched_argument(eta)) < delta*pi]
  707. extra = [Ne(eta, 0), cond_3]
  708. if helper:
  709. extra = []
  710. for t in [tmp1, tmp2, tmp3]:
  711. case1 += [And(*(t + tmp + extra))]
  712. conds += case1
  713. debug(' case 1:', case1)
  714. # case 2
  715. extra = [cond_3]
  716. if helper:
  717. extra = []
  718. case2 = [And(Eq(n, 0), p + 1 <= m, m <= q,
  719. Abs(unbranched_argument(eta)) < delta*pi, *extra)]
  720. conds += case2
  721. debug(' case 2:', case2)
  722. # case 3
  723. extra = [cond_3, cond_4]
  724. if helper:
  725. extra = []
  726. case3 = [And(p < q, 1 <= m, delta > 0, Eq(Abs(unbranched_argument(eta)), delta*pi),
  727. *extra)]
  728. case3 += [And(p <= q - 2, Eq(delta, 0), Eq(Abs(unbranched_argument(eta)), 0), *extra)]
  729. conds += case3
  730. debug(' case 3:', case3)
  731. # TODO altered cases 4-7
  732. # extra case from wofram functions site:
  733. # (reproduced verbatim from Prudnikov, section 2.24.2)
  734. # http://functions.wolfram.com/HypergeometricFunctions/MeijerG/21/02/01/
  735. case_extra = []
  736. case_extra += [Eq(p, q), Eq(delta, 0), Eq(unbranched_argument(eta), 0), Ne(eta, 0)]
  737. if not helper:
  738. case_extra += [cond_3]
  739. s = []
  740. for a, b in zip(g.ap, g.bq):
  741. s += [b - a]
  742. case_extra += [re(Add(*s)) < 0]
  743. case_extra = And(*case_extra)
  744. conds += [case_extra]
  745. debug(' extra case:', [case_extra])
  746. case_extra_2 = [And(delta > 0, Abs(unbranched_argument(eta)) < delta*pi)]
  747. if not helper:
  748. case_extra_2 += [cond_3]
  749. case_extra_2 = And(*case_extra_2)
  750. conds += [case_extra_2]
  751. debug(' second extra case:', [case_extra_2])
  752. # TODO This leaves only one case from the three listed by Prudnikov.
  753. # Investigate if these indeed cover everything; if so, remove the rest.
  754. return Or(*conds)
  755. def _int0oo_1(g, x):
  756. r"""
  757. Evaluate $\int_0^\infty g\, dx$ using G functions,
  758. assuming the necessary conditions are fulfilled.
  759. Examples
  760. ========
  761. >>> from sympy.abc import a, b, c, d, x, y
  762. >>> from sympy import meijerg
  763. >>> from sympy.integrals.meijerint import _int0oo_1
  764. >>> _int0oo_1(meijerg([a], [b], [c], [d], x*y), x)
  765. gamma(-a)*gamma(c + 1)/(y*gamma(-d)*gamma(b + 1))
  766. """
  767. # See [L, section 5.6.1]. Note that s=1.
  768. eta, _ = _get_coeff_exp(g.argument, x)
  769. res = 1/eta
  770. # XXX TODO we should reduce order first
  771. for b in g.bm:
  772. res *= gamma(b + 1)
  773. for a in g.an:
  774. res *= gamma(1 - a - 1)
  775. for b in g.bother:
  776. res /= gamma(1 - b - 1)
  777. for a in g.aother:
  778. res /= gamma(a + 1)
  779. return gammasimp(unpolarify(res))
  780. def _rewrite_saxena(fac, po, g1, g2, x, full_pb=False):
  781. """
  782. Rewrite the integral ``fac*po*g1*g2`` from 0 to oo in terms of G
  783. functions with argument ``c*x``.
  784. Explanation
  785. ===========
  786. Return C, f1, f2 such that integral C f1 f2 from 0 to infinity equals
  787. integral fac ``po``, ``g1``, ``g2`` from 0 to infinity.
  788. Examples
  789. ========
  790. >>> from sympy.integrals.meijerint import _rewrite_saxena
  791. >>> from sympy.abc import s, t, m
  792. >>> from sympy import meijerg
  793. >>> g1 = meijerg([], [], [0], [], s*t)
  794. >>> g2 = meijerg([], [], [m/2], [-m/2], t**2/4)
  795. >>> r = _rewrite_saxena(1, t**0, g1, g2, t)
  796. >>> r[0]
  797. s/(4*sqrt(pi))
  798. >>> r[1]
  799. meijerg(((), ()), ((-1/2, 0), ()), s**2*t/4)
  800. >>> r[2]
  801. meijerg(((), ()), ((m/2,), (-m/2,)), t/4)
  802. """
  803. def pb(g):
  804. a, b = _get_coeff_exp(g.argument, x)
  805. per = g.get_period()
  806. return meijerg(g.an, g.aother, g.bm, g.bother,
  807. _my_principal_branch(a, per, full_pb)*x**b)
  808. _, s = _get_coeff_exp(po, x)
  809. _, b1 = _get_coeff_exp(g1.argument, x)
  810. _, b2 = _get_coeff_exp(g2.argument, x)
  811. if (b1 < 0) == True:
  812. b1 = -b1
  813. g1 = _flip_g(g1)
  814. if (b2 < 0) == True:
  815. b2 = -b2
  816. g2 = _flip_g(g2)
  817. if not b1.is_Rational or not b2.is_Rational:
  818. return
  819. m1, n1 = b1.p, b1.q
  820. m2, n2 = b2.p, b2.q
  821. tau = ilcm(m1*n2, m2*n1)
  822. r1 = tau//(m1*n2)
  823. r2 = tau//(m2*n1)
  824. C1, g1 = _inflate_g(g1, r1)
  825. C2, g2 = _inflate_g(g2, r2)
  826. g1 = pb(g1)
  827. g2 = pb(g2)
  828. fac *= C1*C2
  829. a1, b = _get_coeff_exp(g1.argument, x)
  830. a2, _ = _get_coeff_exp(g2.argument, x)
  831. # arbitrarily tack on the x**s part to g1
  832. # TODO should we try both?
  833. exp = (s + 1)/b - 1
  834. fac = fac/(Abs(b) * a1**exp)
  835. def tr(l):
  836. return [a + exp for a in l]
  837. g1 = meijerg(tr(g1.an), tr(g1.aother), tr(g1.bm), tr(g1.bother), a1*x)
  838. g2 = meijerg(g2.an, g2.aother, g2.bm, g2.bother, a2*x)
  839. return powdenest(fac, polar=True), g1, g2
  840. def _check_antecedents(g1, g2, x):
  841. """ Return a condition under which the integral theorem applies. """
  842. # Yes, this is madness.
  843. # XXX TODO this is a testing *nightmare*
  844. # NOTE if you update these conditions, please update the documentation as well
  845. # The following conditions are found in
  846. # [P], Section 2.24.1
  847. #
  848. # They are also reproduced (verbatim!) at
  849. # http://functions.wolfram.com/HypergeometricFunctions/MeijerG/21/02/03/
  850. #
  851. # Note: k=l=r=alpha=1
  852. sigma, _ = _get_coeff_exp(g1.argument, x)
  853. omega, _ = _get_coeff_exp(g2.argument, x)
  854. s, t, u, v = S([len(g1.bm), len(g1.an), len(g1.ap), len(g1.bq)])
  855. m, n, p, q = S([len(g2.bm), len(g2.an), len(g2.ap), len(g2.bq)])
  856. bstar = s + t - (u + v)/2
  857. cstar = m + n - (p + q)/2
  858. rho = g1.nu + (u - v)/2 + 1
  859. mu = g2.nu + (p - q)/2 + 1
  860. phi = q - p - (v - u)
  861. eta = 1 - (v - u) - mu - rho
  862. psi = (pi*(q - m - n) + Abs(unbranched_argument(omega)))/(q - p)
  863. theta = (pi*(v - s - t) + Abs(unbranched_argument(sigma)))/(v - u)
  864. _debug('Checking antecedents:')
  865. _debug(' sigma=%s, s=%s, t=%s, u=%s, v=%s, b*=%s, rho=%s'
  866. % (sigma, s, t, u, v, bstar, rho))
  867. _debug(' omega=%s, m=%s, n=%s, p=%s, q=%s, c*=%s, mu=%s,'
  868. % (omega, m, n, p, q, cstar, mu))
  869. _debug(' phi=%s, eta=%s, psi=%s, theta=%s' % (phi, eta, psi, theta))
  870. def _c1():
  871. for g in [g1, g2]:
  872. for i in g.an:
  873. for j in g.bm:
  874. diff = i - j
  875. if diff.is_integer and diff.is_positive:
  876. return False
  877. return True
  878. c1 = _c1()
  879. c2 = And(*[re(1 + i + j) > 0 for i in g1.bm for j in g2.bm])
  880. c3 = And(*[re(1 + i + j) < 1 + 1 for i in g1.an for j in g2.an])
  881. c4 = And(*[(p - q)*re(1 + i - 1) - re(mu) > Rational(-3, 2) for i in g1.an])
  882. c5 = And(*[(p - q)*re(1 + i) - re(mu) > Rational(-3, 2) for i in g1.bm])
  883. c6 = And(*[(u - v)*re(1 + i - 1) - re(rho) > Rational(-3, 2) for i in g2.an])
  884. c7 = And(*[(u - v)*re(1 + i) - re(rho) > Rational(-3, 2) for i in g2.bm])
  885. c8 = (Abs(phi) + 2*re((rho - 1)*(q - p) + (v - u)*(q - p) + (mu -
  886. 1)*(v - u)) > 0)
  887. c9 = (Abs(phi) - 2*re((rho - 1)*(q - p) + (v - u)*(q - p) + (mu -
  888. 1)*(v - u)) > 0)
  889. c10 = (Abs(unbranched_argument(sigma)) < bstar*pi)
  890. c11 = Eq(Abs(unbranched_argument(sigma)), bstar*pi)
  891. c12 = (Abs(unbranched_argument(omega)) < cstar*pi)
  892. c13 = Eq(Abs(unbranched_argument(omega)), cstar*pi)
  893. # The following condition is *not* implemented as stated on the wolfram
  894. # function site. In the book of Prudnikov there is an additional part
  895. # (the And involving re()). However, I only have this book in russian, and
  896. # I don't read any russian. The following condition is what other people
  897. # have told me it means.
  898. # Worryingly, it is different from the condition implemented in REDUCE.
  899. # The REDUCE implementation:
  900. # https://reduce-algebra.svn.sourceforge.net/svnroot/reduce-algebra/trunk/packages/defint/definta.red
  901. # (search for tst14)
  902. # The Wolfram alpha version:
  903. # http://functions.wolfram.com/HypergeometricFunctions/MeijerG/21/02/03/03/0014/
  904. z0 = exp(-(bstar + cstar)*pi*S.ImaginaryUnit)
  905. zos = unpolarify(z0*omega/sigma)
  906. zso = unpolarify(z0*sigma/omega)
  907. if zos == 1/zso:
  908. c14 = And(Eq(phi, 0), bstar + cstar <= 1,
  909. Or(Ne(zos, 1), re(mu + rho + v - u) < 1,
  910. re(mu + rho + q - p) < 1))
  911. else:
  912. def _cond(z):
  913. '''Returns True if abs(arg(1-z)) < pi, avoiding arg(0).
  914. Explanation
  915. ===========
  916. If ``z`` is 1 then arg is NaN. This raises a
  917. TypeError on `NaN < pi`. Previously this gave `False` so
  918. this behavior has been hardcoded here but someone should
  919. check if this NaN is more serious! This NaN is triggered by
  920. test_meijerint() in test_meijerint.py:
  921. `meijerint_definite(exp(x), x, 0, I)`
  922. '''
  923. return z != 1 and Abs(arg(1 - z)) < pi
  924. c14 = And(Eq(phi, 0), bstar - 1 + cstar <= 0,
  925. Or(And(Ne(zos, 1), _cond(zos)),
  926. And(re(mu + rho + v - u) < 1, Eq(zos, 1))))
  927. c14_alt = And(Eq(phi, 0), cstar - 1 + bstar <= 0,
  928. Or(And(Ne(zso, 1), _cond(zso)),
  929. And(re(mu + rho + q - p) < 1, Eq(zso, 1))))
  930. # Since r=k=l=1, in our case there is c14_alt which is the same as calling
  931. # us with (g1, g2) = (g2, g1). The conditions below enumerate all cases
  932. # (i.e. we don't have to try arguments reversed by hand), and indeed try
  933. # all symmetric cases. (i.e. whenever there is a condition involving c14,
  934. # there is also a dual condition which is exactly what we would get when g1,
  935. # g2 were interchanged, *but c14 was unaltered*).
  936. # Hence the following seems correct:
  937. c14 = Or(c14, c14_alt)
  938. '''
  939. When `c15` is NaN (e.g. from `psi` being NaN as happens during
  940. 'test_issue_4992' and/or `theta` is NaN as in 'test_issue_6253',
  941. both in `test_integrals.py`) the comparison to 0 formerly gave False
  942. whereas now an error is raised. To keep the old behavior, the value
  943. of NaN is replaced with False but perhaps a closer look at this condition
  944. should be made: XXX how should conditions leading to c15=NaN be handled?
  945. '''
  946. try:
  947. lambda_c = (q - p)*Abs(omega)**(1/(q - p))*cos(psi) \
  948. + (v - u)*Abs(sigma)**(1/(v - u))*cos(theta)
  949. # the TypeError might be raised here, e.g. if lambda_c is NaN
  950. if _eval_cond(lambda_c > 0) != False:
  951. c15 = (lambda_c > 0)
  952. else:
  953. def lambda_s0(c1, c2):
  954. return c1*(q - p)*Abs(omega)**(1/(q - p))*sin(psi) \
  955. + c2*(v - u)*Abs(sigma)**(1/(v - u))*sin(theta)
  956. lambda_s = Piecewise(
  957. ((lambda_s0(+1, +1)*lambda_s0(-1, -1)),
  958. And(Eq(unbranched_argument(sigma), 0), Eq(unbranched_argument(omega), 0))),
  959. (lambda_s0(sign(unbranched_argument(omega)), +1)*lambda_s0(sign(unbranched_argument(omega)), -1),
  960. And(Eq(unbranched_argument(sigma), 0), Ne(unbranched_argument(omega), 0))),
  961. (lambda_s0(+1, sign(unbranched_argument(sigma)))*lambda_s0(-1, sign(unbranched_argument(sigma))),
  962. And(Ne(unbranched_argument(sigma), 0), Eq(unbranched_argument(omega), 0))),
  963. (lambda_s0(sign(unbranched_argument(omega)), sign(unbranched_argument(sigma))), True))
  964. tmp = [lambda_c > 0,
  965. And(Eq(lambda_c, 0), Ne(lambda_s, 0), re(eta) > -1),
  966. And(Eq(lambda_c, 0), Eq(lambda_s, 0), re(eta) > 0)]
  967. c15 = Or(*tmp)
  968. except TypeError:
  969. c15 = False
  970. for cond, i in [(c1, 1), (c2, 2), (c3, 3), (c4, 4), (c5, 5), (c6, 6),
  971. (c7, 7), (c8, 8), (c9, 9), (c10, 10), (c11, 11),
  972. (c12, 12), (c13, 13), (c14, 14), (c15, 15)]:
  973. _debug(' c%s:' % i, cond)
  974. # We will return Or(*conds)
  975. conds = []
  976. def pr(count):
  977. _debug(' case %s:' % count, conds[-1])
  978. conds += [And(m*n*s*t != 0, bstar.is_positive is True, cstar.is_positive is True, c1, c2, c3, c10,
  979. c12)] # 1
  980. pr(1)
  981. conds += [And(Eq(u, v), Eq(bstar, 0), cstar.is_positive is True, sigma.is_positive is True, re(rho) < 1,
  982. c1, c2, c3, c12)] # 2
  983. pr(2)
  984. conds += [And(Eq(p, q), Eq(cstar, 0), bstar.is_positive is True, omega.is_positive is True, re(mu) < 1,
  985. c1, c2, c3, c10)] # 3
  986. pr(3)
  987. conds += [And(Eq(p, q), Eq(u, v), Eq(bstar, 0), Eq(cstar, 0),
  988. sigma.is_positive is True, omega.is_positive is True, re(mu) < 1, re(rho) < 1,
  989. Ne(sigma, omega), c1, c2, c3)] # 4
  990. pr(4)
  991. conds += [And(Eq(p, q), Eq(u, v), Eq(bstar, 0), Eq(cstar, 0),
  992. sigma.is_positive is True, omega.is_positive is True, re(mu + rho) < 1,
  993. Ne(omega, sigma), c1, c2, c3)] # 5
  994. pr(5)
  995. conds += [And(p > q, s.is_positive is True, bstar.is_positive is True, cstar >= 0,
  996. c1, c2, c3, c5, c10, c13)] # 6
  997. pr(6)
  998. conds += [And(p < q, t.is_positive is True, bstar.is_positive is True, cstar >= 0,
  999. c1, c2, c3, c4, c10, c13)] # 7
  1000. pr(7)
  1001. conds += [And(u > v, m.is_positive is True, cstar.is_positive is True, bstar >= 0,
  1002. c1, c2, c3, c7, c11, c12)] # 8
  1003. pr(8)
  1004. conds += [And(u < v, n.is_positive is True, cstar.is_positive is True, bstar >= 0,
  1005. c1, c2, c3, c6, c11, c12)] # 9
  1006. pr(9)
  1007. conds += [And(p > q, Eq(u, v), Eq(bstar, 0), cstar >= 0, sigma.is_positive is True,
  1008. re(rho) < 1, c1, c2, c3, c5, c13)] # 10
  1009. pr(10)
  1010. conds += [And(p < q, Eq(u, v), Eq(bstar, 0), cstar >= 0, sigma.is_positive is True,
  1011. re(rho) < 1, c1, c2, c3, c4, c13)] # 11
  1012. pr(11)
  1013. conds += [And(Eq(p, q), u > v, bstar >= 0, Eq(cstar, 0), omega.is_positive is True,
  1014. re(mu) < 1, c1, c2, c3, c7, c11)] # 12
  1015. pr(12)
  1016. conds += [And(Eq(p, q), u < v, bstar >= 0, Eq(cstar, 0), omega.is_positive is True,
  1017. re(mu) < 1, c1, c2, c3, c6, c11)] # 13
  1018. pr(13)
  1019. conds += [And(p < q, u > v, bstar >= 0, cstar >= 0,
  1020. c1, c2, c3, c4, c7, c11, c13)] # 14
  1021. pr(14)
  1022. conds += [And(p > q, u < v, bstar >= 0, cstar >= 0,
  1023. c1, c2, c3, c5, c6, c11, c13)] # 15
  1024. pr(15)
  1025. conds += [And(p > q, u > v, bstar >= 0, cstar >= 0,
  1026. c1, c2, c3, c5, c7, c8, c11, c13, c14)] # 16
  1027. pr(16)
  1028. conds += [And(p < q, u < v, bstar >= 0, cstar >= 0,
  1029. c1, c2, c3, c4, c6, c9, c11, c13, c14)] # 17
  1030. pr(17)
  1031. conds += [And(Eq(t, 0), s.is_positive is True, bstar.is_positive is True, phi.is_positive is True, c1, c2, c10)] # 18
  1032. pr(18)
  1033. conds += [And(Eq(s, 0), t.is_positive is True, bstar.is_positive is True, phi.is_negative is True, c1, c3, c10)] # 19
  1034. pr(19)
  1035. conds += [And(Eq(n, 0), m.is_positive is True, cstar.is_positive is True, phi.is_negative is True, c1, c2, c12)] # 20
  1036. pr(20)
  1037. conds += [And(Eq(m, 0), n.is_positive is True, cstar.is_positive is True, phi.is_positive is True, c1, c3, c12)] # 21
  1038. pr(21)
  1039. conds += [And(Eq(s*t, 0), bstar.is_positive is True, cstar.is_positive is True,
  1040. c1, c2, c3, c10, c12)] # 22
  1041. pr(22)
  1042. conds += [And(Eq(m*n, 0), bstar.is_positive is True, cstar.is_positive is True,
  1043. c1, c2, c3, c10, c12)] # 23
  1044. pr(23)
  1045. # The following case is from [Luke1969]. As far as I can tell, it is *not*
  1046. # covered by Prudnikov's.
  1047. # Let G1 and G2 be the two G-functions. Suppose the integral exists from
  1048. # 0 to a > 0 (this is easy the easy part), that G1 is exponential decay at
  1049. # infinity, and that the mellin transform of G2 exists.
  1050. # Then the integral exists.
  1051. mt1_exists = _check_antecedents_1(g1, x, helper=True)
  1052. mt2_exists = _check_antecedents_1(g2, x, helper=True)
  1053. conds += [And(mt2_exists, Eq(t, 0), u < s, bstar.is_positive is True, c10, c1, c2, c3)]
  1054. pr('E1')
  1055. conds += [And(mt2_exists, Eq(s, 0), v < t, bstar.is_positive is True, c10, c1, c2, c3)]
  1056. pr('E2')
  1057. conds += [And(mt1_exists, Eq(n, 0), p < m, cstar.is_positive is True, c12, c1, c2, c3)]
  1058. pr('E3')
  1059. conds += [And(mt1_exists, Eq(m, 0), q < n, cstar.is_positive is True, c12, c1, c2, c3)]
  1060. pr('E4')
  1061. # Let's short-circuit if this worked ...
  1062. # the rest is corner-cases and terrible to read.
  1063. r = Or(*conds)
  1064. if _eval_cond(r) != False:
  1065. return r
  1066. conds += [And(m + n > p, Eq(t, 0), Eq(phi, 0), s.is_positive is True, bstar.is_positive is True, cstar.is_negative is True,
  1067. Abs(unbranched_argument(omega)) < (m + n - p + 1)*pi,
  1068. c1, c2, c10, c14, c15)] # 24
  1069. pr(24)
  1070. conds += [And(m + n > q, Eq(s, 0), Eq(phi, 0), t.is_positive is True, bstar.is_positive is True, cstar.is_negative is True,
  1071. Abs(unbranched_argument(omega)) < (m + n - q + 1)*pi,
  1072. c1, c3, c10, c14, c15)] # 25
  1073. pr(25)
  1074. conds += [And(Eq(p, q - 1), Eq(t, 0), Eq(phi, 0), s.is_positive is True, bstar.is_positive is True,
  1075. cstar >= 0, cstar*pi < Abs(unbranched_argument(omega)),
  1076. c1, c2, c10, c14, c15)] # 26
  1077. pr(26)
  1078. conds += [And(Eq(p, q + 1), Eq(s, 0), Eq(phi, 0), t.is_positive is True, bstar.is_positive is True,
  1079. cstar >= 0, cstar*pi < Abs(unbranched_argument(omega)),
  1080. c1, c3, c10, c14, c15)] # 27
  1081. pr(27)
  1082. conds += [And(p < q - 1, Eq(t, 0), Eq(phi, 0), s.is_positive is True, bstar.is_positive is True,
  1083. cstar >= 0, cstar*pi < Abs(unbranched_argument(omega)),
  1084. Abs(unbranched_argument(omega)) < (m + n - p + 1)*pi,
  1085. c1, c2, c10, c14, c15)] # 28
  1086. pr(28)
  1087. conds += [And(
  1088. p > q + 1, Eq(s, 0), Eq(phi, 0), t.is_positive is True, bstar.is_positive is True, cstar >= 0,
  1089. cstar*pi < Abs(unbranched_argument(omega)),
  1090. Abs(unbranched_argument(omega)) < (m + n - q + 1)*pi,
  1091. c1, c3, c10, c14, c15)] # 29
  1092. pr(29)
  1093. conds += [And(Eq(n, 0), Eq(phi, 0), s + t > 0, m.is_positive is True, cstar.is_positive is True, bstar.is_negative is True,
  1094. Abs(unbranched_argument(sigma)) < (s + t - u + 1)*pi,
  1095. c1, c2, c12, c14, c15)] # 30
  1096. pr(30)
  1097. conds += [And(Eq(m, 0), Eq(phi, 0), s + t > v, n.is_positive is True, cstar.is_positive is True, bstar.is_negative is True,
  1098. Abs(unbranched_argument(sigma)) < (s + t - v + 1)*pi,
  1099. c1, c3, c12, c14, c15)] # 31
  1100. pr(31)
  1101. conds += [And(Eq(n, 0), Eq(phi, 0), Eq(u, v - 1), m.is_positive is True, cstar.is_positive is True,
  1102. bstar >= 0, bstar*pi < Abs(unbranched_argument(sigma)),
  1103. Abs(unbranched_argument(sigma)) < (bstar + 1)*pi,
  1104. c1, c2, c12, c14, c15)] # 32
  1105. pr(32)
  1106. conds += [And(Eq(m, 0), Eq(phi, 0), Eq(u, v + 1), n.is_positive is True, cstar.is_positive is True,
  1107. bstar >= 0, bstar*pi < Abs(unbranched_argument(sigma)),
  1108. Abs(unbranched_argument(sigma)) < (bstar + 1)*pi,
  1109. c1, c3, c12, c14, c15)] # 33
  1110. pr(33)
  1111. conds += [And(
  1112. Eq(n, 0), Eq(phi, 0), u < v - 1, m.is_positive is True, cstar.is_positive is True, bstar >= 0,
  1113. bstar*pi < Abs(unbranched_argument(sigma)),
  1114. Abs(unbranched_argument(sigma)) < (s + t - u + 1)*pi,
  1115. c1, c2, c12, c14, c15)] # 34
  1116. pr(34)
  1117. conds += [And(
  1118. Eq(m, 0), Eq(phi, 0), u > v + 1, n.is_positive is True, cstar.is_positive is True, bstar >= 0,
  1119. bstar*pi < Abs(unbranched_argument(sigma)),
  1120. Abs(unbranched_argument(sigma)) < (s + t - v + 1)*pi,
  1121. c1, c3, c12, c14, c15)] # 35
  1122. pr(35)
  1123. return Or(*conds)
  1124. # NOTE An alternative, but as far as I can tell weaker, set of conditions
  1125. # can be found in [L, section 5.6.2].
  1126. def _int0oo(g1, g2, x):
  1127. """
  1128. Express integral from zero to infinity g1*g2 using a G function,
  1129. assuming the necessary conditions are fulfilled.
  1130. Examples
  1131. ========
  1132. >>> from sympy.integrals.meijerint import _int0oo
  1133. >>> from sympy.abc import s, t, m
  1134. >>> from sympy import meijerg, S
  1135. >>> g1 = meijerg([], [], [-S(1)/2, 0], [], s**2*t/4)
  1136. >>> g2 = meijerg([], [], [m/2], [-m/2], t/4)
  1137. >>> _int0oo(g1, g2, t)
  1138. 4*meijerg(((1/2, 0), ()), ((m/2,), (-m/2,)), s**(-2))/s**2
  1139. """
  1140. # See: [L, section 5.6.2, equation (1)]
  1141. eta, _ = _get_coeff_exp(g1.argument, x)
  1142. omega, _ = _get_coeff_exp(g2.argument, x)
  1143. def neg(l):
  1144. return [-x for x in l]
  1145. a1 = neg(g1.bm) + list(g2.an)
  1146. a2 = list(g2.aother) + neg(g1.bother)
  1147. b1 = neg(g1.an) + list(g2.bm)
  1148. b2 = list(g2.bother) + neg(g1.aother)
  1149. return meijerg(a1, a2, b1, b2, omega/eta)/eta
  1150. def _rewrite_inversion(fac, po, g, x):
  1151. """ Absorb ``po`` == x**s into g. """
  1152. _, s = _get_coeff_exp(po, x)
  1153. a, b = _get_coeff_exp(g.argument, x)
  1154. def tr(l):
  1155. return [t + s/b for t in l]
  1156. return (powdenest(fac/a**(s/b), polar=True),
  1157. meijerg(tr(g.an), tr(g.aother), tr(g.bm), tr(g.bother), g.argument))
  1158. def _check_antecedents_inversion(g, x):
  1159. """ Check antecedents for the laplace inversion integral. """
  1160. _debug('Checking antecedents for inversion:')
  1161. z = g.argument
  1162. _, e = _get_coeff_exp(z, x)
  1163. if e < 0:
  1164. _debug(' Flipping G.')
  1165. # We want to assume that argument gets large as |x| -> oo
  1166. return _check_antecedents_inversion(_flip_g(g), x)
  1167. def statement_half(a, b, c, z, plus):
  1168. coeff, exponent = _get_coeff_exp(z, x)
  1169. a *= exponent
  1170. b *= coeff**c
  1171. c *= exponent
  1172. conds = []
  1173. wp = b*exp(S.ImaginaryUnit*re(c)*pi/2)
  1174. wm = b*exp(-S.ImaginaryUnit*re(c)*pi/2)
  1175. if plus:
  1176. w = wp
  1177. else:
  1178. w = wm
  1179. conds += [And(Or(Eq(b, 0), re(c) <= 0), re(a) <= -1)]
  1180. conds += [And(Ne(b, 0), Eq(im(c), 0), re(c) > 0, re(w) < 0)]
  1181. conds += [And(Ne(b, 0), Eq(im(c), 0), re(c) > 0, re(w) <= 0,
  1182. re(a) <= -1)]
  1183. return Or(*conds)
  1184. def statement(a, b, c, z):
  1185. """ Provide a convergence statement for z**a * exp(b*z**c),
  1186. c/f sphinx docs. """
  1187. return And(statement_half(a, b, c, z, True),
  1188. statement_half(a, b, c, z, False))
  1189. # Notations from [L], section 5.7-10
  1190. m, n, p, q = S([len(g.bm), len(g.an), len(g.ap), len(g.bq)])
  1191. tau = m + n - p
  1192. nu = q - m - n
  1193. rho = (tau - nu)/2
  1194. sigma = q - p
  1195. if sigma == 1:
  1196. epsilon = S.Half
  1197. elif sigma > 1:
  1198. epsilon = 1
  1199. else:
  1200. epsilon = S.NaN
  1201. theta = ((1 - sigma)/2 + Add(*g.bq) - Add(*g.ap))/sigma
  1202. delta = g.delta
  1203. _debug(' m=%s, n=%s, p=%s, q=%s, tau=%s, nu=%s, rho=%s, sigma=%s' % (
  1204. m, n, p, q, tau, nu, rho, sigma))
  1205. _debug(' epsilon=%s, theta=%s, delta=%s' % (epsilon, theta, delta))
  1206. # First check if the computation is valid.
  1207. if not (g.delta >= e/2 or (p >= 1 and p >= q)):
  1208. _debug(' Computation not valid for these parameters.')
  1209. return False
  1210. # Now check if the inversion integral exists.
  1211. # Test "condition A"
  1212. for a in g.an:
  1213. for b in g.bm:
  1214. if (a - b).is_integer and a > b:
  1215. _debug(' Not a valid G function.')
  1216. return False
  1217. # There are two cases. If p >= q, we can directly use a slater expansion
  1218. # like [L], 5.2 (11). Note in particular that the asymptotics of such an
  1219. # expansion even hold when some of the parameters differ by integers, i.e.
  1220. # the formula itself would not be valid! (b/c G functions are cts. in their
  1221. # parameters)
  1222. # When p < q, we need to use the theorems of [L], 5.10.
  1223. if p >= q:
  1224. _debug(' Using asymptotic Slater expansion.')
  1225. return And(*[statement(a - 1, 0, 0, z) for a in g.an])
  1226. def E(z):
  1227. return And(*[statement(a - 1, 0, 0, z) for a in g.an])
  1228. def H(z):
  1229. return statement(theta, -sigma, 1/sigma, z)
  1230. def Hp(z):
  1231. return statement_half(theta, -sigma, 1/sigma, z, True)
  1232. def Hm(z):
  1233. return statement_half(theta, -sigma, 1/sigma, z, False)
  1234. # [L], section 5.10
  1235. conds = []
  1236. # Theorem 1 -- p < q from test above
  1237. conds += [And(1 <= n, 1 <= m, rho*pi - delta >= pi/2, delta > 0,
  1238. E(z*exp(S.ImaginaryUnit*pi*(nu + 1))))]
  1239. # Theorem 2, statements (2) and (3)
  1240. conds += [And(p + 1 <= m, m + 1 <= q, delta > 0, delta < pi/2, n == 0,
  1241. (m - p + 1)*pi - delta >= pi/2,
  1242. Hp(z*exp(S.ImaginaryUnit*pi*(q - m))),
  1243. Hm(z*exp(-S.ImaginaryUnit*pi*(q - m))))]
  1244. # Theorem 2, statement (5) -- p < q from test above
  1245. conds += [And(m == q, n == 0, delta > 0,
  1246. (sigma + epsilon)*pi - delta >= pi/2, H(z))]
  1247. # Theorem 3, statements (6) and (7)
  1248. conds += [And(Or(And(p <= q - 2, 1 <= tau, tau <= sigma/2),
  1249. And(p + 1 <= m + n, m + n <= (p + q)/2)),
  1250. delta > 0, delta < pi/2, (tau + 1)*pi - delta >= pi/2,
  1251. Hp(z*exp(S.ImaginaryUnit*pi*nu)),
  1252. Hm(z*exp(-S.ImaginaryUnit*pi*nu)))]
  1253. # Theorem 4, statements (10) and (11) -- p < q from test above
  1254. conds += [And(1 <= m, rho > 0, delta > 0, delta + rho*pi < pi/2,
  1255. (tau + epsilon)*pi - delta >= pi/2,
  1256. Hp(z*exp(S.ImaginaryUnit*pi*nu)),
  1257. Hm(z*exp(-S.ImaginaryUnit*pi*nu)))]
  1258. # Trivial case
  1259. conds += [m == 0]
  1260. # TODO
  1261. # Theorem 5 is quite general
  1262. # Theorem 6 contains special cases for q=p+1
  1263. return Or(*conds)
  1264. def _int_inversion(g, x, t):
  1265. """
  1266. Compute the laplace inversion integral, assuming the formula applies.
  1267. """
  1268. b, a = _get_coeff_exp(g.argument, x)
  1269. C, g = _inflate_fox_h(meijerg(g.an, g.aother, g.bm, g.bother, b/t**a), -a)
  1270. return C/t*g
  1271. ####################################################################
  1272. # Finally, the real meat.
  1273. ####################################################################
  1274. _lookup_table = None
  1275. @cacheit
  1276. @timeit
  1277. def _rewrite_single(f, x, recursive=True):
  1278. """
  1279. Try to rewrite f as a sum of single G functions of the form
  1280. C*x**s*G(a*x**b), where b is a rational number and C is independent of x.
  1281. We guarantee that result.argument.as_coeff_mul(x) returns (a, (x**b,))
  1282. or (a, ()).
  1283. Returns a list of tuples (C, s, G) and a condition cond.
  1284. Returns None on failure.
  1285. """
  1286. from .transforms import (mellin_transform, inverse_mellin_transform,
  1287. IntegralTransformError, MellinTransformStripError)
  1288. global _lookup_table
  1289. if not _lookup_table:
  1290. _lookup_table = {}
  1291. _create_lookup_table(_lookup_table)
  1292. if isinstance(f, meijerg):
  1293. coeff, m = factor(f.argument, x).as_coeff_mul(x)
  1294. if len(m) > 1:
  1295. return None
  1296. m = m[0]
  1297. if m.is_Pow:
  1298. if m.base != x or not m.exp.is_Rational:
  1299. return None
  1300. elif m != x:
  1301. return None
  1302. return [(1, 0, meijerg(f.an, f.aother, f.bm, f.bother, coeff*m))], True
  1303. f_ = f
  1304. f = f.subs(x, z)
  1305. t = _mytype(f, z)
  1306. if t in _lookup_table:
  1307. l = _lookup_table[t]
  1308. for formula, terms, cond, hint in l:
  1309. subs = f.match(formula, old=True)
  1310. if subs:
  1311. subs_ = {}
  1312. for fro, to in subs.items():
  1313. subs_[fro] = unpolarify(polarify(to, lift=True),
  1314. exponents_only=True)
  1315. subs = subs_
  1316. if not isinstance(hint, bool):
  1317. hint = hint.subs(subs)
  1318. if hint == False:
  1319. continue
  1320. if not isinstance(cond, (bool, BooleanAtom)):
  1321. cond = unpolarify(cond.subs(subs))
  1322. if _eval_cond(cond) == False:
  1323. continue
  1324. if not isinstance(terms, list):
  1325. terms = terms(subs)
  1326. res = []
  1327. for fac, g in terms:
  1328. r1 = _get_coeff_exp(unpolarify(fac.subs(subs).subs(z, x),
  1329. exponents_only=True), x)
  1330. try:
  1331. g = g.subs(subs).subs(z, x)
  1332. except ValueError:
  1333. continue
  1334. # NOTE these substitutions can in principle introduce oo,
  1335. # zoo and other absurdities. It shouldn't matter,
  1336. # but better be safe.
  1337. if Tuple(*(r1 + (g,))).has(S.Infinity, S.ComplexInfinity, S.NegativeInfinity):
  1338. continue
  1339. g = meijerg(g.an, g.aother, g.bm, g.bother,
  1340. unpolarify(g.argument, exponents_only=True))
  1341. res.append(r1 + (g,))
  1342. if res:
  1343. return res, cond
  1344. # try recursive mellin transform
  1345. if not recursive:
  1346. return None
  1347. _debug('Trying recursive Mellin transform method.')
  1348. def my_imt(F, s, x, strip):
  1349. """ Calling simplify() all the time is slow and not helpful, since
  1350. most of the time it only factors things in a way that has to be
  1351. un-done anyway. But sometimes it can remove apparent poles. """
  1352. # XXX should this be in inverse_mellin_transform?
  1353. try:
  1354. return inverse_mellin_transform(F, s, x, strip,
  1355. as_meijerg=True, needeval=True)
  1356. except MellinTransformStripError:
  1357. return inverse_mellin_transform(
  1358. simplify(cancel(expand(F))), s, x, strip,
  1359. as_meijerg=True, needeval=True)
  1360. f = f_
  1361. s = _dummy('s', 'rewrite-single', f)
  1362. # to avoid infinite recursion, we have to force the two g functions case
  1363. def my_integrator(f, x):
  1364. r = _meijerint_definite_4(f, x, only_double=True)
  1365. if r is not None:
  1366. res, cond = r
  1367. res = _my_unpolarify(hyperexpand(res, rewrite='nonrepsmall'))
  1368. return Piecewise((res, cond),
  1369. (Integral(f, (x, S.Zero, S.Infinity)), True))
  1370. return Integral(f, (x, S.Zero, S.Infinity))
  1371. try:
  1372. F, strip, _ = mellin_transform(f, x, s, integrator=my_integrator,
  1373. simplify=False, needeval=True)
  1374. g = my_imt(F, s, x, strip)
  1375. except IntegralTransformError:
  1376. g = None
  1377. if g is None:
  1378. # We try to find an expression by analytic continuation.
  1379. # (also if the dummy is already in the expression, there is no point in
  1380. # putting in another one)
  1381. a = _dummy_('a', 'rewrite-single')
  1382. if a not in f.free_symbols and _is_analytic(f, x):
  1383. try:
  1384. F, strip, _ = mellin_transform(f.subs(x, a*x), x, s,
  1385. integrator=my_integrator,
  1386. needeval=True, simplify=False)
  1387. g = my_imt(F, s, x, strip).subs(a, 1)
  1388. except IntegralTransformError:
  1389. g = None
  1390. if g is None or g.has(S.Infinity, S.NaN, S.ComplexInfinity):
  1391. _debug('Recursive Mellin transform failed.')
  1392. return None
  1393. args = Add.make_args(g)
  1394. res = []
  1395. for f in args:
  1396. c, m = f.as_coeff_mul(x)
  1397. if len(m) > 1:
  1398. raise NotImplementedError('Unexpected form...')
  1399. g = m[0]
  1400. a, b = _get_coeff_exp(g.argument, x)
  1401. res += [(c, 0, meijerg(g.an, g.aother, g.bm, g.bother,
  1402. unpolarify(polarify(
  1403. a, lift=True), exponents_only=True)
  1404. *x**b))]
  1405. _debug('Recursive Mellin transform worked:', g)
  1406. return res, True
  1407. def _rewrite1(f, x, recursive=True):
  1408. """
  1409. Try to rewrite ``f`` using a (sum of) single G functions with argument a*x**b.
  1410. Return fac, po, g such that f = fac*po*g, fac is independent of ``x``.
  1411. and po = x**s.
  1412. Here g is a result from _rewrite_single.
  1413. Return None on failure.
  1414. """
  1415. fac, po, g = _split_mul(f, x)
  1416. g = _rewrite_single(g, x, recursive)
  1417. if g:
  1418. return fac, po, g[0], g[1]
  1419. def _rewrite2(f, x):
  1420. """
  1421. Try to rewrite ``f`` as a product of two G functions of arguments a*x**b.
  1422. Return fac, po, g1, g2 such that f = fac*po*g1*g2, where fac is
  1423. independent of x and po is x**s.
  1424. Here g1 and g2 are results of _rewrite_single.
  1425. Returns None on failure.
  1426. """
  1427. fac, po, g = _split_mul(f, x)
  1428. if any(_rewrite_single(expr, x, False) is None for expr in _mul_args(g)):
  1429. return None
  1430. l = _mul_as_two_parts(g)
  1431. if not l:
  1432. return None
  1433. l = list(ordered(l, [
  1434. lambda p: max(len(_exponents(p[0], x)), len(_exponents(p[1], x))),
  1435. lambda p: max(len(_functions(p[0], x)), len(_functions(p[1], x))),
  1436. lambda p: max(len(_find_splitting_points(p[0], x)),
  1437. len(_find_splitting_points(p[1], x)))]))
  1438. for recursive in [False, True]:
  1439. for fac1, fac2 in l:
  1440. g1 = _rewrite_single(fac1, x, recursive)
  1441. g2 = _rewrite_single(fac2, x, recursive)
  1442. if g1 and g2:
  1443. cond = And(g1[1], g2[1])
  1444. if cond != False:
  1445. return fac, po, g1[0], g2[0], cond
  1446. def meijerint_indefinite(f, x):
  1447. """
  1448. Compute an indefinite integral of ``f`` by rewriting it as a G function.
  1449. Examples
  1450. ========
  1451. >>> from sympy.integrals.meijerint import meijerint_indefinite
  1452. >>> from sympy import sin
  1453. >>> from sympy.abc import x
  1454. >>> meijerint_indefinite(sin(x), x)
  1455. -cos(x)
  1456. """
  1457. results = []
  1458. for a in sorted(_find_splitting_points(f, x) | {S.Zero}, key=default_sort_key):
  1459. res = _meijerint_indefinite_1(f.subs(x, x + a), x)
  1460. if not res:
  1461. continue
  1462. res = res.subs(x, x - a)
  1463. if _has(res, hyper, meijerg):
  1464. results.append(res)
  1465. else:
  1466. return res
  1467. if f.has(HyperbolicFunction):
  1468. _debug('Try rewriting hyperbolics in terms of exp.')
  1469. rv = meijerint_indefinite(
  1470. _rewrite_hyperbolics_as_exp(f), x)
  1471. if rv:
  1472. if not isinstance(rv, list):
  1473. return collect(factor_terms(rv), rv.atoms(exp))
  1474. results.extend(rv)
  1475. if results:
  1476. return next(ordered(results))
  1477. def _meijerint_indefinite_1(f, x):
  1478. """ Helper that does not attempt any substitution. """
  1479. _debug('Trying to compute the indefinite integral of', f, 'wrt', x)
  1480. gs = _rewrite1(f, x)
  1481. if gs is None:
  1482. # Note: the code that calls us will do expand() and try again
  1483. return None
  1484. fac, po, gl, cond = gs
  1485. _debug(' could rewrite:', gs)
  1486. res = S.Zero
  1487. for C, s, g in gl:
  1488. a, b = _get_coeff_exp(g.argument, x)
  1489. _, c = _get_coeff_exp(po, x)
  1490. c += s
  1491. # we do a substitution t=a*x**b, get integrand fac*t**rho*g
  1492. fac_ = fac * C / (b*a**((1 + c)/b))
  1493. rho = (c + 1)/b - 1
  1494. # we now use t**rho*G(params, t) = G(params + rho, t)
  1495. # [L, page 150, equation (4)]
  1496. # and integral G(params, t) dt = G(1, params+1, 0, t)
  1497. # (or a similar expression with 1 and 0 exchanged ... pick the one
  1498. # which yields a well-defined function)
  1499. # [R, section 5]
  1500. # (Note that this dummy will immediately go away again, so we
  1501. # can safely pass S.One for ``expr``.)
  1502. t = _dummy('t', 'meijerint-indefinite', S.One)
  1503. def tr(p):
  1504. return [a + rho + 1 for a in p]
  1505. if any(b.is_integer and (b <= 0) == True for b in tr(g.bm)):
  1506. r = -meijerg(
  1507. tr(g.an), tr(g.aother) + [1], tr(g.bm) + [0], tr(g.bother), t)
  1508. else:
  1509. r = meijerg(
  1510. tr(g.an) + [1], tr(g.aother), tr(g.bm), tr(g.bother) + [0], t)
  1511. # The antiderivative is most often expected to be defined
  1512. # in the neighborhood of x = 0.
  1513. if b.is_extended_nonnegative and not f.subs(x, 0).has(S.NaN, S.ComplexInfinity):
  1514. place = 0 # Assume we can expand at zero
  1515. else:
  1516. place = None
  1517. r = hyperexpand(r.subs(t, a*x**b), place=place)
  1518. # now substitute back
  1519. # Note: we really do want the powers of x to combine.
  1520. res += powdenest(fac_*r, polar=True)
  1521. def _clean(res):
  1522. """This multiplies out superfluous powers of x we created, and chops off
  1523. constants:
  1524. >> _clean(x*(exp(x)/x - 1/x) + 3)
  1525. exp(x)
  1526. cancel is used before mul_expand since it is possible for an
  1527. expression to have an additive constant that doesn't become isolated
  1528. with simple expansion. Such a situation was identified in issue 6369:
  1529. Examples
  1530. ========
  1531. >>> from sympy import sqrt, cancel
  1532. >>> from sympy.abc import x
  1533. >>> a = sqrt(2*x + 1)
  1534. >>> bad = (3*x*a**5 + 2*x - a**5 + 1)/a**2
  1535. >>> bad.expand().as_independent(x)[0]
  1536. 0
  1537. >>> cancel(bad).expand().as_independent(x)[0]
  1538. 1
  1539. """
  1540. res = expand_mul(cancel(res), deep=False)
  1541. return Add._from_args(res.as_coeff_add(x)[1])
  1542. res = piecewise_fold(res, evaluate=None)
  1543. if res.is_Piecewise:
  1544. newargs = []
  1545. for e, c in res.args:
  1546. e = _my_unpolarify(_clean(e))
  1547. newargs += [(e, c)]
  1548. res = Piecewise(*newargs, evaluate=False)
  1549. else:
  1550. res = _my_unpolarify(_clean(res))
  1551. return Piecewise((res, _my_unpolarify(cond)), (Integral(f, x), True))
  1552. @timeit
  1553. def meijerint_definite(f, x, a, b):
  1554. """
  1555. Integrate ``f`` over the interval [``a``, ``b``], by rewriting it as a product
  1556. of two G functions, or as a single G function.
  1557. Return res, cond, where cond are convergence conditions.
  1558. Examples
  1559. ========
  1560. >>> from sympy.integrals.meijerint import meijerint_definite
  1561. >>> from sympy import exp, oo
  1562. >>> from sympy.abc import x
  1563. >>> meijerint_definite(exp(-x**2), x, -oo, oo)
  1564. (sqrt(pi), True)
  1565. This function is implemented as a succession of functions
  1566. meijerint_definite, _meijerint_definite_2, _meijerint_definite_3,
  1567. _meijerint_definite_4. Each function in the list calls the next one
  1568. (presumably) several times. This means that calling meijerint_definite
  1569. can be very costly.
  1570. """
  1571. # This consists of three steps:
  1572. # 1) Change the integration limits to 0, oo
  1573. # 2) Rewrite in terms of G functions
  1574. # 3) Evaluate the integral
  1575. #
  1576. # There are usually several ways of doing this, and we want to try all.
  1577. # This function does (1), calls _meijerint_definite_2 for step (2).
  1578. _debug('Integrating', f, 'wrt %s from %s to %s.' % (x, a, b))
  1579. if f.has(DiracDelta):
  1580. _debug('Integrand has DiracDelta terms - giving up.')
  1581. return None
  1582. if f.has(SingularityFunction):
  1583. _debug('Integrand has Singularity Function terms - giving up.')
  1584. return None
  1585. f_, x_, a_, b_ = f, x, a, b
  1586. # Let's use a dummy in case any of the boundaries has x.
  1587. d = Dummy('x')
  1588. f = f.subs(x, d)
  1589. x = d
  1590. if a == b:
  1591. return (S.Zero, True)
  1592. results = []
  1593. if a is S.NegativeInfinity and b is not S.Infinity:
  1594. return meijerint_definite(f.subs(x, -x), x, -b, -a)
  1595. elif a is S.NegativeInfinity:
  1596. # Integrating -oo to oo. We need to find a place to split the integral.
  1597. _debug(' Integrating -oo to +oo.')
  1598. innermost = _find_splitting_points(f, x)
  1599. _debug(' Sensible splitting points:', innermost)
  1600. for c in sorted(innermost, key=default_sort_key, reverse=True) + [S.Zero]:
  1601. _debug(' Trying to split at', c)
  1602. if not c.is_extended_real:
  1603. _debug(' Non-real splitting point.')
  1604. continue
  1605. res1 = _meijerint_definite_2(f.subs(x, x + c), x)
  1606. if res1 is None:
  1607. _debug(' But could not compute first integral.')
  1608. continue
  1609. res2 = _meijerint_definite_2(f.subs(x, c - x), x)
  1610. if res2 is None:
  1611. _debug(' But could not compute second integral.')
  1612. continue
  1613. res1, cond1 = res1
  1614. res2, cond2 = res2
  1615. cond = _condsimp(And(cond1, cond2))
  1616. if cond == False:
  1617. _debug(' But combined condition is always false.')
  1618. continue
  1619. res = res1 + res2
  1620. return res, cond
  1621. elif a is S.Infinity:
  1622. res = meijerint_definite(f, x, b, S.Infinity)
  1623. return -res[0], res[1]
  1624. elif (a, b) == (S.Zero, S.Infinity):
  1625. # This is a common case - try it directly first.
  1626. res = _meijerint_definite_2(f, x)
  1627. if res:
  1628. if _has(res[0], meijerg):
  1629. results.append(res)
  1630. else:
  1631. return res
  1632. else:
  1633. if b is S.Infinity:
  1634. for split in _find_splitting_points(f, x):
  1635. if (a - split >= 0) == True:
  1636. _debug('Trying x -> x + %s' % split)
  1637. res = _meijerint_definite_2(f.subs(x, x + split)
  1638. *Heaviside(x + split - a), x)
  1639. if res:
  1640. if _has(res[0], meijerg):
  1641. results.append(res)
  1642. else:
  1643. return res
  1644. f = f.subs(x, x + a)
  1645. b = b - a
  1646. a = 0
  1647. if b is not S.Infinity:
  1648. phi = exp(S.ImaginaryUnit*arg(b))
  1649. b = Abs(b)
  1650. f = f.subs(x, phi*x)
  1651. f *= Heaviside(b - x)*phi
  1652. b = S.Infinity
  1653. _debug('Changed limits to', a, b)
  1654. _debug('Changed function to', f)
  1655. res = _meijerint_definite_2(f, x)
  1656. if res:
  1657. if _has(res[0], meijerg):
  1658. results.append(res)
  1659. else:
  1660. return res
  1661. if f_.has(HyperbolicFunction):
  1662. _debug('Try rewriting hyperbolics in terms of exp.')
  1663. rv = meijerint_definite(
  1664. _rewrite_hyperbolics_as_exp(f_), x_, a_, b_)
  1665. if rv:
  1666. if not isinstance(rv, list):
  1667. rv = (collect(factor_terms(rv[0]), rv[0].atoms(exp)),) + rv[1:]
  1668. return rv
  1669. results.extend(rv)
  1670. if results:
  1671. return next(ordered(results))
  1672. def _guess_expansion(f, x):
  1673. """ Try to guess sensible rewritings for integrand f(x). """
  1674. res = [(f, 'original integrand')]
  1675. orig = res[-1][0]
  1676. saw = {orig}
  1677. expanded = expand_mul(orig)
  1678. if expanded not in saw:
  1679. res += [(expanded, 'expand_mul')]
  1680. saw.add(expanded)
  1681. expanded = expand(orig)
  1682. if expanded not in saw:
  1683. res += [(expanded, 'expand')]
  1684. saw.add(expanded)
  1685. if orig.has(TrigonometricFunction, HyperbolicFunction):
  1686. expanded = expand_mul(expand_trig(orig))
  1687. if expanded not in saw:
  1688. res += [(expanded, 'expand_trig, expand_mul')]
  1689. saw.add(expanded)
  1690. if orig.has(cos, sin):
  1691. reduced = sincos_to_sum(orig)
  1692. if reduced not in saw:
  1693. res += [(reduced, 'trig power reduction')]
  1694. saw.add(reduced)
  1695. return res
  1696. def _meijerint_definite_2(f, x):
  1697. """
  1698. Try to integrate f dx from zero to infinity.
  1699. The body of this function computes various 'simplifications'
  1700. f1, f2, ... of f (e.g. by calling expand_mul(), trigexpand()
  1701. - see _guess_expansion) and calls _meijerint_definite_3 with each of
  1702. these in succession.
  1703. If _meijerint_definite_3 succeeds with any of the simplified functions,
  1704. returns this result.
  1705. """
  1706. # This function does preparation for (2), calls
  1707. # _meijerint_definite_3 for (2) and (3) combined.
  1708. # use a positive dummy - we integrate from 0 to oo
  1709. # XXX if a nonnegative symbol is used there will be test failures
  1710. dummy = _dummy('x', 'meijerint-definite2', f, positive=True)
  1711. f = f.subs(x, dummy)
  1712. x = dummy
  1713. if f == 0:
  1714. return S.Zero, True
  1715. for g, explanation in _guess_expansion(f, x):
  1716. _debug('Trying', explanation)
  1717. res = _meijerint_definite_3(g, x)
  1718. if res:
  1719. return res
  1720. def _meijerint_definite_3(f, x):
  1721. """
  1722. Try to integrate f dx from zero to infinity.
  1723. This function calls _meijerint_definite_4 to try to compute the
  1724. integral. If this fails, it tries using linearity.
  1725. """
  1726. res = _meijerint_definite_4(f, x)
  1727. if res and res[1] != False:
  1728. return res
  1729. if f.is_Add:
  1730. _debug('Expanding and evaluating all terms.')
  1731. ress = [_meijerint_definite_4(g, x) for g in f.args]
  1732. if all(r is not None for r in ress):
  1733. conds = []
  1734. res = S.Zero
  1735. for r, c in ress:
  1736. res += r
  1737. conds += [c]
  1738. c = And(*conds)
  1739. if c != False:
  1740. return res, c
  1741. def _my_unpolarify(f):
  1742. return _eval_cond(unpolarify(f))
  1743. @timeit
  1744. def _meijerint_definite_4(f, x, only_double=False):
  1745. """
  1746. Try to integrate f dx from zero to infinity.
  1747. Explanation
  1748. ===========
  1749. This function tries to apply the integration theorems found in literature,
  1750. i.e. it tries to rewrite f as either one or a product of two G-functions.
  1751. The parameter ``only_double`` is used internally in the recursive algorithm
  1752. to disable trying to rewrite f as a single G-function.
  1753. """
  1754. # This function does (2) and (3)
  1755. _debug('Integrating', f)
  1756. # Try single G function.
  1757. if not only_double:
  1758. gs = _rewrite1(f, x, recursive=False)
  1759. if gs is not None:
  1760. fac, po, g, cond = gs
  1761. _debug('Could rewrite as single G function:', fac, po, g)
  1762. res = S.Zero
  1763. for C, s, f in g:
  1764. if C == 0:
  1765. continue
  1766. C, f = _rewrite_saxena_1(fac*C, po*x**s, f, x)
  1767. res += C*_int0oo_1(f, x)
  1768. cond = And(cond, _check_antecedents_1(f, x))
  1769. if cond == False:
  1770. break
  1771. cond = _my_unpolarify(cond)
  1772. if cond == False:
  1773. _debug('But cond is always False.')
  1774. else:
  1775. _debug('Result before branch substitutions is:', res)
  1776. return _my_unpolarify(hyperexpand(res)), cond
  1777. # Try two G functions.
  1778. gs = _rewrite2(f, x)
  1779. if gs is not None:
  1780. for full_pb in [False, True]:
  1781. fac, po, g1, g2, cond = gs
  1782. _debug('Could rewrite as two G functions:', fac, po, g1, g2)
  1783. res = S.Zero
  1784. for C1, s1, f1 in g1:
  1785. for C2, s2, f2 in g2:
  1786. r = _rewrite_saxena(fac*C1*C2, po*x**(s1 + s2),
  1787. f1, f2, x, full_pb)
  1788. if r is None:
  1789. _debug('Non-rational exponents.')
  1790. return
  1791. C, f1_, f2_ = r
  1792. _debug('Saxena subst for yielded:', C, f1_, f2_)
  1793. cond = And(cond, _check_antecedents(f1_, f2_, x))
  1794. if cond == False:
  1795. break
  1796. res += C*_int0oo(f1_, f2_, x)
  1797. else:
  1798. continue
  1799. break
  1800. cond = _my_unpolarify(cond)
  1801. if cond == False:
  1802. _debug('But cond is always False (full_pb=%s).' % full_pb)
  1803. else:
  1804. _debug('Result before branch substitutions is:', res)
  1805. if only_double:
  1806. return res, cond
  1807. return _my_unpolarify(hyperexpand(res)), cond
  1808. def meijerint_inversion(f, x, t):
  1809. r"""
  1810. Compute the inverse laplace transform
  1811. $\int_{c+i\infty}^{c-i\infty} f(x) e^{tx}\, dx$,
  1812. for real c larger than the real part of all singularities of ``f``.
  1813. Note that ``t`` is always assumed real and positive.
  1814. Return None if the integral does not exist or could not be evaluated.
  1815. Examples
  1816. ========
  1817. >>> from sympy.abc import x, t
  1818. >>> from sympy.integrals.meijerint import meijerint_inversion
  1819. >>> meijerint_inversion(1/x, x, t)
  1820. Heaviside(t)
  1821. """
  1822. f_ = f
  1823. t_ = t
  1824. t = Dummy('t', polar=True) # We don't want sqrt(t**2) = abs(t) etc
  1825. f = f.subs(t_, t)
  1826. _debug('Laplace-inverting', f)
  1827. if not _is_analytic(f, x):
  1828. _debug('But expression is not analytic.')
  1829. return None
  1830. # Exponentials correspond to shifts; we filter them out and then
  1831. # shift the result later. If we are given an Add this will not
  1832. # work, but the calling code will take care of that.
  1833. shift = S.Zero
  1834. if f.is_Mul:
  1835. args = list(f.args)
  1836. elif isinstance(f, exp):
  1837. args = [f]
  1838. else:
  1839. args = None
  1840. if args:
  1841. newargs = []
  1842. exponentials = []
  1843. while args:
  1844. arg = args.pop()
  1845. if isinstance(arg, exp):
  1846. arg2 = expand(arg)
  1847. if arg2.is_Mul:
  1848. args += arg2.args
  1849. continue
  1850. try:
  1851. a, b = _get_coeff_exp(arg.args[0], x)
  1852. except _CoeffExpValueError:
  1853. b = 0
  1854. if b == 1:
  1855. exponentials.append(a)
  1856. else:
  1857. newargs.append(arg)
  1858. elif arg.is_Pow:
  1859. arg2 = expand(arg)
  1860. if arg2.is_Mul:
  1861. args += arg2.args
  1862. continue
  1863. if x not in arg.base.free_symbols:
  1864. try:
  1865. a, b = _get_coeff_exp(arg.exp, x)
  1866. except _CoeffExpValueError:
  1867. b = 0
  1868. if b == 1:
  1869. exponentials.append(a*log(arg.base))
  1870. newargs.append(arg)
  1871. else:
  1872. newargs.append(arg)
  1873. shift = Add(*exponentials)
  1874. f = Mul(*newargs)
  1875. if x not in f.free_symbols:
  1876. _debug('Expression consists of constant and exp shift:', f, shift)
  1877. cond = Eq(im(shift), 0)
  1878. if cond == False:
  1879. _debug('but shift is nonreal, cannot be a Laplace transform')
  1880. return None
  1881. res = f*DiracDelta(t + shift)
  1882. _debug('Result is a delta function, possibly conditional:', res, cond)
  1883. # cond is True or Eq
  1884. return Piecewise((res.subs(t, t_), cond))
  1885. gs = _rewrite1(f, x)
  1886. if gs is not None:
  1887. fac, po, g, cond = gs
  1888. _debug('Could rewrite as single G function:', fac, po, g)
  1889. res = S.Zero
  1890. for C, s, f in g:
  1891. C, f = _rewrite_inversion(fac*C, po*x**s, f, x)
  1892. res += C*_int_inversion(f, x, t)
  1893. cond = And(cond, _check_antecedents_inversion(f, x))
  1894. if cond == False:
  1895. break
  1896. cond = _my_unpolarify(cond)
  1897. if cond == False:
  1898. _debug('But cond is always False.')
  1899. else:
  1900. _debug('Result before branch substitution:', res)
  1901. res = _my_unpolarify(hyperexpand(res))
  1902. if not res.has(Heaviside):
  1903. res *= Heaviside(t)
  1904. res = res.subs(t, t + shift)
  1905. if not isinstance(cond, bool):
  1906. cond = cond.subs(t, t + shift)
  1907. from .transforms import InverseLaplaceTransform
  1908. return Piecewise((res.subs(t, t_), cond),
  1909. (InverseLaplaceTransform(f_.subs(t, t_), x, t_, None), True))