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

Generated by: LCOV version 1.11