Code coverage tests

This page documents the degree to which the PARI/GP source code is tested by our public test suite, distributed with the source distribution in directory src/test/. This is measured by the gcov utility; we then process gcov output using the lcov frond-end.

We test a few variants depending on Configure flags on the pari.math.u-bordeaux.fr machine (x86_64 architecture), and agregate them in the final report:

The target is 90% coverage for all mathematical modules (given that branches depending on DEBUGLEVEL or DEBUGMEM are not covered). This script is run to produce the results below.

LCOV - code coverage report
Current view: top level - basemath - trans3.c (source / functions) Hit Total Coverage
Test: PARI/GP v2.10.0 lcov report (development 21059-cbe0d6a) Lines: 1800 1884 95.5 %
Date: 2017-09-22 06:24:58 Functions: 114 116 98.3 %
Legend: Lines: hit not hit

          Line data    Source code
       1             : /* Copyright (C) 2000  The PARI group.
       2             : 
       3             : This file is part of the PARI/GP package.
       4             : 
       5             : PARI/GP is free software; you can redistribute it and/or modify it under the
       6             : terms of the GNU General Public License as published by the Free Software
       7             : Foundation. It is distributed in the hope that it will be useful, but WITHOUT
       8             : ANY WARRANTY WHATSOEVER.
       9             : 
      10             : Check the License for details. You should have received a copy of it, along
      11             : with the package; see the file 'COPYING'. If not, write to the Free Software
      12             : Foundation, Inc., 51 Franklin Street, Fifth Floor, Boston, MA 02110-1301 USA. */
      13             : 
      14             : /********************************************************************/
      15             : /**                                                                **/
      16             : /**                   TRANSCENDENTAL FONCTIONS                     **/
      17             : /**                          (part 3)                              **/
      18             : /**                                                                **/
      19             : /********************************************************************/
      20             : #include "pari.h"
      21             : #include "paripriv.h"
      22             : 
      23             : #define HALF_E 1.3591409 /* Exponential / 2 */
      24             : static const long EXTRAPREC =
      25             : #ifdef LONG_IS_64BIT
      26             :   1;
      27             : #else
      28             :   2;
      29             : #endif
      30             : 
      31             : /***********************************************************************/
      32             : /**                                                                   **/
      33             : /**                       BESSEL FUNCTIONS                            **/
      34             : /**                                                                   **/
      35             : /***********************************************************************/
      36             : 
      37             : /* n! sum_{k=0}^m Z^k / (k!*(k+n)!), with Z := (-1)^flag*z^2/4 */
      38             : static GEN
      39        4074 : _jbessel(GEN n, GEN z, long flag, long m)
      40             : {
      41             :   pari_sp av;
      42             :   GEN Z,s;
      43             :   long k;
      44             : 
      45        4074 :   Z = gmul2n(gsqr(z),-2); if (flag & 1) Z = gneg(Z);
      46        4074 :   if (typ(z) == t_SER)
      47             :   {
      48         210 :     long v = valp(z);
      49         210 :     if (v < 0) pari_err_DOMAIN("besselj","valuation", "<", gen_0, z);
      50         203 :     k = lg(Z)-2 - v;
      51         203 :     if (v == 0) pari_err_IMPL("besselj around a!=0");
      52         203 :     if (k <= 0) return scalarser(gen_1, varn(z), 2*v);
      53         203 :     setlg(Z, k+2);
      54             :   }
      55        4067 :   s = gen_1;
      56        4067 :   av = avma;
      57      183341 :   for (k=m; k>=1; k--)
      58             :   {
      59      179274 :     s = gaddsg(1, gdiv(gmul(Z,s), gmulgs(gaddgs(n,k),k)));
      60      179274 :     if (gc_needed(av,1))
      61             :     {
      62           0 :       if (DEBUGMEM>1) pari_warn(warnmem,"besselj");
      63           0 :       s = gerepileupto(av, s);
      64             :     }
      65             :   }
      66        4067 :   return s;
      67             : }
      68             : 
      69             : /* return L * approximate solution to x log x = B */
      70             : static long
      71        6566 : bessel_get_lim(double B, double L)
      72             : {
      73             :   long lim;
      74        6566 :   double x = 1 + B; /* 3 iterations are enough except in pathological cases */
      75        6566 :   x = (x + B)/(log(x)+1);
      76        6566 :   x = (x + B)/(log(x)+1);
      77        6566 :   x = (x + B)/(log(x)+1); x = L*x;
      78        6566 :   lim = (long)x; if (lim < 2) lim = 2;
      79        6566 :   return lim;
      80             : }
      81             : 
      82             : static GEN jbesselintern(GEN n, GEN z, long flag, long prec);
      83             : static GEN kbesselintern(GEN n, GEN z, long flag, long prec);
      84             : static GEN
      85         140 : jbesselvec(GEN n, GEN x, long fl, long prec)
      86             : {
      87             :   long i, l;
      88         140 :   GEN y = cgetg_copy(x, &l);
      89         140 :   for (i=1; i<l; i++) gel(y,i) = jbesselintern(n,gel(x,i),fl,prec);
      90         140 :   return y;
      91             : }
      92             : static GEN
      93          35 : jbesselhvec(GEN n, GEN x, long prec)
      94             : {
      95             :   long i, l;
      96          35 :   GEN y = cgetg_copy(x, &l);
      97          35 :   for (i=1; i<l; i++) gel(y,i) = jbesselh(n,gel(x,i),prec);
      98          35 :   return y;
      99             : }
     100             : static GEN
     101          14 : polylogvec(long m, GEN x, long prec)
     102             : {
     103             :   long l, i;
     104          14 :   GEN y = cgetg_copy(x, &l);
     105          14 :   for (i=1; i<l; i++) gel(y,i) = gpolylog(m,gel(x,i),prec);
     106          14 :   return y;
     107             : }
     108             : static GEN
     109         140 : kbesselvec(GEN n, GEN x, long fl, long prec)
     110             : {
     111             :   long i, l;
     112         140 :   GEN y = cgetg_copy(x, &l);
     113         140 :   for (i=1; i<l; i++) gel(y,i) = kbesselintern(n,gel(x,i),fl,prec);
     114         140 :   return y;
     115             : }
     116             : 
     117             : static GEN
     118        4438 : jbesselintern(GEN n, GEN z, long flag, long prec)
     119             : {
     120             :   long i, ki;
     121        4438 :   pari_sp av = avma;
     122             :   GEN y;
     123             : 
     124        4438 :   switch(typ(z))
     125             :   {
     126             :     case t_INT: case t_FRAC: case t_QUAD:
     127             :     case t_REAL: case t_COMPLEX:
     128             :     {
     129        4004 :       int flz0 = gequal0(z);
     130             :       long lim, k, precnew;
     131             :       GEN p1, p2;
     132             :       double B, L;
     133             : 
     134        4004 :       i = precision(z); if (i) prec = i;
     135        4004 :       if (flz0 && gequal0(n)) return real_1(prec);
     136        4004 :       p2 = gpow(gmul2n(z,-1),n,prec);
     137        3948 :       p2 = gdiv(p2, ggamma(gaddgs(n,1),prec));
     138        3948 :       if (flz0) return gerepileupto(av, p2);
     139        3948 :       L = HALF_E * gtodouble(gabs(gtofp(z,LOWDEFAULTPREC),prec));
     140        3948 :       precnew = prec;
     141        3948 :       if (L >= 1.0)
     142        3493 :         precnew += nbits2extraprec((long)(L/(HALF_E*LOG2) + BITS_IN_LONG));
     143        3948 :       if (issmall(n,&ki)) {
     144        3094 :         k = labs(ki);
     145        3094 :         n = utoi(k);
     146             :       } else {
     147         770 :         i = precision(n);
     148         770 :         if (i && i < precnew) n = gtofp(n,precnew);
     149             :       }
     150        3864 :       z = gtofp(z,precnew);
     151        3864 :       B = prec2nbits_mul(prec, LOG2/2) / L;
     152        3864 :       lim = bessel_get_lim(B, L);
     153        3864 :       p1 = gprec_wtrunc(_jbessel(n,z,flag,lim), prec);
     154        3864 :       return gerepileupto(av, gmul(p2,p1));
     155             :     }
     156             : 
     157             :     case t_VEC: case t_COL: case t_MAT:
     158         112 :       return jbesselvec(n, z, flag, prec);
     159             : 
     160             :     case t_POLMOD:
     161          28 :       y = jbesselvec(n, polmod_to_embed(z, prec), flag, prec);
     162          28 :       return gerepileupto(av,y);
     163             : 
     164          28 :     case t_PADIC: pari_err_IMPL("p-adic jbessel function");
     165             :     default:
     166         266 :       if (!(y = toser_i(z))) break;
     167         238 :       if (issmall(n,&ki)) n = utoi(labs(ki));
     168         210 :       return gerepileupto(av, _jbessel(n,y,flag,lg(y)-2));
     169             :   }
     170          28 :   pari_err_TYPE("jbessel",z);
     171             :   return NULL; /* LCOV_EXCL_LINE */
     172             : }
     173             : GEN
     174         763 : jbessel(GEN n, GEN z, long prec) { return jbesselintern(n,z,1,prec); }
     175             : GEN
     176         252 : ibessel(GEN n, GEN z, long prec) { return jbesselintern(n,z,0,prec); }
     177             : 
     178             : static GEN
     179         140 : _jbesselh(long k, GEN z, long prec)
     180             : {
     181         140 :   GEN s,c,p0,p1,p2, zinv = ginv(z);
     182             :   long i;
     183             : 
     184         140 :   gsincos(z,&s,&c,prec);
     185         140 :   p1 = gmul(zinv,s);
     186         140 :   if (k)
     187             :   {
     188         119 :     p0 = p1; p1 = gmul(zinv,gsub(p0,c));
     189         441 :     for (i=2; i<=k; i++)
     190             :     {
     191         322 :       p2 = gsub(gmul(gmulsg(2*i-1,zinv),p1), p0);
     192         322 :       p0 = p1; p1 = p2;
     193             :     }
     194             :   }
     195         140 :   return p1;
     196             : }
     197             : 
     198             : /* J_{n+1/2}(z) */
     199             : GEN
     200         336 : jbesselh(GEN n, GEN z, long prec)
     201             : {
     202             :   long k, i;
     203             :   pari_sp av;
     204             :   GEN y;
     205             : 
     206         336 :   if (typ(n)!=t_INT) pari_err_TYPE("jbesselh",n);
     207         210 :   k = itos(n);
     208         210 :   if (k < 0) return jbessel(gadd(ghalf,n), z, prec);
     209             : 
     210         210 :   switch(typ(z))
     211             :   {
     212             :     case t_INT: case t_FRAC: case t_QUAD:
     213             :     case t_REAL: case t_COMPLEX:
     214             :     {
     215             :       long bits, precnew, gz, pr;
     216             :       GEN p1;
     217         126 :       if (gequal0(z))
     218             :       {
     219           7 :         av = avma;
     220           7 :         p1 = gmul(gsqrt(gdiv(z,mppi(prec)),prec),gpowgs(z,k));
     221           7 :         p1 = gdiv(p1, mulu_interval(k+1, 2*k+1)); /* x k! / (2k+1)! */
     222           7 :         return gerepileupto(av, gmul2n(p1,2*k));
     223             :       }
     224         119 :       gz = gexpo(z);
     225         119 :       if ( (pr = precision(z)) ) prec = pr;
     226         119 :       y = cgetc(prec);
     227         119 :       bits = -2*k*gz + BITS_IN_LONG;
     228         119 :       av = avma;
     229         119 :       if (bits <= 0)
     230           7 :         precnew = prec;
     231             :       else
     232             :       {
     233         112 :         precnew = prec + nbits2extraprec(bits);
     234         112 :         if (pr) z = gtofp(z, precnew);
     235             :       }
     236         119 :       p1 = gmul(_jbesselh(k,z,prec), gsqrt(gdiv(z,Pi2n(-1,prec)),prec));
     237         119 :       avma = av; return affc_fixlg(p1, y);
     238             :     }
     239             : 
     240             :     case t_VEC: case t_COL: case t_MAT:
     241          28 :       return jbesselhvec(n, z, prec);
     242             : 
     243             :     case t_POLMOD:
     244           7 :       av = avma;
     245           7 :       return gerepileupto(av, jbesselhvec(n, polmod_to_embed(z,prec), prec));
     246             : 
     247           7 :     case t_PADIC: pari_err_IMPL("p-adic jbesselh function");
     248             :     default:
     249             :     {
     250             :       long t, v;
     251          42 :       av = avma; if (!(y = toser_i(z))) break;
     252          35 :       if (gequal0(y)) return gerepileupto(av, gpowgs(y,k));
     253          35 :       v = valp(y);
     254          35 :       if (v < 0) pari_err_DOMAIN("besseljh","valuation", "<", gen_0, y);
     255          28 :       t = lg(y)-2;
     256          28 :       if (v) y = sertoser(y, t + (2*k+1)*v);
     257          28 :       if (!k)
     258           7 :         y = gdiv(gsin(y,prec), y);
     259             :       else
     260             :       {
     261          21 :         GEN T, a = _jbesselh(k, y, prec);
     262          21 :         if (v) y = sertoser(y, t + k*v); /* lower precision */
     263          21 :         y = gdiv(a, gpowgs(y, k));
     264          21 :         T = cgetg(k+1, t_VECSMALL);
     265          21 :         for (i = 1; i <= k; i++) T[i] = 2*i+1;
     266          21 :         y = gmul(y, zv_prod_Z(T));
     267             :       }
     268          28 :       return gerepileupto(av, y);
     269             :     }
     270             :   }
     271           7 :   pari_err_TYPE("besseljh",z);
     272             :   return NULL; /* LCOV_EXCL_LINE */
     273             : }
     274             : 
     275             : static GEN
     276           7 : kbessel2(GEN nu, GEN x, long prec)
     277             : {
     278           7 :   pari_sp av = avma;
     279             :   GEN p1, x2, a;
     280             : 
     281           7 :   if (typ(x)==t_REAL) prec = realprec(x);
     282           7 :   x2 = gshift(x,1);
     283           7 :   a = gtofp(gaddgs(gshift(nu,1), 1), prec);
     284           7 :   p1 = hyperu(gshift(a,-1),a,x2,prec);
     285           7 :   p1 = gmul(gmul(p1,gpow(x2,nu,prec)), sqrtr(mppi(prec)));
     286           7 :   return gerepileupto(av, gmul(p1,gexp(gneg(x),prec)));
     287             : }
     288             : 
     289             : static GEN
     290          28 : kbessel1(GEN nu, GEN gx, long prec)
     291             : {
     292             :   GEN x, y, p1, zf, zz, r, s, t, u, ak, pitemp, nu2;
     293             :   long l, lnew, k, k2, l1, n2, n, ex;
     294             :   pari_sp av;
     295             : 
     296          28 :   if (typ(nu)==t_COMPLEX) return kbessel2(nu,gx,prec);
     297          21 :   l = (typ(gx)==t_REAL)? realprec(gx): prec;
     298          21 :   ex = gexpo(gx);
     299          21 :   if (ex < 0)
     300             :   {
     301           7 :     long rab = nbits2extraprec(-ex);
     302           7 :     lnew = l + rab; prec += rab;
     303             :   }
     304          14 :   else lnew = l;
     305          21 :   y = cgetr(l); l1=lnew+1;
     306          21 :   av = avma; x = gtofp(gx, lnew); nu = gtofp(nu, lnew);
     307          21 :   nu2 = gmul2n(sqrr(nu), 2); togglesign(nu2);
     308          21 :   n = (long) (prec2nbits_mul(l,LOG2) + M_PI*fabs(rtodbl(nu))) / 2;
     309          21 :   n2 = n<<1; pitemp=mppi(l1);
     310          21 :   r = gmul2n(x,1);
     311          21 :   if (cmprs(x, n) < 0)
     312             :   {
     313          14 :     GEN q = stor(n2, l1), v, c, e, f;
     314             :     pari_sp av1, av2;
     315          14 :     u=cgetr(l1); v=cgetr(l1); e=cgetr(l1); f=cgetr(l1);
     316          14 :     av1 = avma;
     317          14 :     zf = sqrtr(divru(pitemp,n2));
     318          14 :     zz = invr(stor(n2<<2, prec));
     319          14 :     s = real_1(prec); t = real_0(prec);
     320        1246 :     for (k=n2,k2=2*n2-1; k > 0; k--,k2-=2)
     321             :     {
     322        1232 :       p1 = addri(nu2, sqrs(k2));
     323        1232 :       ak = divrs(mulrr(p1,zz),-k);
     324        1232 :       s = addsr(1, mulrr(ak,s));
     325        1232 :       t = addsr(k2,mulrr(ak,t));
     326             :     }
     327          14 :     mulrrz(zf, s, u); shiftr_inplace(t, -1);
     328          14 :     divrsz(addrr(mulrr(t,zf),mulrr(u,nu)),-n2,v);
     329         301 :     for(;; avma = av1)
     330             :     {
     331         315 :       GEN d = real_1(l1);
     332         315 :       c = divur(5,q); if (expo(c) >= -1) c = real2n(-1,l1);
     333         315 :       p1 = subsr(1,divrr(r,q)); if (cmprr(c,p1)>0) c = p1;
     334         315 :       togglesign(c);
     335         315 :       affrr(u,e);
     336         315 :       affrr(v,f); av2 = avma;
     337       45267 :       for (k=1;; k++, avma=av2)
     338             :       {
     339       45267 :         GEN w = addrr(gmul2n(mulur(2*k-1,u), -1), mulrr(subrs(q,k),v));
     340       45267 :         w = addrr(w, mulrr(nu, subrr(u,gmul2n(v,1))));
     341       45267 :         divrsz(mulrr(q,v),k,u);
     342       45267 :         divrsz(w,k,v);
     343       45267 :         mulrrz(d,c,d);
     344       45267 :         addrrz(e,mulrr(d,u),e); p1=mulrr(d,v);
     345       45267 :         addrrz(f,p1,f);
     346       45267 :         if (gexpo(p1) - gexpo(f) <= 1-prec2nbits(precision(p1))) break;
     347             : 
     348       44952 :       }
     349         315 :       swap(e, u);
     350         315 :       swap(f, v);
     351         315 :       affrr(mulrr(q,addrs(c,1)), q);
     352         315 :       if (expo(subrr(q,r)) - expo(r) <= 1-prec2nbits(lnew)) break;
     353         301 :     }
     354          14 :     u = mulrr(u, gpow(divru(x,n),nu,prec));
     355             :   }
     356             :   else
     357             :   {
     358           7 :     zf = sqrtr(divrr(pitemp,r));
     359           7 :     zz = ginv(gmul2n(r,2)); s = real_1(prec);
     360         637 :     for (k=n2,k2=2*n2-1; k > 0; k--,k2-=2)
     361             :     {
     362         630 :       p1 = addri(nu2, sqrs(k2));
     363         630 :       ak = divru(mulrr(p1,zz), k);
     364         630 :       s = subsr(1, mulrr(ak,s));
     365             :     }
     366           7 :     u = mulrr(s, zf);
     367             :   }
     368          21 :   affrr(mulrr(u, mpexp(mpneg(x))), y);
     369          21 :   avma = av; return y;
     370             : }
     371             : 
     372             : /*   sum_{k=0}^m Z^k (H(k)+H(k+n)) / (k! (k+n)!)
     373             :  * + sum_{k=0}^{n-1} (-Z)^(k-n) (n-k-1)!/k!   with Z := (-1)^flag*z^2/4.
     374             :  * Warning: contrary to _jbessel, no n! in front.
     375             :  * When flag > 1, compute exactly the H(k) and factorials (slow) */
     376             : static GEN
     377        2800 : _kbessel1(long n, GEN z, long flag, long m, long prec)
     378             : {
     379             :   GEN Z, p1, p2, s, H;
     380             :   pari_sp av;
     381             :   long k;
     382             : 
     383        2800 :   Z = gmul2n(gsqr(z),-2); if (flag & 1) Z = gneg(Z);
     384        2800 :   if (typ(z) == t_SER)
     385             :   {
     386          98 :     long v = valp(z);
     387          98 :     if (v < 0) pari_err_DOMAIN("_kbessel1","valuation", "<", gen_0, z);
     388          84 :     k = lg(Z)-2 - v;
     389          84 :     if (v == 0) pari_err_IMPL("Bessel K around a!=0");
     390          84 :     if (k <= 0) return gadd(gen_1, zeroser(varn(z), 2*v));
     391          84 :     setlg(Z, k+2);
     392             :   }
     393        2786 :   H = cgetg(m+n+2,t_VEC); gel(H,1) = gen_0;
     394        2786 :   if (flag <= 1)
     395             :   {
     396        2702 :     gel(H,2) = s = real_1(prec);
     397        2702 :     for (k=2; k<=m+n; k++) gel(H,k+1) = s = divru(addsr(1,mulur(k,s)),k);
     398             :   }
     399             :   else
     400             :   {
     401          84 :     gel(H,2) = s = gen_1;
     402          84 :     for (k=2; k<=m+n; k++) gel(H,k+1) = s = gdivgs(gaddsg(1,gmulsg(k,s)),k);
     403             :   }
     404        2786 :   s = gadd(gel(H,m+1), gel(H,m+n+1));
     405        2786 :   av = avma;
     406      101066 :   for (k=m; k>0; k--)
     407             :   {
     408       98280 :     s = gadd(gadd(gel(H,k),gel(H,k+n)),gdiv(gmul(Z,s),mulss(k,k+n)));
     409       98280 :     if (gc_needed(av,1))
     410             :     {
     411           0 :       if (DEBUGMEM>1) pari_warn(warnmem,"_kbessel1");
     412           0 :       s = gerepileupto(av, s);
     413             :     }
     414             :   }
     415        2786 :   p1 = (flag <= 1) ? mpfactr(n,prec) : mpfact(n);
     416        2786 :   s = gdiv(s,p1);
     417        2786 :   if (n)
     418             :   {
     419         364 :     Z = gneg(ginv(Z));
     420         364 :     p2 = gmulsg(n, gdiv(Z,p1));
     421         364 :     s = gadd(s,p2);
     422        1372 :     for (k=n-1; k>0; k--)
     423             :     {
     424        1008 :       p2 = gmul(p2, gmul(mulss(k,n-k),Z));
     425        1008 :       s = gadd(s,p2);
     426             :     }
     427             :   }
     428        2786 :   return s;
     429             : }
     430             : 
     431             : /* flag = 0: K / flag = 1: N */
     432             : static GEN
     433        3437 : kbesselintern(GEN n, GEN z, long flag, long prec)
     434             : {
     435        3437 :   const char *f = flag? "besseln": "besselk";
     436        3437 :   const long flK = (flag == 0);
     437             :   long i, k, ki, lim, precnew, fl2, ex;
     438        3437 :   pari_sp av = avma;
     439             :   GEN p1, p2, y, p3, pp, pm, s, c;
     440             :   double B, L;
     441             : 
     442        3437 :   switch(typ(z))
     443             :   {
     444             :     case t_INT: case t_FRAC: case t_QUAD:
     445             :     case t_REAL: case t_COMPLEX:
     446        3059 :       if (gequal0(z)) pari_err_DOMAIN(f, "argument", "=", gen_0, z);
     447        3059 :       i = precision(z); if (i) prec = i;
     448        3059 :       i = precision(n); if (i && prec > i) prec = i;
     449        3059 :       ex = gexpo(z);
     450             :       /* experimental */
     451        3059 :       if (!flag && !gequal0(n) && ex > prec2nbits(prec)/16 + gexpo(n))
     452          28 :         return kbessel1(n,z,prec);
     453        3024 :       L = HALF_E * gtodouble(gabs(z,prec));
     454        3024 :       precnew = prec;
     455        3024 :       if (L >= HALF_E) {
     456        2835 :         long rab = nbits2extraprec((long) (L/(HALF_E*LOG2)));
     457        2835 :         if (flK) rab *= 2;
     458        2835 :          precnew += 1 + rab;
     459             :       }
     460        3024 :       z = gtofp(z, precnew);
     461        3024 :       if (issmall(n,&ki))
     462             :       {
     463        2702 :         GEN z2 = gmul2n(z, -1);
     464        2702 :         k = labs(ki);
     465        2702 :         B = prec2nbits_mul(prec,LOG2/2) / L;
     466        2702 :         if (flK) B += 0.367879;
     467        2702 :         lim = bessel_get_lim(B, L);
     468        2702 :         p1 = gmul(gpowgs(z2,k), _kbessel1(k,z,flag,lim,precnew));
     469        2702 :         p2 = gadd(mpeuler(precnew), glog(z2,precnew));
     470        2702 :         p3 = jbesselintern(stoi(k),z,flag,precnew);
     471        2702 :         p2 = gsub(gmul2n(p1,-1),gmul(p2,p3));
     472        2702 :         p2 = gprec_wtrunc(p2, prec);
     473        2702 :         if (flK) {
     474        2408 :           if (k & 1) p2 = gneg(p2);
     475             :         }
     476             :         else
     477             :         {
     478         294 :           p2 = gdiv(p2, Pi2n(-1,prec));
     479         294 :           if (ki >= 0 || (k&1)==0) p2 = gneg(p2);
     480             :         }
     481        2702 :         return gerepilecopy(av, p2);
     482             :       }
     483             : 
     484         259 :       n = gtofp(n, precnew);
     485         259 :       gsincos(gmul(n,mppi(precnew)), &s,&c,precnew);
     486         259 :       ex = gexpo(s);
     487         259 :       if (ex < 0)
     488             :       {
     489          84 :         long rab = nbits2extraprec(-ex);
     490          84 :         if (flK) rab *= 2;
     491          84 :         precnew += rab;
     492             :       }
     493         259 :       if (i && i < precnew) {
     494          84 :         n = gtofp(n,precnew);
     495          84 :         z = gtofp(z,precnew);
     496          84 :         gsincos(gmul(n,mppi(precnew)), &s,&c,precnew);
     497             :       }
     498             : 
     499         259 :       pp = jbesselintern(n,      z,flag,precnew);
     500         259 :       pm = jbesselintern(gneg(n),z,flag,precnew);
     501         259 :       if (flK)
     502          70 :         p1 = gmul(gsub(pm,pp), Pi2n(-1,precnew));
     503             :       else
     504         189 :         p1 = gsub(gmul(c,pp),pm);
     505         259 :       p1 = gdiv(p1, s);
     506         259 :       return gerepilecopy(av, gprec_wtrunc(p1,prec));
     507             : 
     508             :     case t_VEC: case t_COL: case t_MAT:
     509         112 :       return kbesselvec(n,z,flag,prec);
     510             : 
     511             :     case t_POLMOD:
     512          28 :       y = kbesselvec(n,polmod_to_embed(z,prec),flag,prec);
     513          28 :       return gerepileupto(av, y);
     514             : 
     515          14 :     case t_PADIC: pari_err_IMPL(stack_strcat("p-adic ",f));
     516             :     default:
     517         224 :       if (!(y = toser_i(z))) break;
     518         210 :       if (issmall(n,&ki))
     519             :       {
     520          98 :         k = labs(ki);
     521          98 :         return gerepilecopy(av, _kbessel1(k,y,flag+2,lg(y)-2,prec));
     522             :       }
     523          98 :       if (!issmall(gmul2n(n,1),&ki))
     524          70 :         pari_err_DOMAIN(f, "2n mod Z", "!=", gen_0,n);
     525          28 :       k = labs(ki); n = gmul2n(stoi(k),-1);
     526          28 :       fl2 = (k&3)==1;
     527          28 :       pm = jbesselintern(gneg(n),y,flag,prec);
     528          28 :       if (flK)
     529             :       {
     530           7 :         pp = jbesselintern(n,y,flag,prec);
     531           7 :         p2 = gpowgs(y,-k); if (fl2 == 0) p2 = gneg(p2);
     532           7 :         p3 = gmul2n(diviiexact(mpfact(k + 1),mpfact((k + 1) >> 1)),-(k + 1));
     533           7 :         p3 = gdivgs(gmul2n(gsqr(p3),1),k);
     534           7 :         p2 = gmul(p2,p3);
     535           7 :         p1 = gsub(pp,gmul(p2,pm));
     536             :       }
     537          21 :       else p1 = pm;
     538          28 :       return gerepileupto(av, fl2? gneg(p1): gcopy(p1));
     539             :   }
     540          14 :   pari_err_TYPE(f,z);
     541             :   return NULL; /* LCOV_EXCL_LINE */
     542             : }
     543             : 
     544             : GEN
     545        2611 : kbessel(GEN n, GEN z, long prec) { return kbesselintern(n,z,0,prec); }
     546             : GEN
     547         658 : nbessel(GEN n, GEN z, long prec) { return kbesselintern(n,z,1,prec); }
     548             : /* J + iN */
     549             : GEN
     550         252 : hbessel1(GEN n, GEN z, long prec)
     551             : {
     552         252 :   pari_sp av = avma;
     553         252 :   GEN J = jbessel(n,z,prec);
     554         196 :   GEN N = nbessel(n,z,prec);
     555         182 :   return gerepileupto(av, gadd(J, mulcxI(N)));
     556             : }
     557             : /* J - iN */
     558             : GEN
     559         252 : hbessel2(GEN n, GEN z, long prec)
     560             : {
     561         252 :   pari_sp av = avma;
     562         252 :   GEN J = jbessel(n,z,prec);
     563         196 :   GEN N = nbessel(n,z,prec);
     564         182 :   return gerepileupto(av, gadd(J, mulcxmI(N)));
     565             : }
     566             : 
     567             : /***********************************************************************/
     568             : /*                                                                    **/
     569             : /**                    FONCTION U(a,b,z) GENERALE                     **/
     570             : /**                         ET CAS PARTICULIERS                       **/
     571             : /**                                                                   **/
     572             : /***********************************************************************/
     573             : /* Assume gx > 0 and a,b complex */
     574             : /* This might one day be extended to handle complex gx */
     575             : /* see Temme, N. M. "The numerical computation of the confluent        */
     576             : /* hypergeometric function U(a,b,z)" in Numer. Math. 41 (1983),        */
     577             : /* no. 1, 63--82.                                                      */
     578             : GEN
     579          14 : hyperu(GEN a, GEN b, GEN gx, long prec)
     580             : {
     581          14 :   GEN S, P, T, x, p1, zf, u, a1, mb = gneg(b);
     582          14 :   const int ex = iscomplex(a) || iscomplex(b);
     583          14 :   long k, n, l = (typ(gx)==t_REAL)? realprec(gx): prec, l1 = l+EXTRAPRECWORD;
     584          14 :   GEN y = ex? cgetc(l): cgetr(l);
     585          14 :   pari_sp av = avma;
     586             : 
     587          14 :   if (gsigne(gx) <= 0) pari_err_IMPL("non-positive third argument in hyperu");
     588          14 :   x = gtofp(gx, l);
     589          14 :   a1 = gaddsg(1, gadd(a,mb)); P = gmul(a1, a);
     590          14 :   n = (long)(prec2nbits_mul(l, LOG2) + M_PI*sqrt(dblmodulus(P)));
     591          14 :   S = gadd(a1, a);
     592          14 :   if (cmprs(x,n) < 0)
     593             :   {
     594           7 :     GEN q = stor(n, l1), s = gen_1, t = gen_0, v, c, e, f;
     595             :     pari_sp av1, av2;
     596             : 
     597           7 :     if (ex) { u=cgetc(l1); v=cgetc(l1); e=cgetc(l1); f=cgetc(l1); }
     598           7 :     else    { u=cgetr(l1); v=cgetr(l1); e=cgetr(l1); f=cgetr(l1); }
     599           7 :     av1 = avma;
     600           7 :     zf = gpow(stoi(n),gneg_i(a),l1);
     601           7 :     T = gadd(gadd(P, gmulsg(n-1, S)), sqrs(n-1));
     602         637 :     for (k=n-1; k>=0; k--)
     603             :     {
     604             :       /* T = (a+k)*(a1+k) = a*a1 + k(a+a1) + k^2 = previous(T) - S - 2k + 1 */
     605         637 :       p1 = gdiv(T, mulss(-n, k+1));
     606         637 :       s = gaddgs(gmul(p1,s), 1);
     607         637 :       t = gadd(  gmul(p1,t), gaddgs(a,k));
     608         637 :       if (!k) break;
     609         630 :       T = gsubgs(gsub(T, S), 2*k-1);
     610             :     }
     611           7 :     gmulz(zf, s, u);
     612           7 :     gmulz(zf, gdivgs(t,-n), v);
     613         133 :     for(;; avma = av1)
     614             :     {
     615         140 :       GEN d = real_1(l1), p3 = gadd(q,mb);
     616         140 :       c = divur(5,q); if (expo(c)>= -1) c = real2n(-1, l1);
     617         140 :       p1 = subsr(1,divrr(x,q)); if (cmprr(c,p1)>0) c = p1;
     618         140 :       togglesign(c);
     619         140 :       gaffect(u,e);
     620         140 :       gaffect(v,f); av2 = avma;
     621       11659 :       for(k=1;;k++, avma = av2)
     622             :       {
     623       11659 :         GEN w = gadd(gmul(gaddgs(a,k-1),u), gmul(gaddgs(p3,1-k),v));
     624       11659 :         gmulz(divru(q,k),v, u);
     625       11659 :         gaffect(gdivgs(w,k), v);
     626       11659 :         mulrrz(d,c,d);
     627       11659 :         gaddz(e,gmul(d,u),e); p1=gmul(d,v);
     628       11659 :         gaddz(f,p1,f);
     629       11659 :         if (gequal0(p1) || gexpo(p1) - gexpo(f) <= 1-prec2nbits(precision(p1)))
     630             :           break;
     631       11519 :       }
     632         140 :       swap(e, u);
     633         140 :       swap(f, v);
     634         140 :       affrr(mulrr(q, addrs(c,1)), q);
     635         140 :       if (expo(subrr(q,x)) - expo(x) <= 1-prec2nbits(l)) break;
     636         133 :     }
     637             :   }
     638             :   else
     639             :   {
     640           7 :     GEN zz = invr(x), s = gen_1;
     641           7 :     togglesign(zz); /* -1/x */
     642           7 :     zf = gpow(x,gneg_i(a),l1);
     643           7 :     T = gadd(gadd(P, gmulsg(n-1, S)), sqrs(n-1));
     644         644 :     for (k=n-1; k>=0; k--)
     645             :     {
     646         644 :       p1 = gmul(T,divru(zz,k+1));
     647         644 :       s = gaddsg(1, gmul(p1,s));
     648         644 :       if (!k) break;
     649         637 :       T = gsubgs(gsub(T, S), 2*k-1);
     650             :     }
     651           7 :     u = gmul(s,zf);
     652             :   }
     653          14 :   gaffect(u,y); avma = av; return y;
     654             : }
     655             : 
     656             : /* incgam(0, x, prec) = eint1(x); typ(x) = t_REAL, x > 0 */
     657             : static GEN
     658        6041 : incgam_0(GEN x, GEN expx)
     659             : {
     660             :   pari_sp av;
     661        6041 :   long l = realprec(x), n, i;
     662        6041 :   double mx = rtodbl(x), L = prec2nbits_mul(l,LOG2);
     663             :   GEN z;
     664             : 
     665        6041 :   if (!mx) pari_err_DOMAIN("eint1", "x","=",gen_0, x);
     666        6034 :   if (mx > L)
     667             :   {
     668         280 :     double m = (L + mx)/4;
     669         280 :     n = (long)(1+m*m/mx);
     670         280 :     av = avma;
     671         280 :     z = divsr(-n, addsr(n<<1,x));
     672        7254 :     for (i=n-1; i >= 1; i--)
     673             :     {
     674        6974 :       z = divsr(-i, addrr(addsr(i<<1,x), mulur(i,z))); /* -1 / (2 + z + x/i) */
     675        6974 :       if ((i & 0x1ff) == 0) z = gerepileuptoleaf(av, z);
     676             :     }
     677         280 :     return divrr(addrr(real_1(l),z), mulrr(expx? expx: mpexp(x), x));
     678             :   }
     679             :   else
     680             :   {
     681        5754 :     long prec = nbits2prec(prec2nbits(l) + (mx+log(mx))/LOG2 + 10);
     682        5754 :     GEN S, t, H, run = real_1(prec);
     683        5754 :     n = -prec2nbits(prec);
     684        5754 :     x = rtor(x, prec);
     685        5754 :     av = avma;
     686        5754 :     S = z = t = H = run;
     687      504876 :     for (i = 2; expo(t) - expo(S) >= n; i++)
     688             :     {
     689      499122 :       H = addrr(H, divru(run,i)); /* H = sum_{k<=i} 1/k */
     690      499122 :       z = divru(mulrr(x,z), i);   /* z = x^(i-1)/i! */
     691      499122 :       t = mulrr(z, H); S = addrr(S, t);
     692      499122 :       if ((i & 0x1ff) == 0) gerepileall(av, 4, &z,&t,&S,&H);
     693             :     }
     694        5754 :     return subrr(mulrr(x, divrr(S,expx? expx: mpexp(x))),
     695             :                  addrr(mplog(x), mpeuler(prec)));
     696             :   }
     697             : }
     698             : 
     699             : /* real(z*log(z)-z), z = x+iy */
     700             : static double
     701        7651 : mygamma(double x, double y)
     702             : {
     703        7651 :   if (x == 0.) return -(M_PI/2)*fabs(y);
     704        7651 :   return (x/2)*log(x*x+y*y)-x-y*atan(y/x);
     705             : }
     706             : 
     707             : /* x^s exp(-x) */
     708             : static GEN
     709       10101 : expmx_xs(GEN s, GEN x, GEN logx, long prec)
     710             : {
     711             :   GEN z;
     712       10101 :   long ts = typ(s);
     713       10101 :   if (ts == t_INT || (ts == t_FRAC && absequaliu(gel(s,2), 2)))
     714        5348 :     z = gmul(gexp(gneg(x), prec), gpow(x, s, prec));
     715             :   else
     716        4753 :     z = gexp(gsub(gmul(s, logx? logx: glog(x,prec+EXTRAPREC)), x), prec);
     717       10101 :   return z;
     718             : }
     719             : 
     720             : /* Not yet: doesn't work at low accuracy
     721             : #define INCGAM_CF
     722             : */
     723             : 
     724             : #ifdef INCGAM_CF
     725             : /* Is s very close to a non-positive integer ? */
     726             : static int
     727             : isgammapole(GEN s, long bitprec)
     728             : {
     729             :   pari_sp av = avma;
     730             :   GEN t = imag_i(s);
     731             :   long e, b = bitprec - 10;
     732             : 
     733             :   if (gexpo(t) > - b) return 0;
     734             :   s = real_i(s);
     735             :   if (gsigne(s) > 0 && gexpo(s) > -b) return 0;
     736             :   (void)grndtoi(s, &e); avma = av; return (e < -b);
     737             : }
     738             : 
     739             : /* incgam using the continued fraction. x a t_REAL or t_COMPLEX, mx ~ |x|.
     740             :  * Assume precision(s), precision(x) >= prec */
     741             : static GEN
     742             : incgam_cf(GEN s, GEN x, double mx, long prec)
     743             : {
     744             :   GEN ms, y, S;
     745             :   long n, i, j, LS, bitprec = prec2nbits(prec);
     746             :   double rs, is, m;
     747             : 
     748             :   if (typ(s) == t_COMPLEX)
     749             :   {
     750             :     rs = gtodouble(gel(s,1));
     751             :     is = gtodouble(gel(s,2));
     752             :   }
     753             :   else
     754             :   {
     755             :     rs = gtodouble(s);
     756             :     is = 0.;
     757             :   }
     758             :   if (isgammapole(s, bitprec)) LS = 0;
     759             :   else
     760             :   {
     761             :     double bit,  LGS = mygamma(rs,is);
     762             :     LS = LGS <= 0 ? 0: ceil(LGS);
     763             :     bit = (LGS - (rs-1)*log(mx) + mx)/LOG2;
     764             :     if (bit > 0)
     765             :     {
     766             :       prec += nbits2extraprec((long)bit);
     767             :       x = gtofp(x, prec);
     768             :       if (isinexactreal(s)) s = gtofp(s, prec);
     769             :     }
     770             :   }
     771             :   /* |ln(2*gamma(s)*sin(s*Pi))| <= ln(2) + |lngamma(s)| + |Im(s)*Pi|*/
     772             :   m = bitprec*LOG2 + LS + LOG2 + fabs(is)*M_PI + mx;
     773             :   if (rs < 1) m += (1 - rs)*log(mx);
     774             :   m /= 4;
     775             :   n = (long)(1 + m*m/mx);
     776             :   y = expmx_xs(gsubgs(s,1), x, NULL, prec);
     777             :   if (rs >= 0 && bitprec >= 512)
     778             :   {
     779             :     GEN A = cgetg(n+1, t_VEC), B = cgetg(n+1, t_VEC);
     780             :     ms = gsubsg(1, s);
     781             :     for (j = 1; j <= n; ++j)
     782             :     {
     783             :       gel(A,j) = ms;
     784             :       gel(B,j) = gmulsg(j, gsubgs(s,j));
     785             :       ms = gaddgs(ms, 2);
     786             :     }
     787             :     S = contfraceval_inv(mkvec2(A,B), x, -1);
     788             :   }
     789             :   else
     790             :   {
     791             :     GEN x_s = gsub(x, s);
     792             :     pari_sp av2 = avma;
     793             :     S = gdiv(gsubgs(s,n), gaddgs(x_s,n<<1));
     794             :     for (i=n-1; i >= 1; i--)
     795             :     {
     796             :       S = gdiv(gsubgs(s,i), gadd(gaddgs(x_s,i<<1),gmulsg(i,S)));
     797             :       if (gc_needed(av2,3))
     798             :       {
     799             :         if(DEBUGMEM>1) pari_warn(warnmem,"incgam_cf");
     800             :         S = gerepileupto(av2, S);
     801             :       }
     802             :     }
     803             :     S = gaddgs(S,1);
     804             :   }
     805             :   return gmul(y, S);
     806             : }
     807             : #endif
     808             : 
     809             : static double
     810        5579 : findextraincgam(GEN s, GEN x)
     811             : {
     812        5579 :   double sig = gtodouble(real_i(s)), t = gtodouble(imag_i(s));
     813        5579 :   double xr = gtodouble(real_i(x)), xi = gtodouble(imag_i(x));
     814        5579 :   double exd = 0., Nx = xr*xr + xi*xi, D = Nx - t*t;
     815             :   long n;
     816             : 
     817        5579 :   if (xr < 0)
     818             :   {
     819         819 :     long ex = gexpo(x);
     820         819 :     if (ex > 0 && ex > gexpo(s)) exd = sqrt(Nx)*log(Nx)/2; /* |x| log |x| */
     821             :   }
     822        5579 :   if (D <= 0.) return exd;
     823        4137 :   n = (long)(sqrt(D)-sig);
     824        4137 :   if (n <= 0) return exd;
     825        1484 :   return maxdd(exd, (n*log(Nx)/2 - mygamma(sig+n, t) + mygamma(sig, t)) / LOG2);
     826             : }
     827             : 
     828             : /* use exp(-x) * (x^s/s) * sum_{k >= 0} x^k / prod(i=1, k, s+i) */
     829             : static GEN
     830        5586 : incgamc_i(GEN s, GEN x, long *ptexd, long prec)
     831             : {
     832             :   GEN S, t, y;
     833             :   long l, n, i, exd;
     834        5586 :   pari_sp av = avma, av2;
     835             : 
     836        5586 :   if (gequal0(x))
     837             :   {
     838           7 :     if (ptexd) *ptexd = 0.;
     839           7 :     return gtofp(x, prec);
     840             :   }
     841        5579 :   l = precision(x);
     842        5579 :   if (!l) l = prec;
     843        5579 :   n = -prec2nbits(l)-1;
     844        5579 :   exd = (long)findextraincgam(s, x);
     845        5579 :   if (ptexd) *ptexd = exd;
     846        5579 :   if (exd > 0)
     847             :   {
     848        1372 :     long p = l + nbits2extraprec(exd);
     849        1372 :     x = gtofp(x, p);
     850        1372 :     if (isinexactreal(s)) s = gtofp(s, p);
     851             :   }
     852        4207 :   else x = gtofp(x, l+EXTRAPREC);
     853        5579 :   av2 = avma;
     854        5579 :   S = gdiv(x, gaddsg(1,s));
     855        5579 :   t = gaddsg(1, S);
     856      749012 :   for (i=2; gexpo(S) >= n; i++)
     857             :   {
     858      743433 :     S = gdiv(gmul(x,S), gaddsg(i,s)); /* x^i / ((s+1)...(s+i)) */
     859      743433 :     t = gadd(S,t);
     860      743433 :     if (gc_needed(av2,3))
     861             :     {
     862           0 :       if(DEBUGMEM>1) pari_warn(warnmem,"incgamc");
     863           0 :       gerepileall(av2, 2, &S, &t);
     864             :     }
     865             :   }
     866        5579 :   y = expmx_xs(s, x, NULL, prec);
     867        5579 :   return gerepileupto(av, gmul(gdiv(y,s), t));
     868             : }
     869             : 
     870             : GEN
     871        1414 : incgamc(GEN s, GEN x, long prec)
     872        1414 : { return incgamc_i(s, x, NULL, prec); }
     873             : 
     874             : /* incgamma using asymptotic expansion:
     875             :  *   exp(-x)x^(s-1)(1 + (s-1)/x + (s-1)(s-2)/x^2 + ...) */
     876             : static GEN
     877        2821 : incgam_asymp(GEN s, GEN x, long prec)
     878             : {
     879        2821 :   pari_sp av = avma, av2;
     880             :   GEN S, q, cox, invx;
     881        2821 :   long oldeq = LONG_MAX, eq, esx, j;
     882        2821 :   int flint = (typ(s) == t_INT && signe(s) > 0);
     883             : 
     884        2821 :   x = gtofp(x,prec+EXTRAPREC);
     885        2821 :   invx = ginv(x);
     886        2821 :   esx = -prec2nbits(prec);
     887        2821 :   av2 = avma;
     888        2821 :   q = gmul(gsubgs(s,1), invx);
     889        2821 :   S = gaddgs(q, 1);
     890      129808 :   for (j = 2;; j++)
     891             :   {
     892      129808 :     eq = gexpo(q); if (eq < esx) break;
     893      127106 :     if (!flint && (j & 3) == 0)
     894             :     { /* guard against divergence */
     895       17220 :       if (eq > oldeq) { avma = av; return NULL; } /* regressing, abort */
     896       17101 :       oldeq = eq;
     897             :     }
     898      126987 :     q = gmul(q, gmul(gsubgs(s,j), invx));
     899      126987 :     S = gadd(S, q);
     900      126987 :     if (gc_needed(av2, 1)) gerepileall(av2, 2, &S, &q);
     901      126987 :   }
     902        2702 :   if (DEBUGLEVEL > 2) err_printf("incgam: using asymp\n");
     903        2702 :   cox = expmx_xs(gsubgs(s,1), x, NULL, prec);
     904        2702 :   return gerepileupto(av, gmul(cox, S));
     905             : }
     906             : 
     907             : /* gasx = incgam(s-n,x). Compute incgam(s,x)
     908             :  * = (s-1)(s-2)...(s-n)gasx + exp(-x)x^(s-1) *
     909             :  *   (1 + (s-1)/x + ... + (s-1)(s-2)...(s-n+1)/x^(n-1)) */
     910             : static GEN
     911         546 : incgam_asymp_partial(GEN s, GEN x, GEN gasx, long n, long prec)
     912             : {
     913             :   pari_sp av;
     914         546 :   GEN S, q, cox, invx, s1 = gsubgs(s, 1), sprod;
     915             :   long j;
     916         546 :   cox = expmx_xs(s1, x, NULL, prec);
     917         546 :   if (n == 1) return gadd(cox, gmul(s1, gasx));
     918         546 :   invx = ginv(x);
     919         546 :   av = avma;
     920         546 :   q = gmul(s1, invx);
     921         546 :   S = gaddgs(q, 1);
     922       52164 :   for (j = 2; j < n; j++)
     923             :   {
     924       51618 :     q = gmul(q, gmul(gsubgs(s, j), invx));
     925       51618 :     S = gadd(S, q);
     926       51618 :     if (gc_needed(av, 2)) gerepileall(av, 2, &S, &q);
     927             :   }
     928         546 :   sprod = gmul(gmul(q, gpowgs(x, n-1)), gsubgs(s, n));
     929         546 :   return gadd(gmul(cox, S), gmul(sprod, gasx));
     930             : }
     931             : 
     932             : /* assume s != 0 */
     933             : static GEN
     934        2394 : incgamspec(GEN s, GEN x, GEN g, long prec)
     935             : {
     936        2394 :   GEN q, S, cox = gen_0, P, sk, S1, S2, S3, F3, logx, mx;
     937        2394 :   long n, esk, k = itos(ground(gneg(real_i(s)))), E;
     938             : 
     939        2394 :   if (k && gexpo(x) > 0)
     940             :   {
     941         238 :     GEN xk = gdivgs(x, k);
     942         238 :     long bitprec = prec2nbits(prec);
     943         238 :     double mx = (gexpo(xk) > bitprec)? bitprec*LOG2: log(dblmodulus(xk));
     944         238 :     prec += nbits2extraprec((long)k*(mx + 1)/LOG2);
     945         238 :     if (isinexactreal(s)) s = gtofp(s, prec);
     946             :   }
     947        2394 :   x = gtofp(x, maxss(precision(x), prec) + EXTRAPREC);
     948        2394 :   sk = gaddgs(s, k); /* |Re(sk)| <= 1/2 */
     949        2394 :   logx = glog(x, prec);
     950        2394 :   mx = gneg(x);
     951        2394 :   if (k == 0) { S = gen_0; P = gen_1; }
     952             :   else
     953             :   {
     954             :     long j;
     955         847 :     q = ginv(s); S = q; P = s;
     956        9926 :     for (j = 1; j < k; j++)
     957             :     {
     958        9079 :       GEN sj = gaddgs(s, j);
     959        9079 :       q = gmul(q, gdiv(x, sj));
     960        9079 :       S = gadd(S, q);
     961        9079 :       P = gmul(P, sj);
     962             :     }
     963         847 :     cox = expmx_xs(s, x, logx, prec); /* x^s exp(-x) */
     964         847 :     S = gmul(S, gneg(cox));
     965             :   }
     966        2394 :   if (k && gequal0(sk))
     967         175 :     return gadd(S, gdiv(eint1(x, prec), P));
     968        2219 :   esk = gexpo(sk);
     969        2219 :   if (esk > -7)
     970             :   {
     971        1008 :     GEN a, b, PG = gmul(sk, P);
     972        1008 :     if (g) g = gmul(g, PG);
     973        1008 :     a = incgam0(gaddgs(sk,1), x, g, prec);
     974        1008 :     if (k == 0) cox = expmx_xs(s, x, logx, prec);
     975        1008 :     b = gmul(gpowgs(x, k), cox);
     976        1008 :     return gadd(S, gdiv(gsub(a, b), PG));
     977             :   }
     978        1211 :   E = prec2nbits(prec) + 1;
     979        1211 :   if (gexpo(x) > 0)
     980             :   {
     981         420 :     long X = (long)(dblmodulus(x)/LOG2);
     982         420 :     prec += 2*nbits2extraprec(X);
     983         420 :     x = gtofp(x, prec); mx = gneg(x);
     984         420 :     logx = glog(x, prec); sk = gtofp(sk, prec);
     985         420 :     E += X;
     986             :   }
     987        1211 :   if (isinexactreal(sk)) sk = gtofp(sk, prec+EXTRAPREC);
     988             :   /* |sk| < 2^-7 is small, guard against cancellation */
     989        1211 :   F3 = gexpm1(gmul(sk, logx), prec);
     990             :   /* ( gamma(1+sk) - exp(sk log(x))) ) / sk */
     991        1211 :   S1 = gdiv(gsub(ggamma1m1(sk, prec+EXTRAPREC), F3), sk);
     992        1211 :   q = x; S3 = gdiv(x, gaddsg(1,sk));
     993      255523 :   for (n = 2; gexpo(q) - gexpo(S3) > -E; ++n)
     994             :   {
     995      254312 :     q = gmul(q, gdivgs(mx, n));
     996      254312 :     S3 = gadd(S3, gdiv(q, gaddsg(n, sk)));
     997             :   }
     998        1211 :   S2 = gadd(gadd(S1, S3), gmul(F3, S3));
     999        1211 :   return gadd(S, gdiv(S2, P));
    1000             : }
    1001             : 
    1002             : #if 0
    1003             : static long
    1004             : precision2(GEN x, GEN y)
    1005             : {
    1006             :   long px = precision(x), py = precision(y);
    1007             :   if (!px) return py;
    1008             :   if (!py) return px;
    1009             :   return minss(px, py);
    1010             : }
    1011             : #endif
    1012             : 
    1013             : /* return |x| */
    1014             : double
    1015     4313506 : dblmodulus(GEN x)
    1016             : {
    1017     4313506 :   if (typ(x) == t_COMPLEX)
    1018             :   {
    1019      356860 :     double a = gtodouble(gel(x,1));
    1020      356860 :     double b = gtodouble(gel(x,2));
    1021      356860 :     return sqrt(a*a + b*b);
    1022             :   }
    1023             :   else
    1024     3956646 :     return fabs(gtodouble(x));
    1025             : }
    1026             : 
    1027             : /* Driver routine. If g != NULL, assume that g=gamma(s,prec). */
    1028             : GEN
    1029       11515 : incgam0(GEN s, GEN x, GEN g, long prec)
    1030             : {
    1031             :   pari_sp av;
    1032             :   long E, l, ex;
    1033             :   double mx;
    1034             :   GEN z, rs, is;
    1035             : 
    1036       11515 :   if (gequal0(x)) return g? gcopy(g): ggamma(s,prec);
    1037       11515 :   if (gequal0(s)) return eint1(x, prec);
    1038        9709 :   av = avma;
    1039        9709 :   l = precision(s);
    1040        9709 :   if (!l) l = prec;
    1041        9709 :   E = prec2nbits(l) + 1;
    1042             :   /* avoid overflow in dblmodulus */
    1043        9709 :   ex = gexpo(x);
    1044        9709 :   if (ex > E) mx = E; else mx = dblmodulus(x);
    1045             :   /* use asymptotic expansion */
    1046        9709 :   if (4*mx > 3*E || (typ(s) == t_INT && signe(s) > 0 && ex >= expi(s)))
    1047             :   {
    1048        2716 :     z = incgam_asymp(s, x, l);
    1049        2716 :     if (z) return z;
    1050             :   }
    1051        7112 :   rs = real_i(s);
    1052        7112 :   is = imag_i(s);
    1053             : #ifdef INCGAM_CF
    1054             :   /* Can one use continued fraction ? */
    1055             :   if (gequal0(is) && gequal0(imag_i(x)) && gsigne(x) > 0)
    1056             :   {
    1057             :     double sd = gtodouble(rs), LB, UB;
    1058             :     double xd = gtodouble(real_i(x));
    1059             :     if (sd > 0) {
    1060             :       LB = 15 + 0.1205*E;
    1061             :       UB = 5 + 0.752*E;
    1062             :     } else {
    1063             :       LB = -6 + 0.1205*E;
    1064             :       UB = 5 + 0.752*E + fabs(sd)/54.;
    1065             :     }
    1066             :     if (xd >= LB && xd <= UB)
    1067             :     {
    1068             :       if (DEBUGLEVEL > 2) err_printf("incgam: using continued fraction\n");
    1069             :       return gerepileupto(av, incgam_cf(s, x, xd, prec));
    1070             :     }
    1071             :   }
    1072             : #endif
    1073        7112 :   if (gsigne(rs) > 0 && gexpo(rs) >= -1)
    1074             :   { /* use complementary incomplete gamma */
    1075        4718 :     long n, egs, exd, precg, es = gexpo(s);
    1076        4718 :     if (es < 0) {
    1077         574 :       l += nbits2extraprec(-es) + 1;
    1078         574 :       x = gtofp(x, l);
    1079         574 :       if (isinexactreal(s)) s = gtofp(s, l);
    1080             :     }
    1081        4718 :     n = itos(gceil(rs));
    1082        4718 :     if (n > 100)
    1083             :     {
    1084             :       GEN gasx;
    1085         546 :       n -= 100;
    1086         546 :       if (es > 0)
    1087             :       {
    1088         546 :         es = mygamma(gtodouble(rs) - n, gtodouble(is)) / LOG2;
    1089         546 :         if (es > 0)
    1090             :         {
    1091         546 :           l += nbits2extraprec(es);
    1092         546 :           x = gtofp(x, l);
    1093         546 :           if (isinexactreal(s)) s = gtofp(s, l);
    1094             :         }
    1095             :       }
    1096         546 :       gasx = incgam0(gsubgs(s, n), x, NULL, prec);
    1097         546 :       return gerepileupto(av, incgam_asymp_partial(s, x, gasx, n, prec));
    1098             :     }
    1099        4172 :     if (DEBUGLEVEL > 2) err_printf("incgam: using power series 1\n");
    1100             :     /* egs ~ expo(gamma(s)) */
    1101        4172 :     precg = g? precision(g): 0;
    1102        4172 :     egs = g? gexpo(g): (long)(mygamma(gtodouble(rs), gtodouble(is)) / LOG2);
    1103        4172 :     if (egs > 0) {
    1104        1946 :       l += nbits2extraprec(egs) + 1;
    1105        1946 :       x = gtofp(x, l);
    1106        1946 :       if (isinexactreal(s)) s = gtofp(s, l);
    1107        1946 :       if (precg < l) g = NULL;
    1108             :     }
    1109        4172 :     z = incgamc_i(s, x, &exd, l);
    1110        4172 :     if (exd > 0)
    1111             :     {
    1112         896 :       l += nbits2extraprec(exd);
    1113         896 :       if (isinexactreal(s)) s = gtofp(s, l);
    1114         896 :       if (precg < l) g = NULL;
    1115             :     }
    1116             :     else
    1117             :     { /* gamma(s) negligible ? Compute to lower accuracy */
    1118        3276 :       long e = gexpo(z) - egs;
    1119        3276 :       if (e > 3)
    1120             :       {
    1121         420 :         E -= e;
    1122         420 :         if (E <= 0) g = gen_0; else if (!g) g = ggamma(s, nbits2prec(E));
    1123             :       }
    1124             :     }
    1125             :     /* worry about possible cancellation */
    1126        4172 :     if (!g) g = ggamma(s, maxss(l,precision(z)));
    1127        4172 :     return gerepileupto(av, gsub(g,z));
    1128             :   }
    1129        2394 :   if (DEBUGLEVEL > 2) err_printf("incgam: using power series 2\n");
    1130        2394 :   return gerepileupto(av, incgamspec(s, x, g, l));
    1131             : }
    1132             : 
    1133             : GEN
    1134        1106 : incgam(GEN s, GEN x, long prec) { return incgam0(s, x, NULL, prec); }
    1135             : 
    1136             : /* x >= 0 a t_REAL */
    1137             : GEN
    1138        1974 : mpeint1(GEN x, GEN expx)
    1139             : {
    1140        1974 :   GEN z = cgetr(realprec(x));
    1141        1974 :   pari_sp av = avma;
    1142        1974 :   affrr(incgam_0(x, expx), z);
    1143        1967 :   avma = av; return z;
    1144             : }
    1145             : 
    1146             : static GEN
    1147         357 : cxeint1(GEN x, long prec)
    1148             : {
    1149         357 :   pari_sp av = avma, av2;
    1150             :   GEN q, S3;
    1151             :   GEN run, z, H;
    1152         357 :   long n, E = prec2nbits(prec) + 1, ex = gexpo(x);
    1153             : 
    1154         357 :   if ((ex > E || 4*dblmodulus(x) > 3*E)
    1155         105 :       && (z = incgam_asymp(gen_0, x, prec))) return z;
    1156         252 :   if (ex > 0)
    1157             :   { /* take cancellation into account, log2(\sum |x|^n / n!) = |x| / log(2) */
    1158          42 :     double dbx = dblmodulus(x);
    1159          42 :     long X = (long)((dbx + log(dbx))/LOG2 + 10);
    1160          42 :     prec += nbits2extraprec(X);
    1161          42 :     x = gtofp(x, prec); E += X;
    1162             :   }
    1163         252 :   if (DEBUGLEVEL > 2) err_printf("eint1: using power series\n");
    1164         252 :   run = real_1(prec);
    1165         252 :   av2 = avma;
    1166         252 :   S3 = z = q = H = run;
    1167       48384 :   for (n = 2; gexpo(q) - gexpo(S3) >= -E; n++)
    1168             :   {
    1169       48132 :     H = addrr(H, divru(run, n)); /* H = sum_{k<=n} 1/k */
    1170       48132 :     z = gdivgs(gmul(x,z), n);   /* z = x^(n-1)/n! */
    1171       48132 :     q = gmul(z, H); S3 = gadd(S3, q);
    1172       48132 :     if ((n & 0x1ff) == 0) gerepileall(av2, 4, &z, &q, &S3, &H);
    1173             :   }
    1174         252 :   return gerepileupto(av, gsub(gmul(x, gdiv(S3, gexp(x, prec))), gadd(glog(x, prec), mpeuler(prec))));
    1175             : }
    1176             : 
    1177             : GEN
    1178        2471 : eint1(GEN x, long prec)
    1179             : {
    1180             :   long l, n, i;
    1181             :   pari_sp av;
    1182             :   GEN p1, t, S, y, res;
    1183             : 
    1184        2471 :   switch(typ(x))
    1185             :   {
    1186         357 :     case t_COMPLEX: return cxeint1(x, prec);
    1187        1729 :     case t_REAL: break;
    1188         385 :     default: x = gtofp(x, prec);
    1189             :   }
    1190        2114 :   if (signe(x) >= 0) return mpeint1(x,NULL);
    1191             :   /* rewritten from code contributed by Manfred Radimersky */
    1192         140 :   res = cgetg(3, t_COMPLEX);
    1193         140 :   av = avma;
    1194         140 :   l  = realprec(x);
    1195         140 :   n  = prec2nbits(l);
    1196         140 :   y  = rtor(x, l + EXTRAPREC);
    1197         140 :   setsigne(y,1);
    1198         140 :   if (cmprs(y, (3*n)/4) < 0) {
    1199         126 :     p1 = t = S = y;
    1200       24248 :     for (i = 2; expo(t) - expo(S) >= -n; i++) {
    1201       24122 :       p1 = mulrr(y, divru(p1, i)); /* (-x)^i/i! */
    1202       24122 :       t = divru(p1, i);
    1203       24122 :       S = addrr(S, t);
    1204             :     }
    1205         126 :     y  = addrr(S, addrr(logr_abs(x), mpeuler(l)));
    1206             :   } else { /* ~incgam_asymp: asymptotic expansion */
    1207          14 :     p1 = t = invr(y);
    1208          14 :     S = addrs(t, 1);
    1209         560 :     for (i = 2; expo(t) >= -n; i++) {
    1210         546 :       t = mulrr(p1, mulru(t, i));
    1211         546 :       S = addrr(S, t);
    1212             :     }
    1213          14 :     y  = mulrr(S, mulrr(p1, mpexp(y)));
    1214             :   }
    1215         140 :   gel(res, 1) = gerepileuptoleaf(av, negr(y));
    1216         140 :   y = mppi(prec); setsigne(y, -1);
    1217         140 :   gel(res, 2) = y; return res;
    1218             : }
    1219             : 
    1220             : GEN
    1221          49 : veceint1(GEN C, GEN nmax, long prec)
    1222             : {
    1223          49 :   if (!nmax) return eint1(C,prec);
    1224           7 :   if (typ(nmax) != t_INT) pari_err_TYPE("veceint1",nmax);
    1225           7 :   if (typ(C) != t_REAL) {
    1226           7 :     C = gtofp(C, prec);
    1227           7 :     if (typ(C) != t_REAL) pari_err_TYPE("veceint1",C);
    1228             :   }
    1229           7 :   if (signe(C) <= 0) pari_err_DOMAIN("veceint1", "argument", "<=", gen_0,C);
    1230           7 :   return mpveceint1(C, NULL, itos(nmax));
    1231             : }
    1232             : 
    1233             : /* j > 0, a t_REAL. Return sum_{m >= 0} a^m / j(j+1)...(j+m)).
    1234             :  * Stop when expo(summand) < E; note that s(j-1) = (a s(j) + 1) / (j-1). */
    1235             : static GEN
    1236         252 : mp_sum_j(GEN a, long j, long E, long prec)
    1237             : {
    1238         252 :   pari_sp av = avma;
    1239         252 :   GEN q = divru(real_1(prec), j), s = q;
    1240             :   long m;
    1241        4621 :   for (m = 0;; m++)
    1242             :   {
    1243        4621 :     if (expo(q) < E) break;
    1244        4369 :     q = mulrr(q, divru(a, m+j));
    1245        4369 :     s = addrr(s, q);
    1246        4369 :   }
    1247         252 :   return gerepileuptoleaf(av, s);
    1248             : }
    1249             : /* Return the s_a(j), j <= J */
    1250             : static GEN
    1251         252 : sum_jall(GEN a, long J, long prec)
    1252             : {
    1253         252 :   GEN s = cgetg(J+1, t_VEC);
    1254         252 :   long j, E = -prec2nbits(prec) - 5;
    1255         252 :   gel(s, J) = mp_sum_j(a, J, E, prec);
    1256       10492 :   for (j = J-1; j; j--)
    1257       10240 :     gel(s,j) = divru(addrs(mulrr(a, gel(s,j+1)), 1), j);
    1258         252 :   return s;
    1259             : }
    1260             : 
    1261             : /* T a dense t_POL with t_REAL coeffs. Return T(n) [faster than poleval] */
    1262             : static GEN
    1263      487288 : rX_s_eval(GEN T, long n)
    1264             : {
    1265      487288 :   long i = lg(T)-1;
    1266      487288 :   GEN c = gel(T,i);
    1267      487288 :   for (i--; i>=2; i--) c = gadd(mulrs(c,n),gel(T,i));
    1268      487288 :   return c;
    1269             : }
    1270             : 
    1271             : /* C>0 t_REAL, eC = exp(C). Return eint1(n*C) for 1<=n<=N. Absolute accuracy */
    1272             : GEN
    1273         259 : mpveceint1(GEN C, GEN eC, long N)
    1274             : {
    1275         259 :   const long prec = realprec(C);
    1276         259 :   long Nmin = 15; /* >= 1. E.g. between 10 and 30, but little effect */
    1277         259 :   GEN en, v, w = cgetg(N+1, t_VEC);
    1278             :   pari_sp av0;
    1279             :   double DL;
    1280             :   long n, j, jmax, jmin;
    1281         259 :   if (!N) return w;
    1282         259 :   for (n = 1; n <= N; n++) gel(w,n) = cgetr(prec);
    1283         259 :   av0 = avma;
    1284         259 :   if (N < Nmin) Nmin = N;
    1285         259 :   if (!eC) eC = mpexp(C);
    1286         259 :   en = eC; affrr(incgam_0(C, en), gel(w,1));
    1287        3815 :   for (n = 2; n <= Nmin; n++)
    1288             :   {
    1289             :     pari_sp av2;
    1290        3556 :     en = mulrr(en,eC); /* exp(n C) */
    1291        3556 :     av2 = avma;
    1292        3556 :     affrr(incgam_0(mulru(C,n), en), gel(w,n));
    1293        3556 :     avma = av2;
    1294             :   }
    1295         259 :   if (Nmin == N) { avma = av0; return w; }
    1296             : 
    1297         252 :   DL = prec2nbits_mul(prec, LOG2) + 5;
    1298         252 :   jmin = ceil(DL/log((double)N)) + 1;
    1299         252 :   jmax = ceil(DL/log((double)Nmin)) + 1;
    1300         252 :   v = sum_jall(C, jmax, prec);
    1301         252 :   en = powrs(eC, -N); /* exp(-N C) */
    1302         252 :   affrr(incgam_0(mulru(C,N), invr(en)), gel(w,N));
    1303        6654 :   for (j = jmin, n = N-1; j <= jmax; j++)
    1304             :   {
    1305        6402 :     long limN = maxss((long)ceil(exp(DL/j)), Nmin);
    1306             :     GEN polsh;
    1307        6402 :     setlg(v, j+1);
    1308        6402 :     polsh = RgV_to_RgX_reverse(v, 0);
    1309      493690 :     for (; n >= limN; n--)
    1310             :     {
    1311      487288 :       pari_sp av2 = avma;
    1312      487288 :       GEN S = divri(mulrr(en, rX_s_eval(polsh, -n)), powuu(n,j));
    1313             :       /* w[n+1] - exp(-n C) * polsh(-n) / (-n)^j */
    1314      487288 :       GEN c = odd(j)? addrr(gel(w,n+1), S) : subrr(gel(w,n+1), S);
    1315      487288 :       affrr(c, gel(w,n)); avma = av2;
    1316      487288 :       en = mulrr(en,eC); /* exp(-n C) */
    1317             :     }
    1318             :   }
    1319         252 :   avma = av0; return w;
    1320             : }
    1321             : 
    1322             : /* erfc via numerical integration : assume real(x)>=1 */
    1323             : static GEN
    1324           7 : cxerfc_r1(GEN x, long prec)
    1325             : {
    1326             :   GEN h, h2, eh2, denom, res, lambda;
    1327             :   long u, v;
    1328           7 :   const double D = prec2nbits_mul(prec, LOG2);
    1329           7 :   const long npoints = (long)ceil(D/M_PI)+1;
    1330           7 :   pari_sp av = avma;
    1331             :   {
    1332           7 :     double t = exp(-2*M_PI*M_PI/D); /* ~exp(-2*h^2) */
    1333           7 :     v = 30; /* bits that fit in both long and double mantissa */
    1334           7 :     u = (long)floor(t*(1L<<v));
    1335             :     /* define exp(-2*h^2) to be u*2^(-v) */
    1336             :   }
    1337           7 :   incrprec(prec);
    1338           7 :   x = gtofp(x,prec);
    1339           7 :   eh2 = sqrtr_abs(rtor(shiftr(dbltor(u),-v),prec));
    1340           7 :   h2 = negr(logr_abs(eh2));
    1341           7 :   h = sqrtr_abs(h2);
    1342           7 :   lambda = gdiv(x,h);
    1343           7 :   denom = gsqr(lambda);
    1344             :   { /* res = h/x + 2*x*h*sum(k=1,npoints,exp(-(k*h)^2)/(lambda^2+k^2)); */
    1345             :     GEN Uk; /* = exp(-(kh)^2) */
    1346           7 :     GEN Vk = eh2;/* = exp(-(2k+1)h^2) */
    1347           7 :     pari_sp av2 = avma;
    1348             :     long k;
    1349             :     /* k = 0 moved out for efficiency */
    1350           7 :     denom = gaddsg(1,denom);
    1351           7 :     Uk = Vk;
    1352           7 :     Vk = mulur(u,Vk); shiftr_inplace(Vk, -v);
    1353           7 :     res = gdiv(Uk, denom);
    1354         210 :     for (k = 1; k < npoints; k++)
    1355             :     {
    1356         203 :       if ((k & 255) == 0) gerepileall(av2,4,&denom,&Uk,&Vk,&res);
    1357         203 :       denom = gaddsg(2*k+1,denom);
    1358         203 :       Uk = mpmul(Uk,Vk);
    1359         203 :       Vk = mulur(u,Vk); shiftr_inplace(Vk, -v);
    1360         203 :       res = gadd(res, gdiv(Uk, denom));
    1361             :     }
    1362             :   }
    1363           7 :   res = gmul(res, gshift(lambda,1));
    1364             :   /* 0 term : */
    1365           7 :   res = gadd(res, ginv(lambda));
    1366           7 :   res = gmul(res, gdiv(gexp(gneg(gsqr(x)), prec), mppi(prec)));
    1367           7 :   if (rtodbl(real_i(x)) < sqrt(D))
    1368             :   {
    1369           7 :     GEN t = gmul(divrr(Pi2n(1,prec),h), x);
    1370           7 :     res = gsub(res, gdivsg(2, cxexpm1(t, prec)));
    1371             :   }
    1372           7 :   return gerepileupto(av,res);
    1373             : }
    1374             : 
    1375             : GEN
    1376          28 : gerfc(GEN x, long prec)
    1377             : {
    1378             :   GEN z, xr, xi, res;
    1379             :   pari_sp av;
    1380             : 
    1381          28 :   x = trans_fix_arg(&prec,&x,&xr,&xi, &av,&res);
    1382          28 :   if (signe(xr) >= 0) {
    1383          21 :     if (cmprs(xr, 1) > 0) /* use numerical integration */
    1384           7 :       z = cxerfc_r1(x, prec);
    1385             :     else
    1386             :     { /* erfc(x) = incgam(1/2,x^2)/sqrt(Pi) */
    1387          14 :       GEN sqrtpi = sqrtr(mppi(prec));
    1388          14 :       z = incgam0(ghalf, gsqr(x), sqrtpi, prec);
    1389          14 :       z = gdiv(z, sqrtpi);
    1390             :     }
    1391             :   }
    1392             :   else
    1393             :   { /* erfc(-x)=2-erfc(x) */
    1394             :     /* FIXME could decrease prec
    1395             :     long size = nbits2extraprec(ceil((pow(rtodbl(gimag(x)),2)-pow(rtodbl(greal(x)),2))/LOG2)));
    1396             :     prec = size > 0 ? prec : prec + size;
    1397             :     */
    1398             :     /* NOT gsubsg(2, ...) : would create a result of
    1399             :      * huge accuracy if re(x)>>1, rounded to 2 by subsequent affc_fixlg... */
    1400           7 :     z = gsub(real2n(1,prec+EXTRAPREC), gerfc(gneg(x), prec));
    1401             :   }
    1402          28 :   avma = av; return affc_fixlg(z, res);
    1403             : }
    1404             : 
    1405             : /***********************************************************************/
    1406             : /**                                                                   **/
    1407             : /**                     FONCTION ZETA DE RIEMANN                      **/
    1408             : /**                                                                   **/
    1409             : /***********************************************************************/
    1410             : static const double log2PI = 1.83787706641;
    1411             : 
    1412             : static double
    1413        3822 : get_xinf(double beta)
    1414             : {
    1415        3822 :   const double maxbeta = 0.06415003; /* 3^(-2.5) */
    1416             :   double x0, y0, x1;
    1417             : 
    1418        3822 :   if (beta < maxbeta) return beta + pow(3*beta, 1.0/3.0);
    1419        3822 :   x0 = beta + M_PI/2.0;
    1420             :   for(;;)
    1421             :   {
    1422        6566 :     y0 = x0*x0;
    1423        6566 :     x1 = (beta+atan(x0)) * (1+y0) / y0 - 1/x0;
    1424        6566 :     if (0.99*x0 < x1) return x1;
    1425        2744 :     x0 = x1;
    1426        2744 :   }
    1427             : }
    1428             : /* optimize for zeta( s + it, prec ), assume |s-1| > 0.1
    1429             :  * (if gexpo(u = s-1) < -5, we use the functional equation s->1-s) */
    1430             : static void
    1431        6692 : optim_zeta(GEN S, long prec, long *pp, long *pn)
    1432             : {
    1433             :   double s, t, alpha, beta, n, B;
    1434             :   long p;
    1435        6692 :   if (typ(S) == t_REAL) {
    1436         273 :     s = rtodbl(S);
    1437         273 :     t = 0.;
    1438             :   } else {
    1439        6419 :     s = rtodbl(gel(S,1));
    1440        6419 :     t = fabs( rtodbl(gel(S,2)) );
    1441             :   }
    1442             : 
    1443        6692 :   B = prec2nbits_mul(prec, LOG2);
    1444        6692 :   if (s > 0 && !t) /* positive real input */
    1445             :   {
    1446         273 :     beta = B + 0.61 + s*(log2PI - log(s));
    1447         546 :     if (beta > 0)
    1448             :     {
    1449         266 :       p = (long)ceil(beta / 2.0);
    1450         266 :       n = fabs(s + 2*p-1)/(2*M_PI);
    1451             :     }
    1452             :     else
    1453             :     {
    1454           7 :       p = 0;
    1455           7 :       n = exp((B - LOG2) / s);
    1456             :     }
    1457             :   }
    1458        6419 :   else if (s <= 0 || t < 0.01) /* s < 0 may occur if s ~ 0 */
    1459        2590 :   { /* TODO: the crude bounds below are generally valid. Optimize ? */
    1460        2590 :     double l,l2, la = 1.; /* heuristic */
    1461        2590 :     double rlog, ilog; dcxlog(s-1,t, &rlog,&ilog);
    1462        2590 :     l2 = (s - 0.5)*rlog - t*ilog; /* = Re( (S - 1/2) log (S-1) ) */
    1463        2590 :     l = (B - l2 + s*log2PI) / (2. * (1.+ log((double)la)));
    1464        2590 :     l2 = dabs(s, t)/2;
    1465        2590 :     if (l < l2) l = l2;
    1466        2590 :     p = (long) ceil(l); if (p < 2) p = 2;
    1467        2590 :     n = 1 + dabs(p+s/2.-.25, t/2) * la / M_PI;
    1468             :   }
    1469             :   else
    1470             :   {
    1471        3829 :     double sn = dabs(s, t), L = log(sn/s);
    1472        3829 :     alpha = B - 0.39 + L + s*(log2PI - log(sn));
    1473        3829 :     beta = (alpha+s)/t - atan(s/t);
    1474        3829 :     p = 0;
    1475        3829 :     if (beta > 0)
    1476             :     {
    1477        3822 :       beta = 1.0 - s + t * get_xinf(beta);
    1478        3822 :       if (beta > 0) p = (long)ceil(beta / 2.0);
    1479             :     }
    1480             :     else
    1481           7 :       if (s < 1.0) p = 1;
    1482        3829 :     n = p? dabs(s + 2*p-1, t) / (2*M_PI) : exp((B-LOG2+L) / s);
    1483             :   }
    1484        6692 :   *pp = p;
    1485        6692 :   *pn = (long)ceil(n);
    1486        6692 :   if (*pp < 0 || *pn < 0) pari_err_OVERFLOW("zeta");
    1487        6692 : }
    1488             : 
    1489             : /* 1/zeta(n) using Euler product. Assume n > 0.
    1490             :  * if (lba != 0) it is log(prec2nbits) we _really_ require */
    1491             : GEN
    1492       10211 : inv_szeta_euler(long n, double lba, long prec)
    1493             : {
    1494             :   GEN z, res;
    1495             :   pari_sp av, av2;
    1496             :   double A, D;
    1497             :   ulong p, lim;
    1498             :   forprime_t S;
    1499             : 
    1500       10211 :   if (n > prec2nbits(prec)) return real_1(prec);
    1501             : 
    1502       10204 :   if (!lba) lba = prec2nbits_mul(prec, LOG2);
    1503       10204 :   D = exp((lba - log((double)(n-1))) / (n-1));
    1504       10204 :   lim = 1 + (ulong)ceil(D);
    1505       10204 :   if (lim < 3) return subir(gen_1,real2n(-n,prec));
    1506        9932 :   res = cgetr(prec); incrprec(prec);
    1507        9932 :   av = avma;
    1508        9932 :   z = subir(gen_1, real2n(-n, prec));
    1509             : 
    1510        9932 :   (void)u_forprime_init(&S, 3, lim);
    1511        9932 :   av2 = avma; A = n / LOG2;
    1512       54044 :   while ((p = u_forprime_next(&S)))
    1513             :   {
    1514       34180 :     long l = prec2nbits(prec) - (long)floor(A * log((double)p)) - BITS_IN_LONG;
    1515             :     GEN h;
    1516             : 
    1517       34180 :     if (l < BITS_IN_LONG) l = BITS_IN_LONG;
    1518       34180 :     l = minss(prec, nbits2prec(l));
    1519       34180 :     h = divrr(z, rpowuu(p, (ulong)n, l));
    1520       34180 :     z = subrr(z, h);
    1521       34180 :     if (gc_needed(av,1))
    1522             :     {
    1523           0 :       if (DEBUGMEM>1) pari_warn(warnmem,"inv_szeta_euler, p = %lu/%lu", p,lim);
    1524           0 :       z = gerepileuptoleaf(av2, z);
    1525             :     }
    1526             :   }
    1527        9932 :   affrr(z, res); avma = av; return res;
    1528             : }
    1529             : 
    1530             : /* assume n even > 0, if iz != NULL, assume iz = 1/zeta(n) */
    1531             : GEN
    1532       10033 : bernreal_using_zeta(long n, GEN iz, long prec)
    1533             : {
    1534       10033 :   long l = prec+EXTRAPRECWORD;
    1535             :   GEN z;
    1536             : 
    1537       10033 :   if (!iz) iz = inv_szeta_euler(n, 0., l);
    1538       10033 :   z = divrr(mpfactr(n, l), mulrr(powru(Pi2n(1, l), n), iz));
    1539       10033 :   shiftr_inplace(z, 1); /* 2 * n! * zeta(n) / (2Pi)^n */
    1540       10033 :   if ((n & 3) == 0) setsigne(z, -1);
    1541       10033 :   return z;
    1542             : }
    1543             : 
    1544             : /* assume n even > 0. Faster than standard bernfrac for n >= 6 */
    1545             : GEN
    1546       10026 : bernfrac_using_zeta(long n)
    1547             : {
    1548       10026 :   pari_sp av = avma;
    1549       10026 :   GEN iz, a, d, D = divisorsu(n >> 1);
    1550       10026 :   long i, prec, l = lg(D);
    1551             :   double t, u;
    1552             : 
    1553       10026 :   d = utoipos(6); /* 2 * 3 */
    1554       44347 :   for (i = 2; i < l; i++) /* skip 1 */
    1555             :   { /* Clausen - von Staudt */
    1556       34321 :     ulong p = 2*D[i] + 1;
    1557       34321 :     if (uisprime(p)) d = muliu(d, p);
    1558             :   }
    1559             :   /* 1.712086 = ??? */
    1560       10026 :   t = log( gtodouble(d) ) + (n + 0.5) * log((double)n) - n*(1+log2PI) + 1.712086;
    1561       10026 :   u = t / LOG2; prec = nbits2prec((long)ceil(u) + BITS_IN_LONG);
    1562       10026 :   iz = inv_szeta_euler(n, t, prec);
    1563       10026 :   a = roundr( mulir(d, bernreal_using_zeta(n, iz, prec)) );
    1564       10026 :   return gerepilecopy(av, mkfrac(a, d));
    1565             : }
    1566             : 
    1567             : static int
    1568        2639 : bernreal_use_zeta_i(long n, long prec)
    1569             : {
    1570        2639 :   return (n+0.5) * log((double)n) -n*(1+log2PI) > prec2nbits_mul(prec, LOG2);
    1571             : }
    1572             : static int
    1573        1631 : bernreal_use_zeta(long n, long prec)
    1574             : {
    1575        1631 :   if (bernzone)
    1576             :   {
    1577        1596 :     long k = n >> 1;
    1578        1596 :     if (n+1 < lg(bernzone))
    1579             :     {
    1580         504 :       GEN B = gel(bernzone,k+1);
    1581         504 :       if (typ(B) != t_REAL || realprec(B) >= prec) return 0;
    1582             :     }
    1583             :   }
    1584        1127 :   return bernreal_use_zeta_i(n, prec);
    1585             : }
    1586             : 
    1587             : /* Return B_n */
    1588             : GEN
    1589     8592992 : bernreal(long n, long prec)
    1590             : {
    1591             :   GEN B, storeB;
    1592             :   long k, lbern;
    1593     8592992 :   if (n < 0) pari_err_DOMAIN("bernreal", "index", "<", gen_0, stoi(n));
    1594     8592985 :   if (n == 0) return real_1(prec);
    1595     8592964 :   if (n == 1) return real_m2n(-1,prec); /*-1/2*/
    1596     8592964 :   if (odd(n)) return real_0(prec);
    1597             : 
    1598     8592964 :   k = n >> 1;
    1599     8592964 :   if (!bernzone && k < 100) mpbern(k, prec);
    1600     8592964 :   lbern = bernzone? lg(bernzone): 0;
    1601     8592964 :   if (k < lbern)
    1602             :   {
    1603     8591452 :     B = gel(bernzone,k);
    1604     8591452 :     if (typ(B) != t_REAL) return fractor(B, prec);
    1605     2915307 :     if (realprec(B) >= prec) return rtor(B, prec);
    1606             :   }
    1607             :   /* not cached, must compute */
    1608        1512 :   if (bernreal_use_zeta_i(n, prec))
    1609           7 :     B = storeB = bernreal_using_zeta(n, NULL, prec);
    1610             :   else
    1611             :   {
    1612        1505 :     storeB = bernfrac_using_zeta(n);
    1613        1505 :     B = fractor(storeB, prec);
    1614             :   }
    1615        1512 :   if (k < lbern)
    1616             :   {
    1617           0 :     GEN old = gel(bernzone, k);
    1618           0 :     gel(bernzone, k) = gclone(storeB);
    1619           0 :     gunclone(old);
    1620             :   }
    1621        1512 :   return B;
    1622             : }
    1623             : 
    1624             : /* zeta(a*j+b), j=0..N-1, b>1, using sumalt. Johansonn'b thesis, Algo 4.7.1 */
    1625             : static GEN
    1626        7588 : veczetas(long a, long b, long N, long prec)
    1627             : {
    1628        7588 :   pari_sp av = avma;
    1629        7588 :   const long n = ceil(2 + prec2nbits_mul(prec, LOG2/1.7627));
    1630             :   long j, k;
    1631        7588 :   GEN c, d, z = zerovec(N);
    1632        7588 :   c = d = int2n(2*n-1);
    1633      765374 :   for (k = n; k; k--)
    1634             :   {
    1635      757786 :     GEN u, t = divii(d, powuu(k,b));
    1636      757786 :     if (!odd(k)) t = negi(t);
    1637      757786 :     gel(z,1) = addii(gel(z,1), t);
    1638      757786 :     u = powuu(k,a);
    1639     5176121 :     for (j = 1; j < N; j++)
    1640             :     {
    1641     4426105 :       t = divii(t,u); if (!signe(t)) break;
    1642     4418335 :       gel(z,j+1) = addii(gel(z,j+1), t);
    1643             :     }
    1644      757786 :     c = muluui(k,2*k-1,c);
    1645      757786 :     c = diviuuexact(c, 2*(n-k+1),n+k-1);
    1646      757786 :     d = addii(d,c);
    1647      757786 :     if (gc_needed(av,3))
    1648             :     {
    1649          11 :       if(DEBUGMEM>1) pari_warn(warnmem,"zetaBorwein, k = %ld", k);
    1650          11 :       gerepileall(av, 3, &c,&d,&z);
    1651             :     }
    1652             :   }
    1653       45197 :   for (j = 0; j < N; j++)
    1654             :   {
    1655       37609 :     long u = b+a*j-1;
    1656       37609 :     gel(z,j+1) = rdivii(shifti(gel(z,j+1), u), subii(shifti(d,u), d), prec);
    1657             :   }
    1658        7588 :   return gerepilecopy(av, z);
    1659             : }
    1660             : /* zeta(a*j+b), j=0..N-1, b>1, using sumalt */
    1661             : GEN
    1662        7630 : veczeta(GEN a, GEN b, long N, long prec)
    1663             : {
    1664             :   pari_sp av;
    1665             :   long n, j, k;
    1666             :   GEN L, c, d, z;
    1667        7630 :   if (typ(a) == t_INT && typ(b) == t_INT)
    1668        7588 :     return veczetas(itos(a),  itos(b), N, prec);
    1669          42 :   av = avma; z = zerovec(N);
    1670          42 :   n = ceil(2 + prec2nbits_mul(prec, LOG2/1.7627));
    1671          42 :   c = d = int2n(2*n-1);
    1672        8320 :   for (k = n; k; k--)
    1673             :   {
    1674             :     GEN u, t;
    1675        8278 :     L = logr_abs(utor(k, prec)); /* log(k) */
    1676        8278 :     t = gdiv(d, gexp(gmul(b, L), prec)); /* d / k^b */
    1677        8278 :     if (!odd(k)) t = gneg(t);
    1678        8278 :     gel(z,1) = gadd(gel(z,1), t);
    1679        8278 :     u = gexp(gmul(a, L), prec);
    1680      576704 :     for (j = 1; j < N; j++)
    1681             :     {
    1682      569792 :       t = gdiv(t,u); if (gexpo(t) < 0) break;
    1683      568426 :       gel(z,j+1) = gadd(gel(z,j+1), t);
    1684             :     }
    1685        8278 :     c = muluui(k,2*k-1,c);
    1686        8278 :     c = diviuuexact(c, 2*(n-k+1),n+k-1);
    1687        8278 :     d = addii(d,c);
    1688        8278 :     if (gc_needed(av,3))
    1689             :     {
    1690           0 :       if(DEBUGMEM>1) pari_warn(warnmem,"veczeta, k = %ld", k);
    1691           0 :       gerepileall(av, 3, &c,&d,&z);
    1692             :     }
    1693             :   }
    1694          42 :   L = mplog2(prec);
    1695        2940 :   for (j = 0; j < N; j++)
    1696             :   {
    1697        2898 :     GEN u = gsubgs(gadd(b, gmulgs(a,j)), 1);
    1698        2898 :     GEN w = gexp(gmul(u, L), prec); /* 2^u */
    1699        2898 :     gel(z,j+1) = gdiv(gmul(gel(z,j+1), w), gmul(d,gsubgs(w,1)));
    1700             :   }
    1701          42 :   return gerepilecopy(av, z);
    1702             : }
    1703             : 
    1704             : /* zeta(s) using sumalt, case h=0,N=1. Assume s > 1 */
    1705             : static GEN
    1706         347 : zetaBorwein(long s, long prec)
    1707             : {
    1708         347 :   pari_sp av = avma;
    1709         347 :   const long n = ceil(2 + prec2nbits_mul(prec, LOG2/1.7627));
    1710             :   long k;
    1711         347 :   GEN c, d, z = gen_0;
    1712         347 :   c = d = int2n(2*n-1);
    1713       59250 :   for (k = n; k; k--)
    1714             :   {
    1715       58903 :     GEN t = divii(d, powuu(k,s));
    1716       58903 :     z = odd(k)? addii(z,t): subii(z,t);
    1717       58903 :     c = muluui(k,2*k-1,c);
    1718       58903 :     c = diviuuexact(c, 2*(n-k+1),n+k-1);
    1719       58903 :     d = addii(d,c);
    1720       58903 :     if (gc_needed(av,3))
    1721             :     {
    1722           0 :       if(DEBUGMEM>1) pari_warn(warnmem,"zetaBorwein, k = %ld", k);
    1723           0 :       gerepileall(av, 3, &c,&d,&z);
    1724             :     }
    1725             :   }
    1726         347 :   z = rdivii(shifti(z, s-1), subii(shifti(d,s-1), d), prec);
    1727         347 :   return gerepileuptoleaf(av, z);
    1728             : }
    1729             : 
    1730             : /* assume k != 1 */
    1731             : GEN
    1732        3066 : szeta(long k, long prec)
    1733             : {
    1734        3066 :   pari_sp av = avma;
    1735             :   GEN y;
    1736             :   double p;
    1737             : 
    1738             :   /* treat trivial cases */
    1739        3066 :   if (!k) { y = real2n(-1, prec); setsigne(y,-1); return y; }
    1740        3038 :   if (k < 0)
    1741             :   {
    1742         882 :     if ((k&1) == 0) return gen_0;
    1743             :     /* the one value such that k < 0 and 1 - k < 0, due to overflow */
    1744         882 :     if ((ulong)k == (HIGHBIT | 1))
    1745           0 :       pari_err_OVERFLOW("zeta [large negative argument]");
    1746         882 :     k = 1-k;
    1747         882 :     y = bernreal(k, prec); togglesign(y);
    1748         882 :     return gerepileuptoleaf(av, divru(y, k));
    1749             :   }
    1750        2156 :   if (k > prec2nbits(prec)+1) return real_1(prec);
    1751        2149 :   if ((k&1) == 0)
    1752             :   {
    1753        1631 :     if (bernreal_use_zeta(k, prec))
    1754           7 :       y = invr( inv_szeta_euler(k, 0, prec) );
    1755             :     else
    1756             :     {
    1757        1624 :       y = mulrr(powru(Pi2n(1, prec), k), bernreal(k, prec));
    1758        1624 :       y = divrr(y, mpfactr(k,prec));
    1759        1624 :       setsigne(y, 1);
    1760        1624 :       shiftr_inplace(y, -1);
    1761             :     }
    1762        1631 :     return gerepileuptoleaf(av, y);
    1763             :   }
    1764             :   /* k > 1 odd */
    1765         518 :   p = prec2nbits_mul(prec,0.393); /* 0.393 ~ 1/log_2(3+sqrt(8)) */
    1766         518 :   if (log2(p * log(p))*k > prec2nbits(prec))
    1767         171 :     return gerepileuptoleaf(av, invr( inv_szeta_euler(k, 0, prec) ));
    1768         347 :   return zetaBorwein(k, prec);
    1769             : }
    1770             : 
    1771             : /* s0 a t_INT, t_REAL or t_COMPLEX.
    1772             :  * If a t_INT, assume it's not a trivial case (i.e we have s0 > 1, odd)
    1773             :  * */
    1774             : static GEN
    1775        6706 : czeta(GEN s0, long prec)
    1776             : {
    1777        6706 :   GEN s, u, y, res, tes, sig, tau, invn2, unr, Ns, funeq_factor = NULL;
    1778             :   long i, nn, lim, lim2;
    1779        6706 :   pari_sp av0 = avma, av;
    1780             :   pari_timer T;
    1781             : 
    1782        6706 :   if (DEBUGLEVEL>2) timer_start(&T);
    1783        6706 :   s = trans_fix_arg(&prec,&s0,&sig,&tau,&av,&res);
    1784        6706 :   if (typ(s0) == t_INT) return gerepileupto(av, gzeta(s0, prec));
    1785        6706 :   if (!signe(tau)) /* real */
    1786             :   {
    1787         280 :     long e = expo(sig);
    1788         280 :     if (e >= -5 && (signe(sig) <= 0 || e < -1))
    1789             :     { /* s < 1/2 */
    1790           7 :       s = subsr(1, s);
    1791           7 :       funeq_factor = gen_1;
    1792             :     }
    1793             :   }
    1794             :   else
    1795             :   {
    1796        6426 :     u = gsubsg(1, s); /* temp */
    1797        6426 :     if (gexpo(u) < -5 || ((signe(sig) <= 0 || expo(sig) < -1) && gexpo(s) > -5))
    1798             :     { /* |1-s| < 1/32  || (real(s) < 1/2 && |imag(s)| > 1/32) */
    1799          21 :       s = u;
    1800          21 :       funeq_factor = gen_1;
    1801             :     }
    1802             :   }
    1803             : 
    1804        6706 :   if (funeq_factor)
    1805             :   { /* s <--> 1-s */
    1806             :     GEN t;
    1807          28 :     sig = real_i(s);
    1808             :     /* Gamma(s) (2Pi)^-s 2 cos(Pi s/2) */
    1809          28 :     t = gmul(ggamma(gprec_w(s,prec),prec), gpow(Pi2n(1,prec), gneg(s), prec));
    1810          28 :     funeq_factor = gmul2n(gmul(t, gcos(gmul(Pi2n(-1,prec),s), prec)), 1);
    1811             :   }
    1812        6706 :   if (gcmpgs(sig, prec2nbits(prec) + 1) > 0) { /* zeta(s) = 1 */
    1813          14 :     if (!funeq_factor) { avma = av0; return real_1(prec); }
    1814           7 :     return gerepileupto(av0, funeq_factor);
    1815             :   }
    1816        6692 :   optim_zeta(s, prec, &lim, &nn);
    1817        6692 :   if (DEBUGLEVEL>2) err_printf("lim, nn: [%ld, %ld]\n", lim, nn);
    1818        6692 :   incrprec(prec); unr = real_1(prec); /* one extra word of precision */
    1819             : 
    1820        6692 :   Ns = vecpowug(nn, gneg(s), prec);
    1821        6692 :   y = gmul2n(gel(Ns,nn), -1);
    1822        6692 :   for (i = nn-1; i; i--) y = gadd(y, gel(Ns,i));
    1823        6692 :   if (DEBUGLEVEL>2) timer_printf(&T,"sum from 1 to N-1");
    1824             : 
    1825        6692 :   invn2 = divri(unr, sqru(nn)); lim2 = lim<<1;
    1826        6692 :   mpbern(lim,prec);
    1827        6692 :   tes = bernreal(lim2, prec);
    1828             :   {
    1829             :     GEN s1, s2, s3, s4, s5;
    1830             :     pari_sp av2;
    1831        6692 :     s1 = gsub(gmul2n(s,1), unr);
    1832        6692 :     s2 = gmul(s, gsub(s,unr));
    1833        6692 :     s3 = gmul2n(invn2,3);
    1834        6692 :     av2 = avma;
    1835        6692 :     s4 = gmul(invn2, gmul2n(gaddsg(4*lim-2,s1),1));
    1836        6692 :     s5 = gmul(invn2, gadd(s2, gmulsg(lim2, gaddgs(s1, lim2))));
    1837      666135 :     for (i = lim2-2; i>=2; i -= 2)
    1838             :     {
    1839      659443 :       s5 = gsub(s5, s4);
    1840      659443 :       s4 = gsub(s4, s3);
    1841      659443 :       tes = gadd(bernreal(i,prec), divgunu(gmul(s5,tes), i+1));
    1842      659443 :       if (gc_needed(av2,3))
    1843             :       {
    1844           0 :         if(DEBUGMEM>1) pari_warn(warnmem,"czeta");
    1845           0 :         gerepileall(av2,3, &tes,&s5,&s4);
    1846             :       }
    1847             :     }
    1848        6692 :     u = gmul(gmul(tes,invn2), gmul2n(s2, -1));
    1849        6692 :     tes = gmulsg(nn, gaddsg(1, u));
    1850             :   }
    1851        6692 :   if (DEBUGLEVEL>2) timer_printf(&T,"Bernoulli sum");
    1852             :   /* y += tes n^(-s) / (s-1) */
    1853        6692 :   y = gadd(y, gmul(tes, gdiv(gel(Ns,nn), gsub(s, unr))));
    1854        6692 :   if (funeq_factor) y = gmul(y, funeq_factor);
    1855        6692 :   avma = av; return affc_fixlg(y,res);
    1856             : }
    1857             : 
    1858             : #if 0
    1859             : /* return P mod x^n where P is polynomial in x */
    1860             : static GEN
    1861             : pol_mod_xn(GEN P, long n)
    1862             : {
    1863             :   long j, l = lg(P), N = n+2;
    1864             :   GEN R;
    1865             :   if (l > N) l = N;
    1866             :   R = cgetg(N, t_POL); R[1] = evalvarn(0);
    1867             :   for (j = 2; j < l; j++) gel(R,j) = gel(P,j);
    1868             :   return normalizepol_lg(R, n+2);
    1869             : }
    1870             : 
    1871             : /* compute the values of the twisted partial
    1872             :    zeta function Z_f(a, c, s) for a in va */
    1873             : GEN
    1874             : twistpartialzeta(GEN q, long f, long c, GEN va, GEN cff)
    1875             : {
    1876             :   long j, k, lva = lg(va)-1, N = lg(cff)-1;
    1877             :   pari_sp av, av2;
    1878             :   GEN Ax, Cx, Bx, Dx, x = pol_x(0), y = pol_x(1);
    1879             :   GEN cyc, psm, rep, eta, etaf;
    1880             : 
    1881             :   cyc = gdiv(gsubgs(gpowgs(y, c), 1), gsubgs(y, 1));
    1882             :   psm = polsym(cyc, degpol(cyc) - 1);
    1883             :   eta = mkpolmod(y, cyc);
    1884             :   etaf = gpowgs(eta,f);
    1885             :   av = avma;
    1886             :   Ax  = gsubgs(gpowgs(gaddgs(x, 1), f), 1);
    1887             :   Ax  = gdiv(gmul(Ax, etaf), gsubsg(1, etaf));
    1888             :   Ax  = gerepileupto(av, RgX_to_FqX(Ax, cyc, q));
    1889             :   Cx  = Ax;
    1890             :   Bx  = gen_1;
    1891             :   av  = avma;
    1892             :   for (j = 2; j <= N; j++)
    1893             :   {
    1894             :     Bx = gadd(Bx, Cx);
    1895             :     Bx = FpXQX_red(Bx, cyc, q);
    1896             :     Cx = FpXQX_mul(Cx, Ax, cyc, q);
    1897             :     Cx = pol_mod_xn(Cx, N);
    1898             :     if (gequal0(Cx)) break;
    1899             :     if (gc_needed(av, 1))
    1900             :     {
    1901             :       if(DEBUGMEM>1) pari_warn(warnmem, "twistpartialzeta (1), j = %ld/%ld", j, N);
    1902             :       gerepileall(av, 2, &Cx, &Bx);
    1903             :     }
    1904             :   }
    1905             :   Bx  = lift_shallow(gmul(ginv(gsubsg(1, etaf)), Bx));
    1906             :   Bx  = gerepileupto(av, RgX_to_FqX(Bx, cyc, q));
    1907             :   Cx = lift_shallow(gmul(eta, gaddsg(1, x)));
    1908             :   Dx = pol_1(varn(x));
    1909             :   av2 = avma;
    1910             :   for (j = lva; j > 1; j--)
    1911             :   {
    1912             :     GEN Ex;
    1913             :     long e = va[j] - va[j-1];
    1914             :     if (e == 1)
    1915             :       Ex = Cx;
    1916             :     else
    1917             :       /* e is very small in general and actually very rarely different
    1918             :          to 1, it is always 1 for zetap (so it should be OK not to store
    1919             :          them or to compute them in a smart way) */
    1920             :       Ex = gpowgs(Cx, e);
    1921             :     Dx = gaddsg(1, FpXQX_mul(Dx, Ex, cyc, q));
    1922             :     if (gc_needed(av2, 1))
    1923             :     {
    1924             :       if(DEBUGMEM>1)
    1925             :         pari_warn(warnmem, "twistpartialzeta (2), j = %ld/%ld", lva-j, lva);
    1926             :       Dx = gerepileupto(av2, FpXQX_red(Dx, cyc, q));
    1927             :     }
    1928             :   }
    1929             :   Dx = FpXQX_mul(Dx, Cx, cyc, q); /* va[1] = 1 */
    1930             :   Bx = gerepileupto(av, FpXQX_mul(Dx, Bx, cyc, q));
    1931             :   rep = gen_0;
    1932             :   av2 = avma;
    1933             :   for (k = 1; k <= N; k++)
    1934             :   {
    1935             :     GEN p2, ak = polcoeff_i(Bx, k, 0);
    1936             :     p2  = quicktrace(ak, psm);
    1937             :     rep = modii(addii(rep, mulii(gel(cff, k), p2)), q);
    1938             :     if (gc_needed(av2, 1))
    1939             :     {
    1940             :       if(DEBUGMEM>1) pari_warn(warnmem, "twistpartialzeta (3), j = %ld/%ld", k, N);
    1941             :       rep = gerepileupto(av2, rep);
    1942             :     }
    1943             :   }
    1944             :   return rep;
    1945             : }
    1946             : 
    1947             : /* initialize the roots of unity for the computation
    1948             :    of the Teichmuller character (also the values of f and c) */
    1949             : GEN
    1950             : init_teich(ulong p, GEN q, long prec)
    1951             : {
    1952             :   GEN vz, gp = utoipos(p);
    1953             :   pari_sp av = avma;
    1954             :   long j;
    1955             : 
    1956             :   if (p == 2UL)
    1957             :     return NULL;
    1958             :   else
    1959             :   { /* primitive (p-1)-th root of 1 */
    1960             :     GEN z, z0 = Zp_sqrtnlift(gen_1, utoipos(p-1), pgener_Fp(gp), gp, prec);
    1961             :     z = z0;
    1962             :     vz = cgetg(p, t_VEC);
    1963             :     for (j = 1; j < (long)p-2; j++)
    1964             :     {
    1965             :       gel(vz, umodiu(z, p)) = z; /* z = z0^i */
    1966             :       z = modii(mulii(z, z0), q);
    1967             :     }
    1968             :     gel(vz, umodiu(z, p)) = z; /* z = z0^(p-2) */
    1969             :     gel(vz,1) = gen_1; /* z0^(p-1) */
    1970             :   }
    1971             :   return gerepileupto(av, gcopy(vz));
    1972             : }
    1973             : 
    1974             : /* compute phi^(m)_s(x); s must be an integer */
    1975             : GEN
    1976             : phi_ms(ulong p, GEN q, long m, GEN s, long x, GEN vz)
    1977             : {
    1978             :   long xp = x % p;
    1979             :   GEN p1, p2;
    1980             : 
    1981             :   if (!xp) return gen_0;
    1982             :   if (vz)
    1983             :     p1 =gel(vz,xp); /* vz[x] = Teichmuller(x) */
    1984             :   else
    1985             :     p1 = (x & 2)? gen_m1: gen_1;
    1986             :   p1 = Fp_pow(p1, addis(s, m), q);
    1987             :   p2 = Fp_pow(stoi(x), negi(s), q);
    1988             :   return modii(mulii(p1,p2), q);
    1989             : }
    1990             : 
    1991             : /* compute the first N coefficients of the Mahler expansion
    1992             :    of phi^m_s skipping the first one (which is zero) */
    1993             : GEN
    1994             : coeff_of_phi_ms(ulong p, GEN q, long m, GEN s, long N, GEN vz)
    1995             : {
    1996             :   GEN qs2 = shifti(q, -1), cff = zerovec(N);
    1997             :   pari_sp av;
    1998             :   long k, j;
    1999             : 
    2000             :   av = avma;
    2001             :   for (k = 1; k <= N; k++)
    2002             :   {
    2003             :     gel(cff, k) = phi_ms(p, q, m, s, k, vz);
    2004             :     if (gc_needed(av, 2))
    2005             :     {
    2006             :       if(DEBUGMEM>1)
    2007             :         pari_warn(warnmem, "coeff_of_phi_ms (1), k = %ld/%ld", N-k, N);
    2008             :       cff = gerepileupto(av, gcopy(cff));
    2009             :     }
    2010             :   }
    2011             :   for (j = N; j > 1; j--)
    2012             :   {
    2013             :     GEN b = subii(gel(cff, j), gel(cff, j-1));
    2014             :     gel(cff, j) = centermodii(b, q, qs2);
    2015             :     if (gc_needed(av, 2))
    2016             :     {
    2017             :       if(DEBUGMEM>1)
    2018             :         pari_warn(warnmem, "coeff_of_phi_ms (2), j = %ld/%ld", N-j, N);
    2019             :       cff = gerepileupto(av, gcopy(cff));
    2020             :     }
    2021             :   }
    2022             :   for (k = 1; k < N; k++)
    2023             :     for (j = N; j > k; j--)
    2024             :     {
    2025             :       GEN b = subii(gel(cff, j), gel(cff, j-1));
    2026             :       gel(cff, j) = centermodii(b, q, qs2);
    2027             :       if (gc_needed(av, 2))
    2028             :       {
    2029             :         if(DEBUGMEM>1)
    2030             :           pari_warn(warnmem, "coeff_of_phi_ms (3), (k,j) = (%ld,%ld)/%ld",
    2031             :               k, N-j, N);
    2032             :         cff = gerepileupto(av, gcopy(cff));
    2033             :       }
    2034             :     }
    2035             :   k = N; while(gequal0(gel(cff, k))) k--;
    2036             :   setlg(cff, k+1);
    2037             :   if (DEBUGLEVEL > 2)
    2038             :     err_printf("  coeff_of_phi_ms: %ld coefficients kept out of %ld\n",
    2039             :                k, N);
    2040             :   return gerepileupto(av, cff);
    2041             : }
    2042             : 
    2043             : static long
    2044             : valfact(long N, ulong p)
    2045             : {
    2046             :   long f = 0;
    2047             :   while (N > 1)
    2048             :   {
    2049             :     N /= p;
    2050             :     f += N;
    2051             :   }
    2052             :   return f;
    2053             : }
    2054             : 
    2055             : static long
    2056             : number_of_terms(ulong p, long prec)
    2057             : {
    2058             :   long N, f;
    2059             : 
    2060             :   if (prec == 0) return p;
    2061             :   N = (long)((p-1)*prec + (p>>1)*(log2(prec)/log2(p)));
    2062             :   N = p*(N/p);
    2063             :   f = valfact(N, p);
    2064             :   while (f > prec)
    2065             :   {
    2066             :     N = p*(N/p) - 1;
    2067             :     f -= u_lval(N+1, p);
    2068             :   }
    2069             :   while (f < prec)
    2070             :   {
    2071             :     N = p*(N/p+1);
    2072             :     f += u_lval(N, p);
    2073             :   }
    2074             :   return N;
    2075             : }
    2076             : 
    2077             : static GEN
    2078             : zetap(GEN s)
    2079             : {
    2080             :   ulong p;
    2081             :   long N, f, c, prec = precp(s);
    2082             :   pari_sp av = avma;
    2083             :   GEN gp, q, vz, is, cff, val, va, cft;
    2084             : 
    2085             :   if (valp(s) < 0) pari_err_DOMAIN("zetap", "v_p(s)", "<", gen_0, s);
    2086             :   if (!prec) prec = 1;
    2087             : 
    2088             :   gp = gel(s,2); p = itou(gp);
    2089             :   is = gtrunc(s);  /* make s an integer */
    2090             : 
    2091             :   N  = number_of_terms(p, prec);
    2092             :   q  = powiu(gp, prec);
    2093             : 
    2094             :   /* initialize the roots of unity for the computation
    2095             :      of the Teichmuller character (also the values of f and c) */
    2096             :   if (DEBUGLEVEL > 1) err_printf("zetap: computing (p-1)th roots of 1\n");
    2097             :   vz = init_teich(p, q, prec);
    2098             :   if (p == 2UL) {  f = 4; c = 3; } else { f = (long)p; c = 2; }
    2099             : 
    2100             :   /* compute the first N coefficients of the Mahler expansion
    2101             :      of phi^(-1)_s skipping the first one (which is zero) */
    2102             :   if (DEBUGLEVEL > 1)
    2103             :     err_printf("zetap: computing Mahler expansion of phi^(-1)_s\n");
    2104             :   cff = coeff_of_phi_ms(p, q, -1, is, N, vz);
    2105             : 
    2106             :   /* compute the coefficients of the power series corresponding
    2107             :      to the twisted partial zeta function Z_f(a, c, s) for a in va */
    2108             :   /* The line below looks a bit stupid but it is to keep the
    2109             :      possibility of later adding p-adic Dirichlet L-functions */
    2110             :   va = identity_perm(f - 1);
    2111             :   if (DEBUGLEVEL > 1)
    2112             :     err_printf("zetap: computing values of twisted partial zeta functions\n");
    2113             :   val = twistpartialzeta(q, f, c, va, cff);
    2114             : 
    2115             :   /* sum over all a's the coefficients of the twisted
    2116             :      partial zeta functions and integrate */
    2117             :   if (DEBUGLEVEL > 1)
    2118             :     err_printf("zetap: multiplying by correcting factor\n");
    2119             : 
    2120             :   /* multiply by the corrective factor */
    2121             :   cft = gsubgs(gmulsg(c, phi_ms(p, q, -1, is, c, vz)), 1);
    2122             :   val = gdiv(val, cft);
    2123             : 
    2124             :   /* adjust the precision and return */
    2125             :   return gerepileupto(av, cvtop(val, gp, prec));
    2126             : }
    2127             : #else
    2128             : static GEN
    2129          21 : hurwitz_p(GEN cache, GEN s, GEN x, GEN p, long prec)
    2130             : {
    2131          21 :   GEN S, x2, x2j, s_1 = gsubgs(s,1);
    2132          21 :   long j, J = lg(cache)-2;
    2133          21 :   x = ginv(gadd(x, zeropadic(p, prec)));
    2134          21 :   x2 = gsqr(x);
    2135          21 :   S = gmul2n(gmul(s_1, x), -1);
    2136          21 :   x2j = gen_1;
    2137         126 :   for (j = 0;; j++)
    2138             :   {
    2139         126 :     S = gadd(S, gmul(gel(cache, j+1), x2j));
    2140         126 :     if (j == J) break;
    2141         105 :     x2j = gmul(x2, x2j);
    2142         105 :   }
    2143          21 :   return gmul(gdiv(S, s_1), Qp_exp(gmul(s_1, Qp_log(x))));
    2144             : }
    2145             : 
    2146             : static GEN
    2147          14 : init_cache(long J, GEN s)
    2148             : {
    2149          14 :   GEN C = gen_1, cache = bernvec(J);
    2150             :   long j;
    2151             : 
    2152          77 :   for (j = 1; j <= J; j++)
    2153             :   { /* B_{2j} * binomial(1-s, 2j) */
    2154          63 :     GEN t = gmul(gaddgs(s, 2*j-3), gaddgs(s, 2*j-2));
    2155          63 :     C = gdiv(gmul(C, t), mulss(2*j, 2*j-1));
    2156          63 :     gel(cache, j+1) = gmul(gel(cache, j+1), C);
    2157             :   }
    2158          14 :   return cache;
    2159             : }
    2160             : 
    2161             : static GEN
    2162          14 : zetap(GEN s)
    2163             : {
    2164          14 :   pari_sp av = avma;
    2165          14 :   GEN cache, S, gp = gel(s,2);
    2166          14 :   ulong a, p = itou(gp);
    2167          14 :   long J, prec = valp(s) + precp(s);
    2168             : 
    2169          14 :   if (prec <= 0) prec = 1;
    2170          14 :   if (p == 2) {
    2171           7 :     J = ((long)(1+ceil((prec+1.)/2))) >> 1;
    2172           7 :     cache = init_cache(J, s);
    2173           7 :     S = gmul2n(hurwitz_p(cache, s, gmul2n(gen_1, -2), gen_2, prec), -1);
    2174             :   } else {
    2175           7 :     J = (prec+2) >> 1;
    2176           7 :     cache = init_cache(J, s);
    2177           7 :     S = gen_0;
    2178          21 :     for (a = 1; a <= (p-1)>>1; a++)
    2179          14 :       S = gadd(S, hurwitz_p(cache, s, gdivsg(a, gp), gp, prec));
    2180           7 :     S = gdiv(gmul2n(S, 1), gp);
    2181             :   }
    2182          14 :   return gerepileupto(av, S);
    2183             : }
    2184             : #endif
    2185             : 
    2186             : GEN
    2187        8575 : gzeta(GEN x, long prec)
    2188             : {
    2189        8575 :   if (gequal1(x)) pari_err_DOMAIN("zeta", "argument", "=", gen_1, x);
    2190        8575 :   switch(typ(x))
    2191             :   {
    2192             :     case t_INT:
    2193        1806 :       if (is_bigint(x))
    2194             :       {
    2195          21 :         if (signe(x) > 0) return real_1(prec);
    2196          14 :         if (signe(x) < 0 && mod2(x) == 0) return real_0(prec);
    2197           7 :         pari_err_OVERFLOW("zeta [large negative argument]");
    2198             :       }
    2199        1785 :       return szeta(itos(x),prec);
    2200        6706 :     case t_REAL: case t_COMPLEX: return czeta(x,prec);
    2201          14 :     case t_PADIC: return zetap(x);
    2202           7 :     case t_SER: pari_err_IMPL("zeta(t_SER)");
    2203             :   }
    2204          42 :   return trans_eval("zeta",gzeta,x,prec);
    2205             : }
    2206             : 
    2207             : /********************************************************/
    2208             : /*                   Hurwitz zeta function              */
    2209             : /********************************************************/
    2210             : 
    2211             : /* If s=1 return -psi(x), else zeta(s,x). In particular, if x=1 return zeta(s),
    2212             :  * and Euler's constant if s=1 */
    2213             : GEN
    2214          42 : zetahurwitz(GEN s, GEN x, long prec)
    2215             : {
    2216          42 :   pari_sp av = avma;
    2217             :   GEN al, ral, ral0, Nx, S1, S2, S3, N2;
    2218          42 :   long j, k, m, N, bitprec = prec2nbits(prec);
    2219             : 
    2220          42 :   if (!is_real_t(typ(x))) pari_err_TYPE("zetahurwitz", x);
    2221          42 :   if (gsigne(x) <= 0) pari_err_DOMAIN("zetahurwitz", "x", "<=", gen_0, x);
    2222          42 :   al = gneg(s); ral = greal(al);
    2223          42 :   if (!gequal1(s) && gcmpgs(ral, -1) >= 0)
    2224           0 :     pari_err_IMPL("Re(s) <= 1 in zetahurwitz");
    2225          42 :   ral0 = ground(ral);
    2226          42 :   if (signe(ral0) >= 0 && gexpo(gsub(al, ral0)) < 17-bitprec)
    2227             :   { /* al ~ non negative integer */
    2228           0 :     k = itos(gceil(ral)) + 1;
    2229           0 :     if (odd(k)) k++;
    2230           0 :     N = 1;
    2231             :   }
    2232             :   else
    2233             :   {
    2234             :     GEN C;
    2235          42 :     k = maxss(itos(gceil(gadd(ral, ghalf))) + 1, 50);
    2236          42 :     k = maxss(k, (long)(0.24*bitprec));
    2237          42 :     if (odd(k)) k++;
    2238          42 :     C = gmulsg(2, gmul(binomial(al, k+1), gdivgs(bernfrac(k+2), k+2)));
    2239          42 :     C = gmul2n(gabs(C,LOWDEFAULTPREC), bitprec);
    2240          42 :     C = gadd(gpow(C, ginv(gsubsg(k+1, ral)), LOWDEFAULTPREC), gsubsg(1,x));
    2241          42 :     N = itos(gceil(C));
    2242          42 :     if (N < 1) N = 1;
    2243             :   }
    2244          42 :   S1 = real_0(prec);
    2245          42 :   for (m = 0; m < N; m++) S1 = gadd(S1, gpow(gaddsg(m,x), al, prec));
    2246          42 :   Nx = gaddsg(N-1, x);
    2247          42 :   N2 = ginv(gsqr(Nx));
    2248          42 :   S2 = gen_0;
    2249        1092 :   for (j = k; j >= 2; j -= 2)
    2250             :   {
    2251        1050 :     GEN t = gsubgs(al, j), u = gmul(t, gaddgs(t, 1));
    2252        1050 :     u = gmul(gdivgs(u, j*(j+1)), gmul(S2, N2));
    2253        1050 :     S2 = gadd(gdivgs(bernfrac(j), j), u);
    2254             :   }
    2255          42 :   S2 = gadd(gmul(S2, gdiv(al, Nx)), ghalf);
    2256          42 :   S3 = gpow(Nx, al, prec);
    2257          42 :   if (!gequal1(s))
    2258          35 :     S2 = gmul(S3, gadd(gdiv(Nx, gaddsg(1, al)), S2));
    2259             :   else
    2260           7 :     S2 = gadd(glog(Nx, prec), gmul(S3, S2));
    2261          42 :   return gerepileupto(av, gsub(S1,S2));
    2262             : }
    2263             : 
    2264             : /***********************************************************************/
    2265             : /**                                                                   **/
    2266             : /**                    FONCTIONS POLYLOGARITHME                       **/
    2267             : /**                                                                   **/
    2268             : /***********************************************************************/
    2269             : 
    2270             : /* returns H_n = 1 + 1/2 + ... + 1/n, as a rational number (n "small") */
    2271             : static GEN
    2272          28 : Harmonic(long n)
    2273             : {
    2274          28 :   GEN h = gen_1;
    2275             :   long i;
    2276          28 :   for (i=2; i<=n; i++) h = gadd(h, mkfrac(gen_1, utoipos(i)));
    2277          28 :   return h;
    2278             : }
    2279             : 
    2280             : /* m >= 2. Validity domain contains | log |x| | < 5, best for |x| ~ 1.
    2281             :  * Li_m(x = e^z) = sum_{n >= 0} zeta(m-n) z^n / n!
    2282             :  *    with zeta(1) := H_{m-1} - log(-z) */
    2283             : static GEN
    2284          28 : cxpolylog(long m, GEN x, long prec)
    2285             : {
    2286             :   long li, n;
    2287             :   GEN z, h, q, s;
    2288             :   int real;
    2289             : 
    2290          28 :   if (gequal1(x)) return szeta(m,prec);
    2291             :   /* x real <= 1 ==> Li_m(x) real */
    2292          28 :   real = (typ(x) == t_REAL && (expo(x) < 0 || signe(x) <= 0));
    2293             : 
    2294          28 :   z = glog(x,prec);
    2295             :   /* n = 0 */
    2296          28 :   q = gen_1;
    2297          28 :   s = szeta(m,prec);
    2298          35 :   for (n=1; n < m-1; n++)
    2299             :   {
    2300           7 :     q = gdivgs(gmul(q,z),n);
    2301           7 :     s = gadd(s, gmul(szeta(m-n,prec), real? real_i(q): q));
    2302             :   }
    2303             :   /* n = m-1 */
    2304          28 :     q = gdivgs(gmul(q,z),n); /* multiply by "zeta(1)" */
    2305          28 :     h = gmul(q, gsub(Harmonic(m-1), glog(gneg_i(z),prec)));
    2306          28 :     s = gadd(s, real? real_i(h): h);
    2307             :   /* n = m */
    2308          28 :     q = gdivgs(gmul(q,z),m);
    2309          28 :     s = gadd(s, gmul(szeta(0,prec), real? real_i(q): q));
    2310             :   /* n = m+1 */
    2311          28 :     q = gdivgs(gmul(q,z),m+1);
    2312          28 :     s = gadd(s, gmul(szeta(-1,prec), real? real_i(q): q));
    2313             : 
    2314          28 :   z = gsqr(z); li = -(prec2nbits(prec)+1);
    2315             :   /* n = m+3, m+5, ...; note that zeta(- even integer) = 0 */
    2316         609 :   for(n = m+3;; n += 2)
    2317             :   {
    2318         609 :     GEN zet = szeta(m-n,prec);
    2319         609 :     q = divgunu(gmul(q,z), n-1);
    2320         609 :     s = gadd(s, gmul(zet, real? real_i(q): q));
    2321         609 :     if (gexpo(q) + expo(zet) < li) break;
    2322         581 :   }
    2323          28 :   return s;
    2324             : }
    2325             : 
    2326             : static GEN
    2327         203 : polylog(long m, GEN x, long prec)
    2328             : {
    2329             :   long l, e, i, G, sx;
    2330             :   pari_sp av, av1;
    2331             :   GEN X, Xn, z, p1, p2, y, res;
    2332             : 
    2333         203 :   if (m < 0) pari_err_DOMAIN("polylog", "index", "<", gen_0, stoi(m));
    2334         203 :   if (!m) return mkfrac(gen_m1,gen_2);
    2335         203 :   if (gequal0(x)) return gcopy(x);
    2336         203 :   if (m==1)
    2337             :   {
    2338          35 :     av = avma;
    2339          35 :     return gerepileupto(av, gneg(glog(gsub(gen_1,x), prec)));
    2340             :   }
    2341             : 
    2342         168 :   l = precision(x);
    2343         168 :   if (!l) l = prec; else prec = l;
    2344         168 :   res = cgetc(l); av = avma;
    2345         168 :   x = gtofp(x, l+EXTRAPRECWORD);
    2346         168 :   e = gexpo(gnorm(x));
    2347         168 :   if (!e || e == -1) {
    2348          28 :     y = cxpolylog(m,x,prec);
    2349          28 :     avma = av; return affc_fixlg(y, res);
    2350             :   }
    2351         140 :   X = (e > 0)? ginv(x): x;
    2352         140 :   G = -prec2nbits(l);
    2353         140 :   av1 = avma;
    2354         140 :   y = Xn = X;
    2355       66801 :   for (i=2; ; i++)
    2356             :   {
    2357       66801 :     Xn = gmul(X,Xn); p2 = gdiv(Xn,powuu(i,m));
    2358       66801 :     y = gadd(y,p2);
    2359       66801 :     if (gexpo(p2) <= G) break;
    2360             : 
    2361       66661 :     if (gc_needed(av1,1))
    2362             :     {
    2363           0 :       if(DEBUGMEM>1) pari_warn(warnmem,"polylog");
    2364           0 :       gerepileall(av1,2, &y, &Xn);
    2365             :     }
    2366       66661 :   }
    2367         140 :   if (e < 0) { avma = av; return affc_fixlg(y, res); }
    2368             : 
    2369          28 :   sx = gsigne(imag_i(x));
    2370          28 :   if (!sx)
    2371             :   {
    2372          28 :     if (m&1) sx = gsigne(gsub(gen_1, real_i(x)));
    2373          21 :     else     sx = - gsigne(real_i(x));
    2374             :   }
    2375          28 :   z = divri(mppi(l), mpfact(m-1)); setsigne(z, sx);
    2376          28 :   z = mkcomplex(gen_0, z);
    2377             : 
    2378          28 :   if (m == 2)
    2379             :   { /* same formula as below, written more efficiently */
    2380          21 :     y = gneg_i(y);
    2381          21 :     if (typ(x) == t_REAL && signe(x) < 0)
    2382           7 :       p1 = logr_abs(x);
    2383             :     else
    2384          14 :       p1 = gsub(glog(x,l), z);
    2385          21 :     p1 = gmul2n(gsqr(p1), -1); /* = (log(-x))^2 / 2 */
    2386             : 
    2387          21 :     p1 = gadd(p1, divru(sqrr(mppi(l)), 6));
    2388          21 :     p1 = gneg_i(p1);
    2389             :   }
    2390             :   else
    2391             :   {
    2392           7 :     GEN logx = glog(x,l), logx2 = gsqr(logx);
    2393           7 :     p1 = mkfrac(gen_m1,gen_2);
    2394          14 :     for (i=m-2; i>=0; i-=2)
    2395           7 :       p1 = gadd(szeta(m-i,l), gmul(p1,gdivgs(logx2,(i+1)*(i+2))));
    2396           7 :     if (m&1) p1 = gmul(logx,p1); else y = gneg_i(y);
    2397           7 :     p1 = gadd(gmul2n(p1,1), gmul(z,gpowgs(logx,m-1)));
    2398           7 :     if (typ(x) == t_REAL && signe(x) < 0) p1 = real_i(p1);
    2399             :   }
    2400          28 :   y = gadd(y,p1);
    2401          28 :   avma = av; return affc_fixlg(y, res);
    2402             : }
    2403             : 
    2404             : GEN
    2405          21 : dilog(GEN x, long prec)
    2406             : {
    2407          21 :   return gpolylog(2, x, prec);
    2408             : }
    2409             : 
    2410             : /* x a floating point number, t_REAL or t_COMPLEX of t_REAL */
    2411             : static GEN
    2412          42 : logabs(GEN x)
    2413             : {
    2414             :   GEN y;
    2415          42 :   if (typ(x) == t_COMPLEX)
    2416             :   {
    2417           7 :     y = logr_abs( cxnorm(x) );
    2418           7 :     shiftr_inplace(y, -1);
    2419             :   } else
    2420          35 :     y = logr_abs(x);
    2421          42 :   return y;
    2422             : }
    2423             : 
    2424             : static GEN
    2425          21 : polylogD(long m, GEN x, long flag, long prec)
    2426             : {
    2427             :   long k, l, fl, m2;
    2428             :   pari_sp av;
    2429             :   GEN p1, p2, y;
    2430             : 
    2431          21 :   if (gequal0(x)) return gcopy(x);
    2432          21 :   m2 = m&1;
    2433          21 :   if (gequal1(x) && m>=2) return m2? szeta(m,prec): gen_0;
    2434          21 :   av = avma; l = precision(x);
    2435          21 :   if (!l) { l = prec; x = gtofp(x,l); }
    2436          21 :   p1 = logabs(x);
    2437          21 :   k = signe(p1);
    2438          21 :   if (k > 0) { x = ginv(x); fl = !m2; } else { setabssign(p1); fl = 0; }
    2439             :   /* |x| <= 1, p1 = - log|x| >= 0 */
    2440          21 :   p2 = gen_1;
    2441          21 :   y = polylog(m,x,l);
    2442          21 :   y = m2? real_i(y): imag_i(y);
    2443          84 :   for (k=1; k<m; k++)
    2444             :   {
    2445          63 :     GEN t = polylog(m-k,x,l);
    2446          63 :     p2 = gdivgs(gmul(p2,p1), k); /* (-log|x|)^k / k! */
    2447          63 :     y = gadd(y, gmul(p2, m2? real_i(t): imag_i(t)));
    2448             :   }
    2449          21 :   if (m2)
    2450             :   {
    2451          14 :     if (!flag) p1 = negr( logabs(gsubsg(1,x)) ); else p1 = shiftr(p1,-1);
    2452          14 :     p2 = gdivgs(gmul(p2,p1), -m);
    2453          14 :     y = gadd(y, p2);
    2454             :   }
    2455          21 :   if (fl) y = gneg(y);
    2456          21 :   return gerepileupto(av, y);
    2457             : }
    2458             : 
    2459             : static GEN
    2460          14 : polylogP(long m, GEN x, long prec)
    2461             : {
    2462             :   long k, l, fl, m2;
    2463             :   pari_sp av;
    2464             :   GEN p1,y;
    2465             : 
    2466          14 :   if (gequal0(x)) return gcopy(x);
    2467          14 :   m2 = m&1;
    2468          14 :   if (gequal1(x) && m>=2) return m2? szeta(m,prec): gen_0;
    2469          14 :   av = avma; l = precision(x);
    2470          14 :   if (!l) { l = prec; x = gtofp(x,l); }
    2471          14 :   p1 = logabs(x);
    2472          14 :   if (signe(p1) > 0) { x = ginv(x); fl = !m2; setsigne(p1, -1); } else fl = 0;
    2473             :   /* |x| <= 1 */
    2474          14 :   y = polylog(m,x,l);
    2475          14 :   y = m2? real_i(y): imag_i(y);
    2476          14 :   if (m==1)
    2477             :   {
    2478           7 :     shiftr_inplace(p1, -1); /* log |x| / 2 */
    2479           7 :     y = gadd(y, p1);
    2480             :   }
    2481             :   else
    2482             :   { /* m >= 2, \sum_{0 <= k <= m} 2^k B_k/k! (log |x|)^k Li_{m-k}(x),
    2483             :        with Li_0(x) := -1/2 */
    2484             :     GEN u, t;
    2485           7 :     t = polylog(m-1,x,l);
    2486           7 :     u = gneg_i(p1); /* u = 2 B1 log |x| */
    2487           7 :     y = gadd(y, gmul(u, m2?real_i(t):imag_i(t)));
    2488           7 :     if (m > 2)
    2489             :     {
    2490             :       GEN p2;
    2491           7 :       shiftr_inplace(p1, 1); /* p1 = 2log|x| <= 0 */
    2492           7 :       mpbern(m>>1, l);
    2493           7 :       p1 = sqrr(p1);
    2494           7 :       p2 = shiftr(p1,-1);
    2495          21 :       for (k=2; k<m; k+=2)
    2496             :       {
    2497          14 :         if (k > 2) p2 = divgunu(gmul(p2,p1),k-1);
    2498             :         /* p2 = 2^k/k! log^k |x|*/
    2499          14 :         t = polylog(m-k,x,l);
    2500          14 :         u = gmul(p2, bernreal(k, prec));
    2501          14 :         y = gadd(y, gmul(u, m2?real_i(t):imag_i(t)));
    2502             :       }
    2503             :     }
    2504             :   }
    2505          14 :   if (fl) y = gneg(y);
    2506          14 :   return gerepileupto(av, y);
    2507             : }
    2508             : 
    2509             : GEN
    2510         140 : gpolylog(long m, GEN x, long prec)
    2511             : {
    2512             :   long i, n, v;
    2513         140 :   pari_sp av = avma;
    2514             :   GEN a, y, p1;
    2515             : 
    2516         140 :   if (m <= 0)
    2517             :   {
    2518           7 :     GEN t = mkpoln(2, gen_m1, gen_1); /* 1 - X */
    2519           7 :     p1 = pol_x(0);
    2520          28 :     for (i=2; i <= -m; i++)
    2521          21 :       p1 = RgX_shift_shallow(gadd(gmul(t,ZX_deriv(p1)), gmulsg(i,p1)), 1);
    2522           7 :     p1 = gdiv(p1, gpowgs(t,1-m));
    2523           7 :     return gerepileupto(av, poleval(p1,x));
    2524             :   }
    2525             : 
    2526         133 :   switch(typ(x))
    2527             :   {
    2528             :     case t_INT: case t_REAL: case t_FRAC: case t_COMPLEX: case t_QUAD:
    2529          84 :       return polylog(m,x,prec);
    2530             :     case t_POLMOD:
    2531           7 :       return gerepileupto(av, polylogvec(m, polmod_to_embed(x, prec), prec));
    2532           7 :     case t_INTMOD: case t_PADIC: pari_err_IMPL( "padic polylogarithm");
    2533             :     case t_VEC: case t_COL: case t_MAT:
    2534           7 :       return polylogvec(m, x, prec);
    2535             :     default:
    2536          28 :       av = avma; if (!(y = toser_i(x))) break;
    2537          21 :       if (!m) { avma = av; return mkfrac(gen_m1,gen_2); }
    2538          21 :       if (m==1) return gerepileupto(av, gneg( glog(gsub(gen_1,y),prec) ));
    2539          21 :       if (gequal0(y)) return gerepilecopy(av, y);
    2540          21 :       v = valp(y);
    2541          21 :       if (v < 0) pari_err_DOMAIN("polylog","valuation", "<", gen_0, x);
    2542          14 :       if (v > 0) {
    2543           7 :         n = (lg(y)-3 + v) / v;
    2544           7 :         a = zeroser(varn(y), lg(y)-2);
    2545          35 :         for (i=n; i>=1; i--)
    2546          28 :           a = gmul(y, gadd(a, powis(utoipos(i),-m)));
    2547             :       } else { /* v == 0 */
    2548           7 :         long vy = varn(y);
    2549           7 :         GEN a0 = polcoeff0(y, 0, -1), yprimeovery = gdiv(derivser(y), y);
    2550           7 :         a = gneg( glog(gsub(gen_1,y), prec) );
    2551          14 :         for (i=2; i<=m; i++)
    2552           7 :           a = gadd(gpolylog(i, a0, prec), integ(gmul(yprimeovery, a), vy));
    2553             :       }
    2554          14 :       return gerepileupto(av, a);
    2555             :   }
    2556           7 :   pari_err_TYPE("gpolylog",x);
    2557             :   return NULL; /* LCOV_EXCL_LINE */
    2558             : }
    2559             : 
    2560             : GEN
    2561         126 : polylog0(long m, GEN x, long flag, long prec)
    2562             : {
    2563         126 :   switch(flag)
    2564             :   {
    2565          84 :     case 0: return gpolylog(m,x,prec);
    2566          14 :     case 1: return polylogD(m,x,0,prec);
    2567           7 :     case 2: return polylogD(m,x,1,prec);
    2568          14 :     case 3: return polylogP(m,x,prec);
    2569           7 :     default: pari_err_FLAG("polylog");
    2570             :   }
    2571             :   return NULL; /* LCOV_EXCL_LINE */
    2572             : }
    2573             : 
    2574             : GEN
    2575        8932 : upper_to_cx(GEN x, long *prec)
    2576             : {
    2577        8932 :   long tx = typ(x), l;
    2578        8932 :   if (tx == t_QUAD) { x = quadtofp(x, *prec); tx = typ(x); }
    2579        8932 :   switch(tx)
    2580             :   {
    2581             :     case t_COMPLEX:
    2582        8911 :       if (gsigne(gel(x,2)) > 0) break; /*fall through*/
    2583             :     case t_REAL: case t_INT: case t_FRAC:
    2584          14 :       pari_err_DOMAIN("modular function", "Im(argument)", "<=", gen_0, x);
    2585             :     default:
    2586           7 :       pari_err_TYPE("modular function", x);
    2587             :   }
    2588        8911 :   l = precision(x); if (l) *prec = l;
    2589        8911 :   return x;
    2590             : }
    2591             : 
    2592             : /* sqrt(3)/2 */
    2593             : static GEN
    2594        1477 : sqrt32(long prec) { GEN z = sqrtr_abs(stor(3,prec)); setexpo(z, -1); return z; }
    2595             : /* exp(i x), x = k pi/12 */
    2596             : static GEN
    2597        2093 : e12(ulong k, long prec)
    2598             : {
    2599             :   int s, sPi, sPiov2;
    2600             :   GEN z, t;
    2601        2093 :   k %= 24;
    2602        2093 :   if (!k) return gen_1;
    2603        2093 :   if (k == 12) return gen_m1;
    2604        2093 :   if (k >12) { s = 1; k = 24 - k; } else s = 0; /* x -> 2pi - x */
    2605        2093 :   if (k > 6) { sPi = 1; k = 12 - k; } else sPi = 0; /* x -> pi  - x */
    2606        2093 :   if (k > 3) { sPiov2 = 1; k = 6 - k; } else sPiov2 = 0; /* x -> pi/2 - x */
    2607        2093 :   z = cgetg(3, t_COMPLEX);
    2608        2093 :   switch(k)
    2609             :   {
    2610         161 :     case 0: gel(z,1) = icopy(gen_1); gel(z,2) = gen_0; break;
    2611         728 :     case 1: t = gmul2n(addrs(sqrt32(prec), 1), -1);
    2612         728 :       gel(z,1) = sqrtr(t);
    2613         728 :       gel(z,2) = gmul2n(invr(gel(z,1)), -2); break;
    2614             : 
    2615         749 :     case 2: gel(z,1) = sqrt32(prec);
    2616         749 :             gel(z,2) = real2n(-1, prec); break;
    2617             : 
    2618         455 :     case 3: gel(z,1) = sqrtr_abs(real2n(-1,prec));
    2619         455 :             gel(z,2) = rcopy(gel(z,1)); break;
    2620             :   }
    2621        2093 :   if (sPiov2) swap(gel(z,1), gel(z,2));
    2622        2093 :   if (sPi) togglesign(gel(z,1));
    2623        2093 :   if (s)   togglesign(gel(z,2));
    2624        2093 :   return z;
    2625             : }
    2626             : /* z a t_FRAC */
    2627             : static GEN
    2628       10374 : eiPi_frac(GEN z, long prec)
    2629             : {
    2630             :   GEN n, d;
    2631             :   ulong q, r;
    2632       10374 :   n = gel(z,1);
    2633       10374 :   d = gel(z,2);
    2634       10374 :   q = udivui_rem(12, d, &r);
    2635       10374 :   if (!r) /* relatively frequent case */
    2636        2093 :     return e12(q * umodiu(n, 24), prec);
    2637        8281 :   n = centermodii(n, shifti(d,1), d);
    2638        8281 :   return expIr(divri(mulri(mppi(prec), n), d));
    2639             : }
    2640             : /* exp(i Pi z), z a t_INT or t_FRAC */
    2641             : static GEN
    2642       10990 : exp_IPiQ(GEN z, long prec)
    2643             : {
    2644       10990 :   if (typ(z) == t_INT) return mpodd(z)? gen_m1: gen_1;
    2645        2100 :   return eiPi_frac(z, prec);
    2646             : }
    2647             : /* z a t_COMPLEX */
    2648             : static GEN
    2649       17304 : exp_IPiC(GEN z, long prec)
    2650             : {
    2651       17304 :   GEN r, x = gel(z,1), y = gel(z,2);
    2652       17304 :   GEN pi, mpi = mppi(prec);
    2653       17304 :   togglesign(mpi); /* mpi = -Pi */
    2654       17304 :   r = gexp(gmul(mpi, y), prec);
    2655       17304 :   switch(typ(x))
    2656             :   {
    2657             :     case t_INT:
    2658         252 :       if (mpodd(x)) togglesign(r);
    2659         252 :       return r;
    2660             :     case t_FRAC:
    2661        8274 :       return gmul(r, eiPi_frac(x, prec));
    2662             :     default:
    2663        8778 :       pi = mpi; togglesign(mpi); /* pi = Pi */
    2664        8778 :       return gmul(r, expIr(gmul(pi, x)));
    2665             :   }
    2666             : }
    2667             : 
    2668             : static GEN
    2669          70 : qq(GEN x, long prec)
    2670             : {
    2671          70 :   long tx = typ(x);
    2672             :   GEN y;
    2673             : 
    2674          70 :   if (is_scalar_t(tx))
    2675             :   {
    2676          35 :     if (tx == t_PADIC) return x;
    2677          21 :     x = upper_to_cx(x, &prec);
    2678           7 :     return exp_IPiC(gmul2n(x,1), prec); /* e(x) */
    2679             :   }
    2680          35 :   if (! ( y = toser_i(x)) ) pari_err_TYPE("modular function", x);
    2681          35 :   return y;
    2682             : }
    2683             : 
    2684             : /* return (y * X^d) + x. Assume d > 0, x != 0, valp(x) = 0 */
    2685             : static GEN
    2686          21 : ser_addmulXn(GEN y, GEN x, long d)
    2687             : {
    2688          21 :   long i, lx, ly, l = valp(y) + d; /* > 0 */
    2689             :   GEN z;
    2690             : 
    2691          21 :   lx = lg(x);
    2692          21 :   ly = lg(y) + l; if (lx < ly) ly = lx;
    2693          21 :   if (l > lx-2) return gcopy(x);
    2694          21 :   z = cgetg(ly,t_SER);
    2695          21 :   for (i=2; i<=l+1; i++) gel(z,i) = gel(x,i);
    2696          21 :   for (   ; i < ly; i++) gel(z,i) = gadd(gel(x,i),gel(y,i-l));
    2697          21 :   z[1] = x[1]; return z;
    2698             : }
    2699             : 
    2700             : /* q a t_POL */
    2701             : static GEN
    2702          21 : inteta_pol(GEN q, long v, long l)
    2703             : {
    2704          21 :   pari_sp av = avma;
    2705             :   GEN qn, ps, y;
    2706             :   ulong vps, vqn, n;
    2707             : 
    2708          21 :   y = qn = ps = pol_1(0);
    2709          21 :   vps = vqn = 0;
    2710         238 :   for(n = 0;; n++)
    2711             :   { /* qn = q^n,  ps = (-1)^n q^(n(3n+1)/2),
    2712             :      * vps, vqn valuation of ps, qn HERE */
    2713         238 :     pari_sp av2 = avma;
    2714         238 :     ulong vt = vps + 2*vqn + v; /* valuation of t at END of loop body */
    2715             :     long k1, k2;
    2716             :     GEN t;
    2717         238 :     vqn += v; vps = vt + vqn; /* valuation of qn, ps at END of body */
    2718         238 :     k1 = l-2 + v - vt + 1;
    2719         238 :     k2 = k1 - vqn; /* = l-2 + v - vps + 1 */
    2720         238 :     if (k1 <= 0) break;
    2721         217 :     t = RgX_mul(q, RgX_sqr(qn));
    2722         217 :     t = RgXn_red_shallow(t, k1);
    2723         217 :     t = RgX_mul(ps,t);
    2724         217 :     t = RgXn_red_shallow(t, k1);
    2725         217 :     t = RgX_neg(t); /* t = (-1)^(n+1) q^(n(3n+1)/2 + 2n+1) */
    2726         217 :     t = gerepileupto(av2, t);
    2727         217 :     y = RgX_addmulXn_shallow(t, y, vt);
    2728         217 :     if (k2 <= 0) break;
    2729             : 
    2730         217 :     qn = RgX_mul(qn,q);
    2731         217 :     ps = RgX_mul(t,qn);
    2732         217 :     ps = RgXn_red_shallow(ps, k2);
    2733         217 :     y = RgX_addmulXn_shallow(ps, y, vps);
    2734             : 
    2735         217 :     if (gc_needed(av,1))
    2736             :     {
    2737           0 :       if(DEBUGMEM>1) pari_warn(warnmem,"eta, n = %ld", n);
    2738           0 :       gerepileall(av, 3, &y, &qn, &ps);
    2739             :     }
    2740         217 :   }
    2741          21 :   setvarn(y, varn(q)); return RgX_to_ser(y, l+v);
    2742             : }
    2743             : 
    2744             : static GEN
    2745       16417 : inteta(GEN q)
    2746             : {
    2747       16417 :   long tx = typ(q);
    2748             :   GEN ps, qn, y;
    2749             : 
    2750       16417 :   y = gen_1; qn = gen_1; ps = gen_1;
    2751       16417 :   if (tx==t_PADIC)
    2752             :   {
    2753          28 :     if (valp(q) <= 0) pari_err_DOMAIN("eta", "v_p(q)", "<=",gen_0,q);
    2754             :     for(;;)
    2755             :     {
    2756          77 :       GEN t = gneg_i(gmul(ps,gmul(q,gsqr(qn))));
    2757          77 :       y = gadd(y,t); qn = gmul(qn,q); ps = gmul(t,qn);
    2758          77 :       t = y;
    2759          77 :       y = gadd(y,ps); if (gequal(t,y)) return y;
    2760          56 :     }
    2761             :   }
    2762             : 
    2763       16389 :   if (tx == t_SER)
    2764             :   {
    2765             :     ulong vps, vqn;
    2766          35 :     long l = lg(q), v, n;
    2767             :     pari_sp av;
    2768             : 
    2769          35 :     v = valp(q); /* handle valuation separately to avoid overflow */
    2770          35 :     if (v <= 0) pari_err_DOMAIN("eta", "v_p(q)", "<=",gen_0,q);
    2771          28 :     y = ser2pol_i(q, l); /* t_SER inefficient when input has low degree */
    2772          28 :     n = degpol(y);
    2773          28 :     if (n == 1 || n < (l>>2)) return inteta_pol(y, v, l);
    2774             : 
    2775           7 :     q = leafcopy(q); av = avma;
    2776           7 :     setvalp(q, 0);
    2777           7 :     y = scalarser(gen_1, varn(q), l+v);
    2778           7 :     vps = vqn = 0;
    2779          14 :     for(n = 0;; n++)
    2780             :     { /* qn = q^n,  ps = (-1)^n q^(n(3n+1)/2) */
    2781          14 :       ulong vt = vps + 2*vqn + v;
    2782             :       long k;
    2783             :       GEN t;
    2784          14 :       t = gneg_i(gmul(ps,gmul(q,gsqr(qn))));
    2785             :       /* t = (-1)^(n+1) q^(n(3n+1)/2 + 2n+1) */
    2786          14 :       y = ser_addmulXn(t, y, vt);
    2787          14 :       qn = gmul(qn,q); ps = gmul(t,qn);
    2788          14 :       vqn += v; vps = vt + vqn;
    2789          14 :       k = l+v - vps; if (k <= 2) return y;
    2790             : 
    2791           7 :       y = ser_addmulXn(ps, y, vps);
    2792           7 :       setlg(q, k);
    2793           7 :       setlg(qn, k);
    2794           7 :       setlg(ps, k);
    2795           7 :       if (gc_needed(av,3))
    2796             :       {
    2797           0 :         if(DEBUGMEM>1) pari_warn(warnmem,"eta");
    2798           0 :         gerepileall(av, 3, &y, &qn, &ps);
    2799             :       }
    2800           7 :     }
    2801             :   }
    2802             :   {
    2803             :     long l; /* gcc -Wall */
    2804       16354 :     pari_sp av = avma;
    2805             : 
    2806       16354 :     l = -prec2nbits(precision(q));
    2807             :     for(;;)
    2808             :     {
    2809       61912 :       GEN t = gneg_i(gmul(ps,gmul(q,gsqr(qn))));
    2810             :       /* qn = q^n
    2811             :        * ps = (-1)^n q^(n(3n+1)/2)
    2812             :        * t = (-1)^(n+1) q^(n(3n+1)/2 + 2n+1) */
    2813       61912 :       y = gadd(y,t); qn = gmul(qn,q); ps = gmul(t,qn);
    2814       61912 :       y = gadd(y,ps);
    2815       61912 :       if (gexpo(ps)-gexpo(y) < l) return y;
    2816       45558 :       if (gc_needed(av,3))
    2817             :       {
    2818           0 :         if(DEBUGMEM>1) pari_warn(warnmem,"eta");
    2819           0 :         gerepileall(av, 3, &y, &qn, &ps);
    2820             :       }
    2821       45558 :     }
    2822             :   }
    2823             : }
    2824             : 
    2825             : GEN
    2826          70 : eta(GEN x, long prec)
    2827             : {
    2828          70 :   pari_sp av = avma;
    2829          70 :   GEN z = inteta( qq(x,prec) );
    2830          42 :   if (typ(z) == t_SER) return gerepilecopy(av, z);
    2831          14 :   return gerepileupto(av, z);
    2832             : }
    2833             : 
    2834             : /* s(h,k) = sum(n = 1, k-1, (n/k)*(frac(h*n/k) - 1/2))
    2835             :  * Knuth's algorithm. h integer, k integer > 0, (h,k) = 1 */
    2836             : GEN
    2837        6258 : sumdedekind_coprime(GEN h, GEN k)
    2838             : {
    2839        6258 :   pari_sp av = avma;
    2840             :   GEN s2, s1, p, pp;
    2841             :   long s;
    2842        6258 :   if (lgefint(k) == 3 && uel(k,2) <= (2*(ulong)LONG_MAX) / 3)
    2843             :   {
    2844        6251 :     ulong kk = k[2], hh = umodiu(h, kk);
    2845             :     long s1, s2;
    2846             :     GEN v;
    2847        6251 :     if (signe(k) < 0) { k = negi(k); hh = Fl_neg(hh, kk); }
    2848        6251 :     v = u_sumdedekind_coprime(hh, kk);
    2849        6251 :     s1 = v[1]; s2 = v[2];
    2850        6251 :     return gerepileupto(av, gdiv(addis(mulis(k,s1), s2), muluu(12, kk)));
    2851             :   }
    2852           7 :   s = 1;
    2853           7 :   s1 = gen_0; p = gen_1; pp = gen_0;
    2854           7 :   s2 = h = modii(h, k);
    2855          42 :   while (signe(h)) {
    2856          28 :     GEN r, nexth, a = dvmdii(k, h, &nexth);
    2857          28 :     if (is_pm1(h)) s2 = s == 1? addii(s2, p): subii(s2, p);
    2858          28 :     s1 = s == 1? addii(s1, a): subii(s1, a);
    2859          28 :     s = -s;
    2860          28 :     k = h; h = nexth;
    2861          28 :     r = addii(mulii(a,p), pp); pp = p; p = r;
    2862             :   }
    2863             :   /* at this point p = original k */
    2864           7 :   if (s == -1) s1 = subiu(s1, 3);
    2865           7 :   return gerepileupto(av, gdiv(addii(mulii(p,s1), s2), muliu(p,12)));
    2866             : }
    2867             : /* as above, for ulong arguments.
    2868             :  * k integer > 0, 0 <= h < k, (h,k) = 1. Returns [s1,s2] such that
    2869             :  * s(h,k) = (s2 + k s1) / (12k). Requires max(h + k/2, k) < LONG_MAX
    2870             :  * to avoid overflow, in particular k <= LONG_MAX * 2/3 is fine */
    2871             : GEN
    2872        6251 : u_sumdedekind_coprime(long h, long k)
    2873             : {
    2874        6251 :   long s = 1, s1 = 0, s2 = h, p = 1, pp = 0;
    2875       17675 :   while (h) {
    2876        5173 :     long r, nexth = k % h, a = k / h; /* a >= 1, a >= 2 if h = 1 */
    2877        5173 :     if (h == 1) s2 += p * s; /* occurs exactly once, last step */
    2878        5173 :     s1 += a * s;
    2879        5173 :     s = -s;
    2880        5173 :     k = h; h = nexth;
    2881        5173 :     r = a*p + pp; pp = p; p = r; /* p >= pp >= 0 */
    2882             :   }
    2883             :   /* in the above computation, p increases from 1 to original k,
    2884             :    * -k/2 <= s2 <= h + k/2, and |s1| <= k */
    2885        6251 :   if (s < 0) s1 -= 3; /* |s1| <= k+3 ? */
    2886             :   /* But in fact, |s2 + p s1| <= k^2 + 1/2 - 3k; if (s < 0), we have
    2887             :    * |s2| <= k/2 and it follows that |s1| < k here as well */
    2888             :   /* p = k; s(h,k) = (s2 + p s1)/12p. */
    2889        6251 :   return mkvecsmall2(s1, s2);
    2890             : }
    2891             : GEN
    2892          28 : sumdedekind(GEN h, GEN k)
    2893             : {
    2894          28 :   pari_sp av = avma;
    2895             :   GEN d;
    2896          28 :   if (typ(h) != t_INT) pari_err_TYPE("sumdedekind",h);
    2897          28 :   if (typ(k) != t_INT) pari_err_TYPE("sumdedekind",k);
    2898          28 :   d = gcdii(h,k);
    2899          28 :   if (!is_pm1(d))
    2900             :   {
    2901           7 :     h = diviiexact(h, d);
    2902           7 :     k = diviiexact(k, d);
    2903             :   }
    2904          28 :   return gerepileupto(av, sumdedekind_coprime(h,k));
    2905             : }
    2906             : 
    2907             : /* eta(x); assume Im x >> 0 (e.g. x in SL2's standard fundamental domain) */
    2908             : static GEN
    2909       17290 : eta_reduced(GEN x, long prec)
    2910             : {
    2911       17290 :   GEN z = exp_IPiC(gdivgs(x, 12), prec); /* e(x/24) */
    2912       17290 :   if (24 * gexpo(z) >= -prec2nbits(prec))
    2913       16333 :     z = gmul(z, inteta( gpowgs(z,24) ));
    2914       17290 :   return z;
    2915             : }
    2916             : 
    2917             : /* x = U.z (flag = 1), or x = U^(-1).z (flag = 0)
    2918             :  * Return [s,t] such that eta(z) = eta(x) * sqrt(s) * exp(I Pi t) */
    2919             : static GEN
    2920       17304 : eta_correction(GEN x, GEN U, long flag)
    2921             : {
    2922             :   GEN a,b,c,d, s,t;
    2923             :   long sc;
    2924       17304 :   a = gcoeff(U,1,1);
    2925       17304 :   b = gcoeff(U,1,2);
    2926       17304 :   c = gcoeff(U,2,1);
    2927       17304 :   d = gcoeff(U,2,2);
    2928             :   /* replace U by U^(-1) */
    2929       17304 :   if (flag) {
    2930        8932 :     swap(a,d);
    2931        8932 :     togglesign_safe(&b);
    2932        8932 :     togglesign_safe(&c);
    2933             :   }
    2934       17304 :   sc = signe(c);
    2935       17304 :   if (!sc) {
    2936       11074 :     if (signe(d) < 0) togglesign_safe(&b);
    2937       11074 :     s = gen_1;
    2938       11074 :     t = gdivgs(utoi(umodiu(b, 24)), 12);
    2939             :   } else {
    2940        6230 :     if (sc < 0) {
    2941        1722 :       togglesign_safe(&a);
    2942        1722 :       togglesign_safe(&b);
    2943        1722 :       togglesign_safe(&c);
    2944        1722 :       togglesign_safe(&d);
    2945             :     } /* now c > 0 */
    2946        6230 :     s = mulcxmI(gadd(gmul(c,x), d));
    2947        6230 :     t = gadd(gdiv(addii(a,d),muliu(c,12)), sumdedekind_coprime(negi(d),c));
    2948             :     /* correction : exp(I Pi (((a+d)/12c) + s(-d,c)) ) sqrt(-i(cx+d))  */
    2949             :   }
    2950       17304 :   return mkvec2(s, t);
    2951             : }
    2952             : 
    2953             : /* returns the true value of eta(x) for Im(x) > 0, using reduction to
    2954             :  * standard fundamental domain */
    2955             : GEN
    2956        8862 : trueeta(GEN x, long prec)
    2957             : {
    2958        8862 :   pari_sp av = avma;
    2959             :   GEN U, st, s, t;
    2960             : 
    2961        8862 :   if (!is_scalar_t(typ(x))) pari_err_TYPE("trueeta",x);
    2962        8862 :   x = upper_to_cx(x, &prec);
    2963        8862 :   x = cxredsl2(x, &U);
    2964        8862 :   st = eta_correction(x, U, 1);
    2965        8862 :   x = eta_reduced(x, prec);
    2966        8862 :   s = gel(st, 1);
    2967        8862 :   t = gel(st, 2);
    2968        8862 :   x = gmul(x, exp_IPiQ(t, prec));
    2969        8862 :   if (s != gen_1) x = gmul(x, gsqrt(s, prec));
    2970        8862 :   return gerepileupto(av, x);
    2971             : }
    2972             : 
    2973             : GEN
    2974          70 : eta0(GEN x, long flag,long prec)
    2975          70 : { return flag? trueeta(x,prec): eta(x,prec); }
    2976             : 
    2977             : /* eta(q) = 1 + \sum_{n>0} (-1)^n * (q^(n(3n-1)/2) + q^(n(3n+1)/2)) */
    2978             : static GEN
    2979           7 : ser_eta(long prec)
    2980             : {
    2981           7 :   GEN e = cgetg(prec+2, t_SER), ed = e+2;
    2982             :   long n, j;
    2983           7 :   e[1] = evalsigne(1)|_evalvalp(0)|evalvarn(0);
    2984           7 :   gel(ed,0) = gen_1;
    2985           7 :   for (n = 1; n < prec; n++) gel(ed,n) = gen_0;
    2986          49 :   for (n = 1, j = 0; n < prec; n++)
    2987             :   {
    2988             :     GEN s;
    2989          49 :     j += 3*n-2; /* = n*(3*n-1) / 2 */;
    2990          49 :     if (j >= prec) break;
    2991          42 :     s = odd(n)? gen_m1: gen_1;
    2992          42 :     gel(ed, j) = s;
    2993          42 :     if (j+n >= prec) break;
    2994          42 :     gel(ed, j+n) = s;
    2995             :   }
    2996           7 :   return e;
    2997             : }
    2998             : 
    2999             : static GEN
    3000         476 : coeffEu(GEN fa)
    3001             : {
    3002         476 :   pari_sp av = avma;
    3003         476 :   return gerepileuptoint(av, mului(65520, usumdivk_fact(fa,11)));
    3004             : }
    3005             : /* E12 = 1 + q*E/691 */
    3006             : static GEN
    3007           7 : ser_E(long prec)
    3008             : {
    3009           7 :   GEN e = cgetg(prec+2, t_SER), ed = e+2;
    3010           7 :   GEN F = vecfactoru_i(2, prec); /* F[n] = factoru(n+1) */
    3011             :   long n;
    3012           7 :   e[1] = evalsigne(1)|_evalvalp(0)|evalvarn(0);
    3013           7 :   gel(ed,0) = utoipos(65520);
    3014           7 :   for (n = 1; n < prec; n++) gel(ed,n) = coeffEu(gel(F,n));
    3015           7 :   return e;
    3016             : }
    3017             : /* j = E12/Delta + 432000/691, E12 = 1 + q*E/691 */
    3018             : static GEN
    3019           7 : ser_j2(long prec, long v)
    3020             : {
    3021           7 :   pari_sp av = avma;
    3022           7 :   GEN iD = gpowgs(ginv(ser_eta(prec)), 24); /* q/Delta */
    3023           7 :   GEN J = gmul(ser_E(prec), iD);
    3024           7 :   setvalp(iD,-1); /* now 1/Delta */
    3025           7 :   J = gadd(gdivgs(J, 691), iD);
    3026           7 :   J = gerepileupto(av, J);
    3027           7 :   if (prec > 1) gel(J,3) = utoipos(744);
    3028           7 :   setvarn(J,v); return J;
    3029             : }
    3030             : 
    3031             : /* j(q) = \sum_{n >= -1} c(n)q^n,
    3032             :  * \sum_{n = -1}^{N-1} c(n) (-10n \sigma_3(N-n) + 21 \sigma_5(N-n))
    3033             :  * = c(N) (N+1)/24 */
    3034             : static GEN
    3035          14 : ser_j(long prec, long v)
    3036             : {
    3037             :   GEN j, J, S3, S5, F;
    3038             :   long i, n;
    3039          14 :   if (prec > 64) return ser_j2(prec, v);
    3040           7 :   S3 = cgetg(prec+1, t_VEC);
    3041           7 :   S5 = cgetg(prec+1,t_VEC);
    3042           7 :   F = vecfactoru_i(1, prec);
    3043          35 :   for (n = 1; n <= prec; n++)
    3044             :   {
    3045          28 :     GEN fa = gel(F,n);
    3046          28 :     gel(S3,n) = mului(10, usumdivk_fact(fa,3));
    3047          28 :     gel(S5,n) = mului(21, usumdivk_fact(fa,5));
    3048             :   }
    3049           7 :   J = cgetg(prec+2, t_SER),
    3050           7 :   J[1] = evalvarn(v)|evalsigne(1)|evalvalp(-1);
    3051           7 :   j = J+3;
    3052           7 :   gel(j,-1) = gen_1;
    3053           7 :   gel(j,0) = utoipos(744);
    3054           7 :   gel(j,1) = utoipos(196884);
    3055          21 :   for (n = 2; n < prec; n++)
    3056             :   {
    3057          14 :     pari_sp av = avma;
    3058          14 :     GEN c, s3 = gel(S3,n+1), s5 = gel(S5,n+1);
    3059          14 :     c = addii(s3, s5);
    3060          49 :     for (i = 0; i < n; i++)
    3061             :     {
    3062          35 :       s3 = gel(S3,n-i); s5 = gel(S5,n-i);
    3063          35 :       c = addii(c, mulii(gel(j,i), subii(s5, mului(i,s3))));
    3064             :     }
    3065          14 :     gel(j,n) = gerepileuptoint(av, diviuexact(muliu(c,24), n+1));
    3066             :   }
    3067           7 :   return J;
    3068             : }
    3069             : 
    3070             : GEN
    3071          42 : jell(GEN x, long prec)
    3072             : {
    3073          42 :   long tx = typ(x);
    3074          42 :   pari_sp av = avma;
    3075             :   GEN q, h, U;
    3076             : 
    3077          42 :   if (!is_scalar_t(tx))
    3078             :   {
    3079             :     long v;
    3080          21 :     if (gequalX(x)) return ser_j(precdl, varn(x));
    3081          21 :     q = toser_i(x); if (!q) pari_err_TYPE("ellj",x);
    3082          14 :     v = fetch_var_higher();
    3083          14 :     h = ser_j(lg(q)-2, v);
    3084          14 :     h = gsubst(h, v, q);
    3085          14 :     delete_var(); return gerepileupto(av, h);
    3086             :   }
    3087          21 :   if (tx == t_PADIC)
    3088             :   {
    3089           7 :     GEN p2, p1 = gdiv(inteta(gsqr(x)), inteta(x));
    3090           7 :     p1 = gmul2n(gsqr(p1),1);
    3091           7 :     p1 = gmul(x,gpowgs(p1,12));
    3092           7 :     p2 = gaddsg(768,gadd(gsqr(p1),gdivsg(4096,p1)));
    3093           7 :     p1 = gmulsg(48,p1);
    3094           7 :     return gerepileupto(av, gadd(p2,p1));
    3095             :   }
    3096             :   /* Let h = Delta(2x) / Delta(x), then j(x) = (1 + 256h)^3 / h */
    3097          14 :   x = upper_to_cx(x, &prec);
    3098           7 :   x = cxredsl2(x, &U); /* forget about Ua : j has weight 0 */
    3099             :   { /* cf eta_reduced, raised to power 24
    3100             :      * Compute
    3101             :      *   t = (inteta(q(2x)) / inteta(q(x))) ^ 24;
    3102             :      * then
    3103             :      *   h = t * (q(2x) / q(x) = t * q(x);
    3104             :      * but inteta(q) costly and useless if expo(q) << 1  => inteta(q) = 1.
    3105             :      * log_2 ( exp(-2Pi Im tau) ) < -prec2nbits(prec)
    3106             :      * <=> Im tau > prec2nbits(prec) * log(2) / 2Pi */
    3107           7 :     long C = (long)prec2nbits_mul(prec, LOG2/(2*M_PI));
    3108           7 :     q = exp_IPiC(gmul2n(x,1), prec); /* e(x) */
    3109           7 :     if (gcmpgs(gel(x,2), C) > 0) /* eta(q(x)) = 1 : no need to compute q(2x) */
    3110           0 :       h = q;
    3111             :     else
    3112             :     {
    3113           7 :       GEN t = gdiv(inteta(gsqr(q)), inteta(q));
    3114           7 :       h = gmul(q, gpowgs(t, 24));
    3115             :     }
    3116             :   }
    3117             :   /* real_1 important ! gaddgs(, 1) could increase the accuracy ! */
    3118           7 :   return gerepileupto(av, gdiv(gpowgs(gadd(gmul2n(h,8), real_1(prec)), 3), h));
    3119             : }
    3120             : 
    3121             : static GEN
    3122        8372 : to_form(GEN a, GEN w, GEN C)
    3123        8372 : { return mkvec3(a, w, diviiexact(C, a)); }
    3124             : static GEN
    3125        8372 : form_to_quad(GEN f, GEN sqrtD)
    3126             : {
    3127        8372 :   long a = itos(gel(f,1)), a2 = a << 1;
    3128        8372 :   GEN b = gel(f,2);
    3129        8372 :   return mkcomplex(gdivgs(b, -a2), gdivgs(sqrtD, a2));
    3130             : }
    3131             : static GEN
    3132        8372 : eta_form(GEN f, GEN sqrtD, GEN *s_t, long prec)
    3133             : {
    3134        8372 :   GEN U, t = form_to_quad(redimagsl2(f, &U), sqrtD);
    3135        8372 :   *s_t = eta_correction(t, U, 0);
    3136        8372 :   return eta_reduced(t, prec);
    3137             : }
    3138             : 
    3139             : /* eta(t/p)eta(t/q) / (eta(t)eta(t/pq)), t = (-w + sqrt(D)) / 2a */
    3140             : GEN
    3141        2093 : double_eta_quotient(GEN a, GEN w, GEN D, long p, long q, GEN pq, GEN sqrtD)
    3142             : {
    3143        2093 :   GEN C = shifti(subii(sqri(w), D), -2);
    3144             :   GEN d, t, z, zp, zq, zpq, s_t, s_tp, s_tpq, s, sp, spq;
    3145        2093 :   long prec = realprec(sqrtD);
    3146             : 
    3147        2093 :   z = eta_form(to_form(a, w, C), sqrtD, &s_t, prec);
    3148        2093 :   s = gel(s_t, 1);
    3149        2093 :   zp = eta_form(to_form(mului(p, a), w, C), sqrtD, &s_tp, prec);
    3150        2093 :   sp = gel(s_tp, 1);
    3151        2093 :   zpq = eta_form(to_form(mulii(pq, a), w, C), sqrtD, &s_tpq, prec);
    3152        2093 :   spq = gel(s_tpq, 1);
    3153        2093 :   if (p == q) {
    3154           0 :     z = gdiv(gsqr(zp), gmul(z, zpq));
    3155           0 :     t = gsub(gmul2n(gel(s_tp,2), 1),
    3156           0 :              gadd(gel(s_t,2), gel(s_tpq,2)));
    3157           0 :     if (sp != gen_1) z = gmul(z, sp);
    3158             :   } else {
    3159             :     GEN s_tq, sq;
    3160        2093 :     zq = eta_form(to_form(mului(q, a), w, C), sqrtD, &s_tq, prec);
    3161        2093 :     sq = gel(s_tq, 1);
    3162        2093 :     z = gdiv(gmul(zp, zq), gmul(z, zpq));
    3163        4186 :     t = gsub(gadd(gel(s_tp,2), gel(s_tq,2)),
    3164        4186 :              gadd(gel(s_t,2), gel(s_tpq,2)));
    3165        2093 :     if (sp != gen_1) z = gmul(z, gsqrt(sp, prec));
    3166        2093 :     if (sq != gen_1) z = gmul(z, gsqrt(sq, prec));
    3167             :   }
    3168        2093 :   d = NULL;
    3169        2093 :   if (s != gen_1) d = gsqrt(s, prec);
    3170        2093 :   if (spq != gen_1) {
    3171        2065 :     GEN x = gsqrt(spq, prec);
    3172        2065 :     d = d? gmul(d, x): x;
    3173             :   }
    3174        2093 :   if (d) z = gdiv(z, d);
    3175        2093 :   return gmul(z, exp_IPiQ(t, prec));
    3176             : }
    3177             : 
    3178             : typedef struct { GEN u; long v, t; } cxanalyze_t;
    3179             : 
    3180             : /* Check whether a t_COMPLEX, t_REAL or t_INT z != 0 can be written as
    3181             :  * z = u * 2^(v/2) * exp(I Pi/4 t), u > 0, v = 0,1 and -3 <= t <= 4.
    3182             :  * Allow z t_INT/t_REAL to simplify handling of eta_correction() output */
    3183             : static int
    3184          70 : cxanalyze(cxanalyze_t *T, GEN z)
    3185             : {
    3186             :   GEN a, b;
    3187             :   long ta, tb;
    3188             : 
    3189          70 :   T->v = 0;
    3190          70 :   if (is_intreal_t(typ(z)))
    3191             :   {
    3192          63 :     T->u = mpabs_shallow(z);
    3193          63 :     T->t = signe(z) < 0? 4: 0;
    3194          63 :     return 1;
    3195             :   }
    3196           7 :   a = gel(z,1); ta = typ(a);
    3197           7 :   b = gel(z,2); tb = typ(b);
    3198             : 
    3199           7 :   T->t = 0;
    3200           7 :   if (ta == t_INT && !signe(a))
    3201             :   {
    3202           0 :     T->u = R_abs_shallow(b);
    3203           0 :     T->t = gsigne(b) < 0? -2: 2;
    3204           0 :     return 1;
    3205             :   }
    3206           7 :   if (tb == t_INT && !signe(b))
    3207             :   {
    3208           0 :     T->u = R_abs_shallow(a);
    3209           0 :     T->t = gsigne(a) < 0? 4: 0;
    3210           0 :     return 1;
    3211             :   }
    3212           7 :   if (ta != tb || ta == t_REAL) { T->u = z; return 0; }
    3213             :   /* a,b both non zero, both t_INT or t_FRAC */
    3214           7 :   if (ta == t_INT)
    3215             :   {
    3216           7 :     if (!absequalii(a, b)) return 0;
    3217           7 :     T->u = absi_shallow(a);
    3218           7 :     T->v = 1;
    3219           7 :     if (signe(a) == signe(b))
    3220           0 :     { T->t = signe(a) < 0? -3: 1; }
    3221             :     else
    3222           7 :     { T->t = signe(a) < 0? 3: -1; }
    3223             :   }
    3224             :   else
    3225             :   {
    3226           0 :     if (!absequalii(gel(a,2), gel(b,2)) || !absequalii(gel(a,1),gel(b,1)))
    3227           0 :       return 0;
    3228           0 :     T->u = absfrac_shallow(a);
    3229           0 :     T->v = 1;
    3230           0 :     a = gel(a,1);
    3231           0 :     b = gel(b,1);
    3232           0 :     if (signe(a) == signe(b))
    3233           0 :     { T->t = signe(a) < 0? -3: 1; }
    3234             :     else
    3235           0 :     { T->t = signe(a) < 0? 3: -1; }
    3236             :   }
    3237           7 :   return 1;
    3238             : }
    3239             : 
    3240             : /* z * sqrt(st_b) / sqrt(st_a) exp(I Pi (t + t0)). Assume that
    3241             :  * sqrt2 = gsqrt(gen_2, prec) or NULL */
    3242             : static GEN
    3243          35 : apply_eta_correction(GEN z, GEN st_a, GEN st_b, GEN t0, GEN sqrt2, long prec)
    3244             : {
    3245          35 :   GEN t, s_a = gel(st_a, 1), s_b = gel(st_b, 1);
    3246             :   cxanalyze_t Ta, Tb;
    3247             :   int ca, cb;
    3248             : 
    3249          35 :   t = gsub(gel(st_b,2), gel(st_a,2));
    3250          35 :   if (t0 != gen_0) t = gadd(t, t0);
    3251          35 :   ca = cxanalyze(&Ta, s_a);
    3252          35 :   cb = cxanalyze(&Tb, s_b);
    3253          35 :   if (ca || cb)
    3254          35 :   { /* compute sqrt(s_b) / sqrt(s_a) in a more efficient way:
    3255             :      * sb = ub sqrt(2)^vb exp(i Pi/4 tb) */
    3256          35 :     GEN u = gdiv(Tb.u,Ta.u);
    3257          35 :     switch(Tb.v - Ta.v)
    3258             :     {
    3259           0 :       case -1: u = gmul2n(u,-1); /* fall through: write 1/sqrt2 = sqrt2/2 */
    3260           7 :       case 1: u = gmul(u, sqrt2? sqrt2: sqrtr_abs(real2n(1, prec)));
    3261             :     }
    3262          35 :     if (!isint1(u)) z = gmul(z, gsqrt(u, prec));
    3263          35 :     t = gadd(t, gmul2n(stoi(Tb.t - Ta.t), -3));
    3264             :   }
    3265             :   else
    3266             :   {
    3267           0 :     z = gmul(z, gsqrt(s_b, prec));
    3268           0 :     z = gdiv(z, gsqrt(s_a, prec));
    3269             :   }
    3270          35 :   return gmul(z, exp_IPiQ(t, prec));
    3271             : }
    3272             : 
    3273             : /* sqrt(2) eta(2x) / eta(x) */
    3274             : GEN
    3275           7 : weberf2(GEN x, long prec)
    3276             : {
    3277           7 :   pari_sp av = avma;
    3278             :   GEN z, sqrt2, a,b, Ua,Ub, st_a,st_b;
    3279             : 
    3280           7 :   x = upper_to_cx(x, &prec);
    3281           7 :   a = cxredsl2(x, &Ua);
    3282           7 :   b = cxredsl2(gmul2n(x,1), &Ub);
    3283           7 :   if (gequal(a,b)) /* not infrequent */
    3284           0 :     z = gen_1;
    3285             :   else
    3286           7 :     z = gdiv(eta_reduced(b,prec), eta_reduced(a,prec));
    3287           7 :   st_a = eta_correction(a, Ua, 1);
    3288           7 :   st_b = eta_correction(b, Ub, 1);
    3289           7 :   sqrt2 = sqrtr_abs(real2n(1, prec));
    3290           7 :   z = apply_eta_correction(z, st_a, st_b, gen_0, sqrt2, prec);
    3291           7 :   return gerepileupto(av, gmul(z, sqrt2));
    3292             : }
    3293             : 
    3294             : /* eta(x/2) / eta(x) */
    3295             : GEN
    3296          14 : weberf1(GEN x, long prec)
    3297             : {
    3298          14 :   pari_sp av = avma;
    3299             :   GEN z, a,b, Ua,Ub, st_a,st_b;
    3300             : 
    3301          14 :   x = upper_to_cx(x, &prec);
    3302          14 :   a = cxredsl2(x, &Ua);
    3303          14 :   b = cxredsl2(gmul2n(x,-1), &Ub);
    3304          14 :   if (gequal(a,b)) /* not infrequent */
    3305           0 :     z = gen_1;
    3306             :   else
    3307          14 :     z = gdiv(eta_reduced(b,prec), eta_reduced(a,prec));
    3308          14 :   st_a = eta_correction(a, Ua, 1);
    3309          14 :   st_b = eta_correction(b, Ub, 1);
    3310          14 :   z = apply_eta_correction(z, st_a, st_b, gen_0, NULL, prec);
    3311          14 :   return gerepileupto(av, z);
    3312             : }
    3313             : /* exp(-I*Pi/24) * eta((x+1)/2) / eta(x) */
    3314             : GEN
    3315          14 : weberf(GEN x, long prec)
    3316             : {
    3317          14 :   pari_sp av = avma;
    3318             :   GEN z, t0, a,b, Ua,Ub, st_a,st_b;
    3319          14 :   x = upper_to_cx(x, &prec);
    3320          14 :   a = cxredsl2(x, &Ua);
    3321          14 :   b = cxredsl2(gmul2n(gaddgs(x,1),-1), &Ub);
    3322          14 :   if (gequal(a,b)) /* not infrequent */
    3323           7 :     z = gen_1;
    3324             :   else
    3325           7 :     z = gdiv(eta_reduced(b,prec), eta_reduced(a,prec));
    3326          14 :   st_a = eta_correction(a, Ua, 1);
    3327          14 :   st_b = eta_correction(b, Ub, 1);
    3328          14 :   t0 = mkfrac(gen_m1, utoipos(24));
    3329          14 :   z = apply_eta_correction(z, st_a, st_b, t0, NULL, prec);
    3330          14 :   if (typ(z) == t_COMPLEX && isexactzero(real_i(x)))
    3331           0 :     z = gerepilecopy(av, gel(z,1));
    3332             :   else
    3333          14 :     z = gerepileupto(av, z);
    3334          14 :   return z;
    3335             : }
    3336             : GEN
    3337          35 : weber0(GEN x, long flag,long prec)
    3338             : {
    3339          35 :   switch(flag)
    3340             :   {
    3341          14 :     case 0: return weberf(x,prec);
    3342          14 :     case 1: return weberf1(x,prec);
    3343           7 :     case 2: return weberf2(x,prec);
    3344           0 :     default: pari_err_FLAG("weber");
    3345             :   }
    3346             :   return NULL; /* LCOV_EXCL_LINE */
    3347             : }
    3348             : 
    3349             : /* check |q| < 1 */
    3350             : static GEN
    3351          21 : check_unit_disc(const char *fun, GEN q, long prec)
    3352             : {
    3353          21 :   GEN Q = gtofp(q, prec), Qlow;
    3354          21 :   Qlow = (prec > LOWDEFAULTPREC)? gtofp(Q,LOWDEFAULTPREC): Q;
    3355          21 :   if (gcmp(gnorm(Qlow), gen_1) >= 0)
    3356           0 :     pari_err_DOMAIN(fun, "abs(q)", ">=", gen_1, q);
    3357          21 :   return Q;
    3358             : }
    3359             : 
    3360             : GEN
    3361          14 : theta(GEN q, GEN z, long prec)
    3362             : {
    3363             :   long l, n;
    3364          14 :   pari_sp av = avma, av2;
    3365             :   GEN s, c, snz, cnz, s2z, c2z, ps, qn, y, zy, ps2, k, zold;
    3366             : 
    3367          14 :   l = precision(q);
    3368          14 :   n = precision(z); if (n && n < l) l = n;
    3369          14 :   if (l) prec = l;
    3370          14 :   z = gtofp(z, prec);
    3371          14 :   q = check_unit_disc("theta", q, prec);
    3372          14 :   zold = NULL; /* gcc -Wall */
    3373          14 :   zy = imag_i(z);
    3374          14 :   if (gequal0(zy)) k = gen_0;
    3375             :   else
    3376             :   {
    3377           7 :     GEN lq = glog(q,prec); k = roundr(divrr(zy, real_i(lq)));
    3378           7 :     if (signe(k)) { zold = z; z = gadd(z, mulcxmI(gmul(lq,k))); }
    3379             :   }
    3380          14 :   qn = gen_1;
    3381          14 :   ps2 = gsqr(q);
    3382          14 :   ps = gneg_i(ps2);
    3383          14 :   gsincos(z, &s, &c, prec);
    3384          14 :   s2z = gmul2n(gmul(s,c),1); /* sin 2z */
    3385          14 :   c2z = gsubgs(gmul2n(gsqr(c),1), 1); /* cos 2z */
    3386          14 :   snz = s;
    3387          14 :   cnz = c; y = s;
    3388          14 :   av2 = avma;
    3389         273 :   for (n = 3;; n += 2)
    3390             :   {
    3391             :     long e;
    3392         273 :     s = gadd(gmul(snz, c2z), gmul(cnz,s2z));
    3393         273 :     qn = gmul(qn,ps);
    3394         273 :     y = gadd(y, gmul(qn, s));
    3395         273 :     e = gexpo(s); if (e < 0) e = 0;
    3396         273 :     if (gexpo(qn) + e < -prec2nbits(prec)) break;
    3397             : 
    3398         259 :     ps = gmul(ps,ps2);
    3399         259 :     c = gsub(gmul(cnz, c2z), gmul(snz,s2z));
    3400         259 :     snz = s; /* sin nz */
    3401         259 :     cnz = c; /* cos nz */
    3402         259 :     if (gc_needed(av,2))
    3403             :     {
    3404           0 :       if (DEBUGMEM>1) pari_warn(warnmem,"theta (n = %ld)", n);
    3405           0 :       gerepileall(av2, 5, &snz, &cnz, &ps, &qn, &y);
    3406             :     }
    3407         259 :   }
    3408          14 :   if (signe(k))
    3409             :   {
    3410           7 :     y = gmul(y, gmul(powgi(q,sqri(k)),
    3411             :                      gexp(gmul(mulcxI(zold),shifti(k,1)), prec)));
    3412           7 :     if (mod2(k)) y = gneg_i(y);
    3413             :   }
    3414          14 :   return gerepileupto(av, gmul(y, gmul2n(gsqrt(gsqrt(q,prec),prec),1)));
    3415             : }
    3416             : 
    3417             : GEN
    3418           7 : thetanullk(GEN q, long k, long prec)
    3419             : {
    3420             :   long l, n;
    3421           7 :   pari_sp av = avma;
    3422             :   GEN p1, ps, qn, y, ps2;
    3423             : 
    3424           7 :   if (k < 0)
    3425           0 :     pari_err_DOMAIN("thetanullk", "k", "<", gen_0, stoi(k));
    3426           7 :   l = precision(q);
    3427           7 :   if (l) prec = l;
    3428           7 :   q = check_unit_disc("thetanullk", q, prec);
    3429             : 
    3430           7 :   if (!(k&1)) { avma = av; return gen_0; }
    3431           7 :   qn = gen_1;
    3432           7 :   ps2 = gsqr(q);
    3433           7 :   ps = gneg_i(ps2);
    3434           7 :   y = gen_1;
    3435         287 :   for (n = 3;; n += 2)
    3436             :   {
    3437             :     GEN t;
    3438         287 :     qn = gmul(qn,ps);
    3439         287 :     ps = gmul(ps,ps2);
    3440         287 :     t = gmul(qn, powuu(n, k)); y = gadd(y, t);
    3441         287 :     if (gexpo(t) < -prec2nbits(prec)) break;
    3442         280 :   }
    3443           7 :   p1 = gmul2n(gsqrt(gsqrt(q,prec),prec),1);
    3444           7 :   if (k&2) y = gneg_i(y);
    3445           7 :   return gerepileupto(av, gmul(p1, y));
    3446             : }
    3447             : 
    3448             : /* q2 = q^2 */
    3449             : static GEN
    3450       13440 : vecthetanullk_loop(GEN q2, long k, long prec)
    3451             : {
    3452       13440 :   GEN ps, qn = gen_1, y = const_vec(k, gen_1);
    3453       13440 :   pari_sp av = avma;
    3454       13440 :   const long bit = prec2nbits(prec);
    3455             :   long i, n;
    3456             : 
    3457       13440 :   ps = gneg_i(q2);
    3458       87444 :   for (n = 3;; n += 2)
    3459             :   {
    3460       87444 :     GEN t = NULL/*-Wall*/, P = utoipos(n), N2 = sqru(n);
    3461       87444 :     qn = gmul(qn,ps);
    3462       87444 :     ps = gmul(ps,q2);
    3463      349776 :     for (i = 1; i <= k; i++)
    3464             :     {
    3465      262332 :       t = gmul(qn, P); gel(y,i) = gadd(gel(y,i), t);
    3466      262332 :       P = mulii(P, N2);
    3467             :     }
    3468      100884 :     if (gexpo(t) < -bit) return y;
    3469       74004 :     if (gc_needed(av,2))
    3470             :     {
    3471           0 :       if (DEBUGMEM>1) pari_warn(warnmem,"vecthetanullk_loop, n = %ld",n);
    3472           0 :       gerepileall(av, 3, &qn, &ps, &y);
    3473             :     }
    3474       74004 :   }
    3475             : }
    3476             : /* [d^i theta/dz^i(q, 0), i = 1, 3, .., 2*k - 1] */
    3477             : GEN
    3478           0 : vecthetanullk(GEN q, long k, long prec)
    3479             : {
    3480           0 :   long i, l = precision(q);
    3481           0 :   pari_sp av = avma;
    3482             :   GEN p1, y;
    3483             : 
    3484           0 :   if (l) prec = l;
    3485           0 :   q = check_unit_disc("vecthetanullk", q, prec);
    3486           0 :   y = vecthetanullk_loop(gsqr(q), k, prec);
    3487             : 
    3488           0 :   p1 = gmul2n(gsqrt(gsqrt(q,prec),prec),1);
    3489           0 :   for (i = 2; i <= k; i+=2) gel(y,i) = gneg_i(gel(y,i));
    3490           0 :   return gerepileupto(av, gmul(p1, y));
    3491             : }
    3492             : 
    3493             : /* [d^i theta/dz^i(q, 0), i = 1, 3, .., 2*k - 1], q = exp(2iPi tau) */
    3494             : GEN
    3495           0 : vecthetanullk_tau(GEN tau, long k, long prec)
    3496             : {
    3497           0 :   long i, l = precision(tau);
    3498           0 :   pari_sp av = avma;
    3499             :   GEN p1, q4, y;
    3500             : 
    3501           0 :   if (l) prec = l;
    3502           0 :   if (typ(tau) != t_COMPLEX || gsigne(gel(tau,2)) <= 0)
    3503           0 :     pari_err_DOMAIN("vecthetanullk_tau", "imag(tau)", "<=", gen_0, tau);
    3504             : 
    3505           0 :   q4 = expIxy(Pi2n(-1, prec), tau, prec); /* q^(1/4) */
    3506           0 :   y = vecthetanullk_loop(gpowgs(q4,8), k, prec);
    3507           0 :   p1 = gmul2n(q4,1);
    3508           0 :   for (i = 2; i <= k; i+=2) gel(y,i) = gneg_i(gel(y,i));
    3509           0 :   return gerepileupto(av, gmul(p1, y));
    3510             : }
    3511             : 
    3512             : /* -theta^(3)(tau/2) / theta^(1)(tau/2). Assume that Im tau > 0 */
    3513             : GEN
    3514       13440 : trueE2(GEN tau, long prec)
    3515             : {
    3516       13440 :   long l = precision(tau);
    3517       13440 :   pari_sp av = avma;
    3518             :   GEN q2, y;
    3519       13440 :   if (l) prec = l;
    3520       13440 :   q2 = expIxy(Pi2n(1, prec), tau, prec);
    3521       13440 :   y = vecthetanullk_loop(q2, 3, prec);
    3522       13440 :   return gerepileupto(av, gdiv(gel(y,2), gel(y,1)));
    3523             : }
    3524             : 
    3525             : /* Lambert W function : solution x of x*exp(x)=y, using Newton. y >= 0 t_REAL.
    3526             :  * Good for low accuracy: the precision remains constant. Not memory-clean */
    3527             : static GEN
    3528         196 : mplambertW0(GEN y)
    3529             : {
    3530         196 :   long bitprec = bit_prec(y) - 2;
    3531             :   GEN x, tmp;
    3532             : 
    3533         196 :   x = mplog(addrs(y,1));
    3534             :   do {
    3535        1135 :    tmp = x;
    3536             :    /* f(x) = log(x)+x-log(y), f'(x) = (1+1/x)
    3537             :     * x *= (1-log(x/y))/(x+1) */
    3538        1135 :     x = mulrr(x, divrr(subsr(1, mplog(divrr(x,y))), addrs(x,1)));
    3539        1135 :   } while (expo(tmp) - expo(subrr(x,tmp)) < bitprec);
    3540         196 :   return x;
    3541             : }
    3542             : 
    3543             : /* Lambert W function using Newton, increasing prec */
    3544             : GEN
    3545         203 : mplambertW(GEN y)
    3546             : {
    3547         203 :   pari_sp av = avma;
    3548             :   GEN x;
    3549         203 :   long p = 1, s = signe(y);
    3550         203 :   ulong mask = quadratic_prec_mask(realprec(y)-1);
    3551             : 
    3552         203 :   if (s<0) pari_err_DOMAIN("Lw", "y", "<", gen_0, y);
    3553         196 :   if(s==0) return rcopy(y);
    3554         196 :   x = mplambertW0(rtor(y, LOWDEFAULTPREC));
    3555         690 :   while (mask!=1)
    3556             :   {
    3557         298 :     p <<= 1; if (mask & 1) p--;
    3558         298 :     mask >>= 1;
    3559         298 :     x = rtor(x, p+2);
    3560         298 :     x = mulrr(x, divrr(subsr(1, mplog(divrr(x,y))),addrs(x,1)));
    3561             :   }
    3562         196 :   return gerepileuptoleaf(av,x);
    3563             : }
    3564             : 
    3565             : /* exp(t (1 + O(t^n))), n >= 1 */
    3566             : static GEN
    3567          84 : serexp0(long v, long n)
    3568             : {
    3569          84 :   long i, l = n+3;
    3570          84 :   GEN y = cgetg(l, t_SER), t;
    3571          84 :   y[1] = evalsigne(1) | evalvarn(v) | evalvalp(0);
    3572          84 :   gel(y,2) = gen_1; t = gen_1;
    3573        1148 :   for (i = 3; i < l; i++)
    3574             :   {
    3575        1064 :     t = muliu(t, i-2);
    3576        1064 :     gel(y,i) = mkfrac(gen_1, t);
    3577             :   }
    3578          84 :   return y;
    3579             : }
    3580             : 
    3581             : static GEN
    3582          84 : reverse(GEN y)
    3583             : {
    3584          84 :   GEN z = ser2rfrac_i(y);
    3585          84 :   long l = lg(z);
    3586          84 :   return RgX_to_ser(RgXn_reverse(z, l-2), l);
    3587             : }
    3588             : static GEN
    3589         126 : serlambertW(GEN y, long prec)
    3590             : {
    3591             :   GEN x, t, y0;
    3592             :   long n, l, vy, val, v;
    3593             : 
    3594         126 :   if (!signe(y)) return gcopy(y);
    3595         112 :   v = valp(y);
    3596         112 :   vy = varn(y);
    3597         112 :   n = lg(y)-3;
    3598         112 :   y0 = gel(y,2);
    3599         504 :   for (val = 1; val < n; val++)
    3600         476 :     if (!gequal0(polcoeff0(y, val, vy))) break;
    3601         112 :   if (v < 0) pari_err_DOMAIN("lambertw","valuation", "<", gen_0, y);
    3602         105 :   if (val >= n)
    3603             :   {
    3604          21 :     if (v) return zeroser(vy, n);
    3605          21 :     x = glambertW(y0,prec);
    3606          21 :     return scalarser(x, vy, n+1);
    3607             :   }
    3608          84 :   l = 3 + n/val;
    3609          84 :   if (v)
    3610             :   {
    3611          42 :     t = serexp0(vy, l-3);
    3612          42 :     setvalp(t, 1); /* t exp(t) */
    3613          42 :     t = reverse(t);
    3614             :   }
    3615             :   else
    3616             :   {
    3617          42 :     y = serchop0(y);
    3618          42 :     x = glambertW(y0, prec);
    3619             :     /* (x + t) exp(x + t) = (y0 + t y0/x) * exp(t) */
    3620          42 :     t = gmul(deg1pol_shallow(gdiv(y0,x), y0, vy), serexp0(vy, l-3));
    3621          42 :     t = gadd(x, reverse(serchop0(t)));
    3622             :   }
    3623          84 :   t = gsubst(t, vy, y);
    3624          84 :   return normalize(t);
    3625             : }
    3626             : 
    3627             : GEN
    3628         532 : glambertW(GEN y, long prec)
    3629             : {
    3630             :   pari_sp av;
    3631             :   GEN z;
    3632         532 :   switch(typ(y))
    3633             :   {
    3634         203 :     case t_REAL: return mplambertW(y);
    3635           7 :     case t_COMPLEX: pari_err_IMPL("lambert(t_COMPLEX)");
    3636             :     default:
    3637         322 :       av = avma; if (!(z = toser_i(y))) break;
    3638         126 :       return gerepileupto(av, serlambertW(z, prec));
    3639             :   }
    3640         196 :   return trans_eval("lambert",glambertW,y,prec);
    3641             : }
    3642             : 
    3643             : #if 0
    3644             : /* Another Lambert-like function: solution of exp(x)/x=y, y >= e t_REAL,
    3645             :  * using Newton with constant prec. Good for low accuracy */
    3646             : GEN
    3647             : mplambertX(GEN y)
    3648             : {
    3649             :   long bitprec = bit_prec(y)-2;
    3650             :   GEN tmp, x = mplog(y);
    3651             :   if (typ(x) != t_REAL || signe(subrs(x,1))<0)
    3652             :     pari_err(e_MISC,"Lx : argument less than e");
    3653             :   do {
    3654             :     tmp = x;
    3655             :    /* f(x)=x-log(xy), f'(x)=1-1/x */
    3656             :    /* x *= (log(x*y)-1)/(x-1) */
    3657             :     x = mulrr(x, divrr(subrs(mplog(mulrr(x,y)),1), subrs(x,1)));
    3658             :   } while(expo(tmp)-expo(subrr(x,tmp)) < bitprec);
    3659             :   return x;
    3660             : }
    3661             : #endif

Generated by: LCOV version 1.11