meijerint.py 79 KB

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