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 - language - intnum.c (source / functions) Hit Total Coverage
Test: PARI/GP v2.10.0 lcov report (development 21196-f12677d) Lines: 1387 1418 97.8 %
Date: 2017-10-22 06:23:24 Functions: 115 116 99.1 %
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             : #include "pari.h"
      15             : #include "paripriv.h"
      16             : 
      17             : static const long EXTRAPREC =
      18             : #ifdef LONG_IS_64BIT
      19             :   1;
      20             : #else
      21             :   2;
      22             : #endif
      23             : 
      24             : /********************************************************************/
      25             : /**                NUMERICAL INTEGRATION (Romberg)                 **/
      26             : /********************************************************************/
      27             : typedef struct {
      28             :   void *E;
      29             :   GEN (*f)(void *E, GEN);
      30             : } invfun;
      31             : 
      32             : /* 1/x^2 f(1/x) */
      33             : static GEN
      34       12474 : _invf(void *E, GEN x)
      35             : {
      36       12474 :   invfun *S = (invfun*)E;
      37       12474 :   GEN y = ginv(x);
      38       12474 :   return gmul(S->f(S->E, y), gsqr(y));
      39             : }
      40             : 
      41             : /* h and s are arrays of the same length L > D. The h[i] are (decreasing)
      42             :  * step sizes, s[i] is the computed Riemann sum for step size h[i].
      43             :  * Interpolate the last D+1 values so that s ~ polynomial in h of degree D.
      44             :  * Guess that limit_{h->0} = s(0) */
      45             : static GEN
      46         105 : interp(GEN h, GEN s, long L, long bit, long D)
      47             : {
      48         105 :   pari_sp av = avma;
      49             :   long e1,e2;
      50         105 :   GEN dss, ss = polint_i(h + L-D,s + L-D, gen_0, D+1, &dss);
      51             : 
      52         105 :   e1 = gexpo(ss);
      53         105 :   e2 = gexpo(dss);
      54         105 :   if (DEBUGLEVEL>2)
      55             :   {
      56           0 :     err_printf("romb: iteration %ld, guess: %Ps\n", L,ss);
      57           0 :     err_printf("romb: relative error < 2^-%ld [target %ld bits]\n",e1-e2,bit);
      58             :   }
      59         105 :   if (e1-e2 <= bit && (L <= 10 || e1 >= -bit)) { avma = av; return NULL; }
      60          70 :   if (typ(ss) == t_COMPLEX && gequal0(gel(ss,2))) ss = gel(ss,1);
      61          70 :   return ss;
      62             : }
      63             : 
      64             : static GEN
      65           7 : qrom3(void *E, GEN (*eval)(void *, GEN), GEN a, GEN b, long bit)
      66             : {
      67           7 :   const long JMAX = 25, KLOC = 4;
      68             :   GEN ss,s,h,p1,p2,qlint,del,x,sum;
      69           7 :   long j, j1, it, sig, prec = nbits2prec(bit);
      70             : 
      71           7 :   a = gtofp(a,prec);
      72           7 :   b = gtofp(b,prec);
      73           7 :   qlint = subrr(b,a); sig = signe(qlint);
      74           7 :   if (!sig) return gen_0;
      75           7 :   if (sig < 0) { setabssign(qlint); swap(a,b); }
      76             : 
      77           7 :   s = new_chunk(JMAX+KLOC-1);
      78           7 :   h = new_chunk(JMAX+KLOC-1);
      79           7 :   gel(h,0) = real_1(prec);
      80             : 
      81           7 :   p1 = eval(E, a); if (p1 == a) p1 = rcopy(p1);
      82           7 :   p2 = eval(E, b);
      83           7 :   gel(s,0) = gmul2n(gmul(qlint,gadd(p1,p2)),-1);
      84          28 :   for (it=1,j=1; j<JMAX; j++, it<<=1) /* it = 2^(j-1) */
      85             :   {
      86             :     pari_sp av, av2;
      87          28 :     gel(h,j) = real2n(-2*j, prec); /* 2^(-2j) */
      88          28 :     av = avma; del = divru(qlint,it);
      89          28 :     x = addrr(a, shiftr(del,-1));
      90          28 :     av2 = avma;
      91         133 :     for (sum = gen_0, j1 = 1; j1 <= it; j1++, x = addrr(x,del))
      92             :     {
      93         105 :       sum = gadd(sum, eval(E, x));
      94         105 :       if ((j1 & 0x1ff) == 0) gerepileall(av2, 2, &sum,&x);
      95             :     }
      96          28 :     sum = gmul(sum,del);
      97          28 :     gel(s,j) = gerepileupto(av, gmul2n(gadd(gel(s,j-1), sum), -1));
      98          28 :     if (j >= KLOC && (ss = interp(h, s, j, bit-j-6, KLOC)))
      99           7 :       return gmulsg(sig,ss);
     100             :   }
     101           0 :   pari_err_IMPL("intnumromb recovery [too many iterations]");
     102             :   return NULL; /* LCOV_EXCL_LINE */
     103             : }
     104             : 
     105             : static GEN
     106          63 : qrom2(void *E, GEN (*eval)(void *, GEN), GEN a, GEN b, long bit)
     107             : {
     108          63 :   const long JMAX = 16, KLOC = 4;
     109             :   GEN ss,s,h,p1,qlint,del,ddel,x,sum;
     110          63 :   long j, j1, it, sig, prec = nbits2prec(bit);
     111             : 
     112          63 :   a = gtofp(a, prec);
     113          63 :   b = gtofp(b, prec);
     114          63 :   qlint = subrr(b,a); sig = signe(qlint);
     115          63 :   if (!sig)  return gen_0;
     116          63 :   if (sig < 0) { setabssign(qlint); swap(a,b); }
     117             : 
     118          63 :   s = new_chunk(JMAX+KLOC-1);
     119          63 :   h = new_chunk(JMAX+KLOC-1);
     120          63 :   gel(h,0) = real_1(prec);
     121             : 
     122          63 :   p1 = shiftr(addrr(a,b),-1);
     123          63 :   gel(s,0) = gmul(qlint, eval(E, p1));
     124         287 :   for (it=1, j=1; j<JMAX; j++, it*=3) /* it = 3^(j-1) */
     125             :   {
     126             :     pari_sp av, av2;
     127         287 :     gel(h,j) = divru(gel(h,j-1), 9); /* 3^(-2j) */
     128         287 :     av = avma; del = divru(qlint,3*it); ddel = shiftr(del,1);
     129         287 :     x = addrr(a, shiftr(del,-1));
     130         287 :     av2 = avma;
     131        7910 :     for (sum = gen_0, j1 = 1; j1 <= it; j1++)
     132             :     {
     133        7623 :       sum = gadd(sum, eval(E, x)); x = addrr(x,ddel);
     134        7623 :       sum = gadd(sum, eval(E, x)); x = addrr(x,del);
     135        7623 :       if ((j1 & 0x1ff) == 0) gerepileall(av2, 2, &sum,&x);
     136             :     }
     137         287 :     sum = gmul(sum,del); p1 = gdivgs(gel(s,j-1),3);
     138         287 :     gel(s,j) = gerepileupto(av, gadd(p1,sum));
     139         287 :     if (j >= KLOC && (ss = interp(h, s, j, bit-(3*j/2)+3, KLOC)))
     140          63 :       return gmulsg(sig, ss);
     141             :   }
     142           0 :   pari_err_IMPL("intnumromb recovery [too many iterations]");
     143             :   return NULL; /* LCOV_EXCL_LINE */
     144             : }
     145             : 
     146             : /* integrate after change of variables x --> 1/x */
     147             : static GEN
     148          28 : qromi(void *E, GEN (*eval)(void*, GEN), GEN a, GEN b, long bit)
     149             : {
     150          28 :   GEN A = ginv(b), B = ginv(a);
     151             :   invfun S;
     152          28 :   S.f = eval;
     153          28 :   S.E = E; return qrom2(&S, &_invf, A, B, bit);
     154             : }
     155             : 
     156             : /* a < b, assume b "small" (< 100 say) */
     157             : static GEN
     158          28 : rom_bsmall(void *E, GEN (*eval)(void*, GEN), GEN a, GEN b, long bit)
     159             : {
     160          28 :   if (gcmpgs(a,-100) >= 0) return qrom2(E,eval,a,b,bit);
     161           7 :   if (gcmpgs(b, -1) < 0)   return qromi(E,eval,a,b,bit); /* a<-100, b<-1 */
     162             :   /* a<-100, b>=-1, split at -1 */
     163           7 :   return gadd(qromi(E,eval,a,gen_m1,bit),
     164             :               qrom2(E,eval,gen_m1,b,bit));
     165             : }
     166             : 
     167             : static GEN
     168          35 : rombint(void *E, GEN (*eval)(void*, GEN), GEN a, GEN b, long bit)
     169             : {
     170          35 :   long l = gcmp(b,a);
     171             :   GEN z;
     172             : 
     173          35 :   if (!l) return gen_0;
     174          35 :   if (l < 0) swap(a,b);
     175          35 :   if (gcmpgs(b,100) >= 0)
     176             :   {
     177          14 :     if (gcmpgs(a,1) >= 0)
     178           7 :       z = qromi(E,eval,a,b,bit);
     179             :     else /* split at 1 */
     180           7 :       z = gadd(rom_bsmall(E,eval,a,gen_1,bit), qromi(E,eval,gen_1,b,bit));
     181             :   }
     182             :   else
     183          21 :     z = rom_bsmall(E,eval,a,b,bit);
     184          35 :   if (l < 0) z = gneg(z);
     185          35 :   return z;
     186             : }
     187             : 
     188             : GEN
     189          56 : intnumromb_bitprec(void *E, GEN (*f)(void *, GEN), GEN a,GEN b, long fl, long B)
     190             : {
     191          56 :   pari_sp av = avma;
     192             :   GEN z;
     193          56 :   switch(fl)
     194             :   {
     195           7 :     case 0: z = qrom3  (E, f, a, b, B); break;
     196          35 :     case 1: z = rombint(E, f, a, b, B); break;
     197           7 :     case 2: z = qromi  (E, f, a, b, B); break;
     198           7 :     case 3: z = qrom2  (E, f, a, b, B); break;
     199             :     default: pari_err_FLAG("intnumromb"); return NULL; /* LCOV_EXCL_LINE */
     200             :   }
     201          56 :   return gerepileupto(av, z);
     202             : }
     203             : GEN
     204           0 : intnumromb(void *E, GEN (*f)(void *, GEN), GEN a, GEN b, long flag, long prec)
     205           0 : { return intnumromb_bitprec(E,f,a,b,flag,prec2nbits(prec));}
     206             : GEN
     207          56 : intnumromb0_bitprec(GEN a, GEN b, GEN code, long flag, long bit)
     208          56 : { EXPR_WRAP(code, intnumromb_bitprec(EXPR_ARG, a, b, flag, bit)); }
     209             : 
     210             : /********************************************************************/
     211             : /**             NUMERICAL INTEGRATION (Gauss-Legendre)             **/
     212             : /********************************************************************/
     213             : GEN
     214          35 : intnumgaussinit(long n, long prec)
     215             : {
     216          35 :   pari_sp ltop = avma;
     217             :   GEN L, dp1, p1, p2, R, W;
     218          35 :   long prec0 = prec + EXTRAPRECWORD;
     219          35 :   long bitprec = prec2nbits(prec), i, d1;
     220          35 :   if (n <= 0) n = (long)(bitprec*0.2258);
     221          35 :   if (odd(n)) n++;
     222          35 :   if (n == 2) n = 4;
     223             :   /* n even >= 4, p1 is even */
     224          35 :   prec = nbits2prec(3*bitprec/2 + 32);
     225          35 :   L = pollegendre(n, 0); /* L_n = p1(x^2) */
     226          35 :   p1 = Q_remove_denom(RgX_deflate(L, 2), &dp1);
     227          35 :   d1 = vali(dp1);
     228          35 :   p2 = ZX_deriv(p1); /* L_n' = 2x p2(x^2) / 2^d1 */
     229          35 :   R = ZX_Uspensky(p1, gen_0, 1, 3*bitprec/2 + 32); /* positive roots of p1 */
     230          35 :   n >>= 1;
     231          35 :   W = cgetg(n+1, t_VEC);
     232         973 :   for (i = 1; i <= n; ++i)
     233             :   {
     234         938 :     GEN t, r2 = gel(R,i);
     235         938 :     if (typ(r2) != t_REAL) r2 = gtofp(r2, prec);
     236         938 :     gel(R,i) = sqrtr_abs(r2); /* positive root of L_n */
     237             :     /* 2 / (L'(r)^2(1-r^2)) =  2^(2d1 - 1) / (1-r2)r2 (p2(r2))^2 */
     238         938 :     t = mulrr(subrr(r2, sqrr(r2)), sqrr(poleval(p2, r2)));
     239         938 :     shiftr_inplace(t,1-2*d1);
     240         938 :     gel(W,i) = invr(t);
     241             :   }
     242          35 :   R = gprec_wtrunc(R,prec0);
     243          35 :   W = gprec_wtrunc(W,prec0);
     244          35 :   return gerepilecopy(ltop, mkvec2(R,W));
     245             : }
     246             : 
     247             : GEN
     248          56 : intnumgauss(void *E, GEN (*eval)(void*, GEN), GEN a, GEN b, GEN tab, long prec)
     249             : {
     250          56 :   pari_sp ltop = avma;
     251             :   GEN R, W, bma, bpa, S;
     252             :   long n, i;
     253          56 :   if (!tab)
     254           7 :     tab = intnumgaussinit(0,prec);
     255          49 :   else if (typ(tab) != t_INT)
     256             :   {
     257          35 :     if (typ(tab) != t_VEC || lg(tab) != 3)
     258           7 :       pari_err_TYPE("intnumgauss",tab);
     259             :   }
     260             :   else
     261          14 :     tab = intnumgaussinit(itos(tab),prec);
     262             : 
     263          49 :   R = gel(tab,1); n = lg(R)-1;
     264          49 :   W = gel(tab,2);
     265          49 :   a = gprec_w(a, prec+EXTRAPREC);
     266          49 :   b = gprec_w(b, prec+EXTRAPREC);
     267          49 :   bma = gmul2n(gsub(b,a), -1); /* (b-a)/2 */
     268          49 :   bpa = gadd(bma, a); /* (b+a)/2 */
     269          49 :   S = gen_0;
     270        1491 :   for (i = 1; i <= n; ++i)
     271             :   {
     272        1442 :     GEN r = gel(R,i);
     273        1442 :     GEN P = eval(E, gadd(bpa, gmul(bma, r)));
     274        1442 :     GEN M = eval(E, gsub(bpa, gmul(bma, r)));
     275        1442 :     S = gadd(S, gmul(gel(W,i), gadd(P,M)));
     276             :   }
     277          49 :   return gerepilecopy(ltop, gprec_wtrunc(gmul(bma,S), prec));
     278             : }
     279             : 
     280             : GEN
     281          56 : intnumgauss0(GEN a, GEN b, GEN code, GEN tab, long prec)
     282          56 : { EXPR_WRAP(code, intnumgauss(EXPR_ARG, a, b, tab, prec)); }
     283             : 
     284             : /********************************************************************/
     285             : /**                DOUBLE EXPONENTIAL INTEGRATION                  **/
     286             : /********************************************************************/
     287             : 
     288             : typedef struct _intdata {
     289             :   long eps;  /* bit accuracy of current precision */
     290             :   long l; /* table lengths */
     291             :   GEN tabx0; /* abscissa phi(0) for t = 0 */
     292             :   GEN tabw0; /* weight phi'(0) for t = 0 */
     293             :   GEN tabxp; /* table of abscissas phi(kh) for k > 0 */
     294             :   GEN tabwp; /* table of weights phi'(kh) for k > 0 */
     295             :   GEN tabxm; /* table of abscissas phi(kh) for k < 0, possibly empty */
     296             :   GEN tabwm; /* table of weights phi'(kh) for k < 0, possibly empty */
     297             :   GEN h; /* integration step */
     298             : } intdata;
     299             : 
     300             : static const long LGTAB = 8;
     301             : #define TABh(v) gel(v,1)
     302             : #define TABx0(v) gel(v,2)
     303             : #define TABw0(v) gel(v,3)
     304             : #define TABxp(v) gel(v,4)
     305             : #define TABwp(v) gel(v,5)
     306             : #define TABxm(v) gel(v,6)
     307             : #define TABwm(v) gel(v,7)
     308             : 
     309             : static int
     310       11137 : isinR(GEN z) { return is_real_t(typ(z)); }
     311             : static int
     312       10787 : isinC(GEN z)
     313       10787 : { return (typ(z) == t_COMPLEX)? isinR(gel(z,1)) && isinR(gel(z,2)): isinR(z); }
     314             : 
     315             : static int
     316        8820 : checktabsimp(GEN tab)
     317             : {
     318             :   long L, LN, LW;
     319        8820 :   if (!tab || typ(tab) != t_VEC) return 0;
     320        8820 :   if (lg(tab) != LGTAB) return 0;
     321        8820 :   if (typ(TABxp(tab)) != t_VEC) return 0;
     322        8820 :   if (typ(TABwp(tab)) != t_VEC) return 0;
     323        8820 :   if (typ(TABxm(tab)) != t_VEC) return 0;
     324        8820 :   if (typ(TABwm(tab)) != t_VEC) return 0;
     325        8820 :   L = lg(TABxp(tab)); if (lg(TABwp(tab)) != L) return 0;
     326        8820 :   LN = lg(TABxm(tab)); if (LN != 1 && LN != L) return 0;
     327        8820 :   LW = lg(TABwm(tab)); if (LW != 1 && LW != L) return 0;
     328        8820 :   return 1;
     329             : }
     330             : 
     331             : static int
     332         476 : checktabdoub(GEN tab)
     333             : {
     334             :   long L;
     335         476 :   if (typ(tab) != t_VEC) return 0;
     336         476 :   if (lg(tab) != LGTAB) return 0;
     337         476 :   L = lg(TABxp(tab));
     338         476 :   if (lg(TABwp(tab)) != L) return 0;
     339         476 :   if (lg(TABxm(tab)) != L) return 0;
     340         476 :   if (lg(TABwm(tab)) != L) return 0;
     341         476 :   return 1;
     342             : }
     343             : 
     344             : static int
     345        4445 : checktab(GEN tab)
     346             : {
     347        4445 :   if (typ(tab) != t_VEC) return 0;
     348        4445 :   if (lg(tab) != 3) return checktabsimp(tab);
     349          14 :   return checktabsimp(gel(tab,1))
     350           7 :       && checktabsimp(gel(tab,2));
     351             : }
     352             : 
     353             : /* the TUNE parameter is heuristic */
     354             : static void
     355         693 : intinit_start(intdata *D, long m, double TUNE, long prec)
     356             : {
     357         693 :   long l, n, bitprec = prec2nbits(prec);
     358         693 :   double d = bitprec*LOG10_2;
     359         693 :   GEN h, nh, pi = mppi(prec);
     360             : 
     361         693 :   n = (long)ceil(d*log(d) / TUNE); /* heuristic */
     362             :   /* nh ~ log(2npi/log(n)) */
     363         693 :   nh = logr_abs(divrr(mulur(2*n, pi), logr_abs(utor(n,prec))));
     364         693 :   h = divru(nh, n);
     365         693 :   if (m > 0) { h = gmul2n(h,-m); n <<= m; }
     366         693 :   D->h = h;
     367         693 :   D->eps = bitprec;
     368         693 :   D->l = l = n+1;
     369         693 :   D->tabxp = cgetg(l, t_VEC);
     370         693 :   D->tabwp = cgetg(l, t_VEC);
     371         693 :   D->tabxm = cgetg(l, t_VEC);
     372         693 :   D->tabwm = cgetg(l, t_VEC);
     373         693 : }
     374             : 
     375             : static GEN
     376         693 : intinit_end(intdata *D, long pnt, long mnt)
     377             : {
     378         693 :   GEN v = cgetg(LGTAB, t_VEC);
     379         693 :   if (pnt < 0) pari_err_DOMAIN("intnuminit","table length","<",gen_0,stoi(pnt));
     380         693 :   TABx0(v) = D->tabx0;
     381         693 :   TABw0(v) = D->tabw0;
     382         693 :   TABh(v) = D->h;
     383         693 :   TABxp(v) = D->tabxp; setlg(D->tabxp, pnt+1);
     384         693 :   TABwp(v) = D->tabwp; setlg(D->tabwp, pnt+1);
     385         693 :   TABxm(v) = D->tabxm; setlg(D->tabxm, mnt+1);
     386         693 :   TABwm(v) = D->tabwm; setlg(D->tabwm, mnt+1); return v;
     387             : }
     388             : 
     389             : /* divide by 2 in place */
     390             : static GEN
     391      261378 : divr2_ip(GEN x) { shiftr_inplace(x, -1); return x; }
     392             : 
     393             : /* phi(t)=tanh((Pi/2)sinh(t)): from -1 to 1, hence also from a to b compact
     394             :  * interval */
     395             : static GEN
     396         119 : inittanhsinh(long m, long prec)
     397             : {
     398         119 :   GEN et, ex, pi = mppi(prec);
     399         119 :   long k, nt = -1;
     400             :   intdata D;
     401             : 
     402         119 :   intinit_start(&D, m, 1.86, prec);
     403         119 :   D.tabx0 = real_0(prec);
     404         119 :   D.tabw0 = Pi2n(-1,prec);
     405         119 :   et = ex = mpexp(D.h);
     406       27602 :   for (k = 1; k < D.l; k++)
     407             :   {
     408             :     GEN xp, wp, ct, st, z;
     409             :     pari_sp av;
     410       27602 :     gel(D.tabxp,k) = cgetr(prec);
     411       27602 :     gel(D.tabwp,k) = cgetr(prec); av = avma;
     412       27602 :     ct = divr2_ip(addrr(et, invr(et))); /* ch(kh) */
     413       27602 :     st = subrr(et, ct); /* sh(kh) */
     414       27602 :     z = invr( addrs(mpexp(mulrr(pi, st)), 1) );
     415       27602 :     shiftr_inplace(z, 1);
     416       27602 :     xp = subsr(1, z);
     417       27602 :     wp = divr2_ip(mulrr(mulrr(pi,ct), mulrr(z, subsr(2, z))));
     418       27602 :     if (expo(wp) < -D.eps) { nt = k-1; break; }
     419       27601 :     affrr(xp, gel(D.tabxp,k));
     420       27601 :     if (absrnz_equal1(gel(D.tabxp,k))) { nt = k-1; break; }
     421       27483 :     affrr(wp, gel(D.tabwp,k)); et = gerepileuptoleaf(av, mulrr(et, ex));
     422             :   }
     423         119 :   return intinit_end(&D, nt, 0);
     424             : }
     425             : 
     426             : /* phi(t)=sinh(sinh(t)): from -oo to oo, slowly decreasing, at least
     427             :  * as 1/x^2. */
     428             : static GEN
     429          14 : initsinhsinh(long m, long prec)
     430             : {
     431             :   pari_sp av;
     432             :   GEN et, ct, st, ex;
     433          14 :   long k, nt = -1;
     434             :   intdata D;
     435             : 
     436          14 :   intinit_start(&D, m, 0.666, prec);
     437          14 :   D.tabx0 = real_0(prec);
     438          14 :   D.tabw0 = real_1(prec);
     439          14 :   et = ex = mpexp(D.h);
     440        8184 :   for (k = 1; k < D.l; k++)
     441             :   {
     442             :     GEN xp, wp, ext, exu;
     443        8184 :     gel(D.tabxp,k) = cgetr(prec);
     444        8184 :     gel(D.tabwp,k) = cgetr(prec); av = avma;
     445        8184 :     ct = divr2_ip(addrr(et, invr(et)));
     446        8184 :     st = subrr(et, ct);
     447        8184 :     ext = mpexp(st);
     448        8184 :     exu = invr(ext);
     449        8184 :     xp = divr2_ip(subrr(ext, exu));
     450        8184 :     wp = divr2_ip(mulrr(ct, addrr(ext, exu)));
     451        8184 :     if (expo(wp) - 2*expo(xp) < -D.eps) { nt = k-1; break; }
     452        8170 :     affrr(xp, gel(D.tabxp,k));
     453        8170 :     affrr(wp, gel(D.tabwp,k)); et = gerepileuptoleaf(av, mulrr(et, ex));
     454             :   }
     455          14 :   return intinit_end(&D, nt, 0);
     456             : }
     457             : 
     458             : /* phi(t)=2sinh(t): from -oo to oo, exponentially decreasing as exp(-x) */
     459             : static GEN
     460         126 : initsinh(long m, long prec)
     461             : {
     462             :   pari_sp av;
     463             :   GEN et, ex, eti, xp, wp;
     464         126 :   long k, nt = -1;
     465             :   intdata D;
     466             : 
     467         126 :   intinit_start(&D, m, 1.0, prec);
     468         126 :   D.tabx0 = real_0(prec);
     469         126 :   D.tabw0 = real2n(1, prec);
     470         126 :   et = ex = mpexp(D.h);
     471       38136 :   for (k = 1; k < D.l; k++)
     472             :   {
     473       38136 :     gel(D.tabxp,k) = cgetr(prec);
     474       38136 :     gel(D.tabwp,k) = cgetr(prec); av = avma;
     475       38136 :     eti = invr(et);
     476       38136 :     xp = subrr(et, eti);
     477       38136 :     wp = addrr(et, eti);
     478       38136 :     if (cmprs(xp, (long)(LOG2*(expo(wp)+D.eps) + 1)) > 0) { nt = k-1; break; }
     479       38010 :     affrr(xp, gel(D.tabxp,k));
     480       38010 :     affrr(wp, gel(D.tabwp,k)); et = gerepileuptoleaf(av, mulrr(et, ex));
     481             :   }
     482         126 :   return intinit_end(&D, nt, 0);
     483             : }
     484             : 
     485             : /* phi(t)=exp(2sinh(t)): from 0 to oo, slowly decreasing at least as 1/x^2 */
     486             : static GEN
     487         238 : initexpsinh(long m, long prec)
     488             : {
     489             :   GEN et, ex;
     490         238 :   long k, nt = -1;
     491             :   intdata D;
     492             : 
     493         238 :   intinit_start(&D, m, 1.05, prec);
     494         238 :   D.tabx0 = real_1(prec);
     495         238 :   D.tabw0 = real2n(1, prec);
     496         238 :   ex = mpexp(D.h);
     497         238 :   et = real_1(prec);
     498      111791 :   for (k = 1; k < D.l; k++)
     499             :   {
     500             :     GEN t, eti, xp;
     501      111791 :     et = mulrr(et, ex);
     502      111791 :     eti = invr(et); t = addrr(et, eti);
     503      111791 :     xp = mpexp(subrr(et, eti));
     504      111791 :     gel(D.tabxp,k) = xp;
     505      111791 :     gel(D.tabwp,k) = mulrr(xp, t);
     506      111791 :     gel(D.tabxm,k) = invr(xp);
     507      111791 :     gel(D.tabwm,k) = mulrr(gel(D.tabxm,k), t);
     508      111791 :     if (expo(gel(D.tabxm,k)) < -D.eps) { nt = k-1; break; }
     509             :   }
     510         238 :   return intinit_end(&D, nt, nt);
     511             : }
     512             : 
     513             : /* phi(t)=exp(t-exp(-t)) : from 0 to \infty, exponentially decreasing. */
     514             : static GEN
     515          84 : initexpexp(long m, long prec)
     516             : {
     517             :   pari_sp av;
     518             :   GEN et, ex;
     519          84 :   long k, nt = -1;
     520             :   intdata D;
     521             : 
     522          84 :   intinit_start(&D, m, 1.76, prec);
     523          84 :   D.tabx0 = mpexp(real_m1(prec));
     524          84 :   D.tabw0 = gmul2n(D.tabx0, 1);
     525          84 :   et = ex = mpexp(negr(D.h));
     526       35644 :   for (k = 1; k < D.l; k++)
     527             :   {
     528             :     GEN xp, xm, wp, wm, eti, kh;
     529       35644 :     gel(D.tabxp,k) = cgetr(prec);
     530       35644 :     gel(D.tabwp,k) = cgetr(prec);
     531       35644 :     gel(D.tabxm,k) = cgetr(prec);
     532       35644 :     gel(D.tabwm,k) = cgetr(prec); av = avma;
     533       35644 :     eti = invr(et); kh = mulur(k,D.h);
     534       35644 :     xp = mpexp(subrr(kh, et));
     535       35644 :     xm = mpexp(negr(addrr(kh, eti)));
     536       35644 :     wp = mulrr(xp, addsr(1, et));
     537       35644 :     wm = mulrr(xm, addsr(1, eti));
     538       35644 :     if (expo(xm) < -D.eps && cmprs(xp, (long)(LOG2*(expo(wp)+D.eps) + 1)) > 0) { nt = k-1; break; }
     539       35560 :     affrr(xp, gel(D.tabxp,k));
     540       35560 :     affrr(wp, gel(D.tabwp,k));
     541       35560 :     affrr(xm, gel(D.tabxm,k));
     542       35560 :     affrr(wm, gel(D.tabwm,k)); et = gerepileuptoleaf(av, mulrr(et, ex));
     543             :   }
     544          84 :   return intinit_end(&D, nt, nt);
     545             : }
     546             : 
     547             : /* phi(t)=(Pi/h)*t/(1-exp(-sinh(t))) from 0 to oo, sine oscillation */
     548             : static GEN
     549         112 : initnumsine(long m, long prec)
     550             : {
     551             :   pari_sp av;
     552         112 :   GEN invh, et, eti, ex, pi = mppi(prec);
     553         112 :   long exh, k, nt = -1;
     554             :   intdata D;
     555             : 
     556         112 :   intinit_start(&D, m, 0.666, prec);
     557         112 :   invh = invr(D.h);
     558         112 :   D.tabx0 = mulrr(pi, invh);
     559         112 :   D.tabw0 = gmul2n(D.tabx0,-1);
     560         112 :   exh = expo(invh); /*  expo(1/h) */
     561         112 :   et = ex = mpexp(D.h);
     562       90811 :   for (k = 1; k < D.l; k++)
     563             :   {
     564             :     GEN xp,xm, wp,wm, ct,st, extp,extp1,extp2, extm,extm1,extm2, kct, kpi;
     565       90811 :     gel(D.tabxp,k) = cgetr(prec);
     566       90811 :     gel(D.tabwp,k) = cgetr(prec);
     567       90811 :     gel(D.tabxm,k) = cgetr(prec);
     568       90811 :     gel(D.tabwm,k) = cgetr(prec); av = avma;
     569       90811 :     eti = invr(et); /* exp(-kh) */
     570       90811 :     ct = divr2_ip(addrr(et, eti)); /* ch(kh) */
     571       90811 :     st = divr2_ip(subrr(et, eti)); /* sh(kh) */
     572       90811 :     extp = mpexp(st);  extp1 = subsr(1, extp);
     573       90811 :     extp2 = invr(extp1); /* 1/(1-exp(sh(kh))) */
     574       90811 :     extm = invr(extp); extm1 = subsr(1, extm);
     575       90811 :     extm2 = invr(extm1);/* 1/(1-exp(sh(-kh))) */
     576       90811 :     kpi = mulur(k, pi);
     577       90811 :     kct = mulur(k, ct);
     578       90811 :     extm1 = mulrr(extm1, invh);
     579       90811 :     extp1 = mulrr(extp1, invh);
     580       90811 :     xp = mulrr(kpi, extm2); /* phi(kh) */
     581       90811 :     wp = mulrr(subrr(extm1, mulrr(kct, extm)), mulrr(pi, sqrr(extm2)));
     582       90811 :     xm = mulrr(negr(kpi), extp2); /* phi(-kh) */
     583       90811 :     wm = mulrr(addrr(extp1, mulrr(kct, extp)), mulrr(pi, sqrr(extp2)));
     584       90811 :     if (expo(wm) < -D.eps && expo(extm) + exh + expu(10 * k) < -D.eps) { nt = k-1; break; }
     585       90699 :     affrr(xp, gel(D.tabxp,k));
     586       90699 :     affrr(wp, gel(D.tabwp,k));
     587       90699 :     affrr(xm, gel(D.tabxm,k));
     588       90699 :     affrr(wm, gel(D.tabwm,k)); et = gerepileuptoleaf(av, mulrr(et, ex));
     589             :   }
     590         112 :   return intinit_end(&D, nt, nt);
     591             : }
     592             : 
     593             : /* End of initialization functions. These functions can be executed once
     594             :  * and for all for a given accuracy and type of integral ([a,b], [a,oo[ or
     595             :  * ]-oo,a], ]-oo,oo[) */
     596             : 
     597             : /* The numbers below can be changed, but NOT the ordering */
     598             : enum {
     599             :   f_REG    = 0, /* regular function */
     600             :   f_SING   = 1, /* algebraic singularity */
     601             :   f_YSLOW  = 2, /* +\infty, slowly decreasing, at least x^(-2)  */
     602             :   f_YVSLO  = 3, /* +\infty, very slowly decreasing, worse than x^(-2) */
     603             :   f_YFAST  = 4, /* +\infty, exponentially decreasing */
     604             :   f_YOSCS  = 5, /* +\infty, sine oscillating */
     605             :   f_YOSCC  = 6  /* +\infty, cosine oscillating */
     606             : };
     607             : /* is finite ? */
     608             : static int
     609         966 : is_fin_f(long c) { return c == f_REG || c == f_SING; }
     610             : /* is oscillatory ? */
     611             : static int
     612         140 : is_osc(long c) { long a = labs(c); return a == f_YOSCC|| a == f_YOSCS; }
     613             : 
     614             : /* All inner functions such as intn, etc... must be called with a
     615             :  * valid 'tab' table. The wrapper intnum provides a higher level interface */
     616             : 
     617             : /* compute \int_a^b f(t)dt with [a,b] compact and f nonsingular. */
     618             : static GEN
     619        3745 : intn(void *E, GEN (*eval)(void*, GEN), GEN a, GEN b, GEN tab)
     620             : {
     621             :   GEN tabx0, tabw0, tabxp, tabwp;
     622             :   GEN bpa, bma, bmb, S;
     623             :   long i;
     624        3745 :   pari_sp ltop = avma, av;
     625             : 
     626        3745 :   if (!checktabsimp(tab)) pari_err_TYPE("intnum",tab);
     627        3745 :   tabx0 = TABx0(tab); tabw0 = TABw0(tab);
     628        3745 :   tabxp = TABxp(tab); tabwp = TABwp(tab);
     629        3745 :   bpa = gmul2n(gadd(b, a), -1); /* (b+a)/2 */
     630        3745 :   bma = gsub(bpa, a); /* (b-a)/2 */
     631        3745 :   av = avma;
     632        3745 :   bmb = gmul(bma, tabx0); /* (b-a)/2 phi(0) */
     633             :   /* phi'(0) f( (b+a)/2 + (b-a)/2 * phi(0) ) */
     634        3745 :   S = gmul(tabw0, eval(E, gadd(bpa, bmb)));
     635      982639 :   for (i = lg(tabxp)-1; i > 0; i--)
     636             :   {
     637             :     GEN SP, SM;
     638      978894 :     bmb = gmul(bma, gel(tabxp,i));
     639      978894 :     SP = eval(E, gsub(bpa, bmb));
     640      978894 :     SM = eval(E, gadd(bpa, bmb));
     641      978894 :     S = gadd(S, gmul(gel(tabwp,i), gadd(SP, SM)));
     642      978894 :     if ((i & 0x7f) == 1) S = gerepileupto(av, S);
     643             :   }
     644        3745 :   return gerepileupto(ltop, gmul(S, gmul(bma, TABh(tab))));
     645             : }
     646             : 
     647             : /* compute \int_a^b f(t)dt with [a,b] compact, possible singularity with
     648             :  * exponent a[2] at lower extremity, b regular. Use tanh(sinh(t)). */
     649             : static GEN
     650          42 : intnsing(void *E, GEN (*eval)(void*, GEN), GEN a, GEN b, GEN tab, long prec)
     651             : {
     652             :   GEN tabx0, tabw0, tabxp, tabwp, ea, ba, S;
     653             :   long i;
     654          42 :   pari_sp ltop = avma, av;
     655             : 
     656          42 :   if (!checktabsimp(tab)) pari_err_TYPE("intnum",tab);
     657          42 :   tabx0 = TABx0(tab); tabw0 = TABw0(tab);
     658          42 :   tabxp = TABxp(tab); tabwp = TABwp(tab);
     659          42 :   ea = ginv(gaddsg(1, gel(a,2)));
     660          42 :   a = gel(a,1);
     661          42 :   ba = gdiv(gsub(b, a), gpow(gen_2, ea, prec));
     662          42 :   av = avma;
     663          42 :   S = gmul(gmul(tabw0, ba), eval(E, gadd(gmul(ba, addsr(1, tabx0)), a)));
     664       11046 :   for (i = lg(tabxp)-1; i > 0; i--)
     665             :   {
     666       11004 :     GEN p = addsr(1, gel(tabxp,i));
     667       11004 :     GEN m = subsr(1, gel(tabxp,i));
     668       11004 :     GEN bp = gmul(ba, gpow(p, ea, prec));
     669       11004 :     GEN bm = gmul(ba, gpow(m, ea, prec));
     670       11004 :     GEN SP = gmul(gdiv(bp, p), eval(E, gadd(bp, a)));
     671       11004 :     GEN SM = gmul(gdiv(bm, m), eval(E, gadd(bm, a)));
     672       11004 :     S = gadd(S, gmul(gel(tabwp,i), gadd(SP, SM)));
     673       11004 :     if ((i & 0x7f) == 1) S = gerepileupto(av, S);
     674             :   }
     675          42 :   return gerepileupto(ltop, gmul(gmul(S, TABh(tab)), ea));
     676             : }
     677             : 
     678      170450 : static GEN id(GEN x) { return x; }
     679             : 
     680             : /* compute  \int_a^oo f(t)dt if si>0 or \int_{-oo}^a f(t)dt if si<0$.
     681             :  * Use exp(2sinh(t)) for slowly decreasing functions, exp(1+t-exp(-t)) for
     682             :  * exponentially decreasing functions, and (pi/h)t/(1-exp(-sinh(t))) for
     683             :  * oscillating functions. */
     684             : static GEN
     685         476 : intninfpm(void *E, GEN (*eval)(void*, GEN), GEN a, long sb, GEN tab)
     686             : {
     687             :   GEN tabx0, tabw0, tabxp, tabwp, tabxm, tabwm;
     688             :   GEN S;
     689             :   long L, i;
     690         476 :   pari_sp av = avma;
     691             : 
     692         476 :   if (!checktabdoub(tab)) pari_err_TYPE("intnum",tab);
     693         476 :   tabx0 = TABx0(tab); tabw0 = TABw0(tab);
     694         476 :   tabxp = TABxp(tab); tabwp = TABwp(tab); L = lg(tabxp);
     695         476 :   tabxm = TABxm(tab); tabwm = TABwm(tab);
     696         476 :   if (gequal0(a))
     697             :   {
     698         238 :     GEN (*NEG)(GEN) = sb > 0? id: gneg;
     699         238 :     S = gmul(tabw0, eval(E, NEG(tabx0)));
     700      128814 :     for (i = 1; i < L; i++)
     701             :     {
     702      128576 :       GEN SP = eval(E, NEG(gel(tabxp,i)));
     703      128576 :       GEN SM = eval(E, NEG(gel(tabxm,i)));
     704      128576 :       S = gadd(S, gadd(gmul(gel(tabwp,i), SP), gmul(gel(tabwm,i), SM)));
     705      128576 :       if ((i & 0x7f) == 1) S = gerepileupto(av, S);
     706             :     }
     707             :   }
     708         238 :   else if (gexpo(a) <= 0 || is_osc(sb))
     709         105 :   { /* a small */
     710         105 :     GEN (*ADD)(GEN,GEN) = sb > 0? gadd: gsub;
     711         105 :     S = gmul(tabw0, eval(E, ADD(a, tabx0)));
     712       61460 :     for (i = 1; i < L; i++)
     713             :     {
     714       61355 :       GEN SP = eval(E, ADD(a, gel(tabxp,i)));
     715       61355 :       GEN SM = eval(E, ADD(a, gel(tabxm,i)));
     716       61355 :       S = gadd(S, gadd(gmul(gel(tabwp,i), SP), gmul(gel(tabwm,i), SM)));
     717       61355 :       if ((i & 0x7f) == 1) S = gerepileupto(av, S);
     718             :     }
     719             :   }
     720             :   else
     721             :   { /* a large, |a|*\int_sgn(a)^{oo} f(|a|*x)dx (sb > 0)*/
     722         133 :     GEN (*ADD)(long,GEN) = sb > 0? addsr: subsr;
     723         133 :     long sa = gsigne(a);
     724         133 :     GEN A = sa > 0? a: gneg(a);
     725         133 :     pari_sp av2 = avma;
     726         133 :     S = gmul(tabw0, eval(E, gmul(A, ADD(sa, tabx0))));
     727       88705 :     for (i = 1; i < L; i++)
     728             :     {
     729       88572 :       GEN SP = eval(E, gmul(A, ADD(sa, gel(tabxp,i))));
     730       88572 :       GEN SM = eval(E, gmul(A, ADD(sa, gel(tabxm,i))));
     731       88572 :       S = gadd(S, gadd(gmul(gel(tabwp,i), SP), gmul(gel(tabwm,i), SM)));
     732       88572 :       if ((i & 0x7f) == 1) S = gerepileupto(av2, S);
     733             :     }
     734         133 :     S = gmul(S,A);
     735             :   }
     736         476 :   return gerepileupto(av, gmul(S, TABh(tab)));
     737             : }
     738             : 
     739             : /* Compute  \int_{-oo}^oo f(t)dt
     740             :  * use sinh(sinh(t)) for slowly decreasing functions and sinh(t) for
     741             :  * exponentially decreasing functions.
     742             :  * HACK: in case TABwm(tab) contains something, assume function to be integrated
     743             :  * satisfies f(-x) = conj(f(x)).
     744             :  */
     745             : static GEN
     746         581 : intninfinf(void *E, GEN (*eval)(void*, GEN), GEN tab)
     747             : {
     748             :   GEN tabx0, tabw0, tabxp, tabwp, tabwm;
     749             :   GEN S;
     750             :   long L, i, spf;
     751         581 :   pari_sp ltop = avma;
     752             : 
     753         581 :   if (!checktabsimp(tab)) pari_err_TYPE("intnum",tab);
     754         581 :   tabx0 = TABx0(tab); tabw0 = TABw0(tab);
     755         581 :   tabxp = TABxp(tab); tabwp = TABwp(tab); L = lg(tabxp);
     756         581 :   tabwm = TABwm(tab);
     757         581 :   spf = (lg(tabwm) == lg(tabwp));
     758         581 :   S = gmul(tabw0, eval(E, tabx0));
     759         581 :   if (spf) S = gmul2n(real_i(S), -1);
     760      176099 :   for (i = L-1; i > 0; i--)
     761             :   {
     762      175518 :     GEN SP = eval(E, gel(tabxp,i));
     763      175518 :     if (spf)
     764      170044 :       S = gadd(S, real_i(gmul(gel(tabwp,i), SP)));
     765             :     else
     766             :     {
     767        5474 :       GEN SM = eval(E, negr(gel(tabxp,i)));
     768        5474 :       S = gadd(S, gmul(gel(tabwp,i), gadd(SP,SM)));
     769             :     }
     770      175518 :     if ((i & 0x7f) == 1) S = gerepileupto(ltop, S);
     771             :   }
     772         581 :   if (spf) S = gmul2n(S,1);
     773         581 :   return gerepileupto(ltop, gmul(S, TABh(tab)));
     774             : }
     775             : 
     776             : /* general num integration routine int_a^b f(t)dt, where a and b are as follows:
     777             :  - a scalar : the scalar, no singularity worse than logarithmic at a.
     778             :  - [a, e] : the scalar a, singularity exponent -1 < e <= 0.
     779             :  - +oo: slowly decreasing function (at least O(t^-2))
     780             :  - [[+oo], a], a nonnegative real : +oo, function behaving like exp(-a|t|)
     781             :  - [[+oo], e], e < -1 : +oo, function behaving like t^e
     782             :  - [[+oo], a*I], a > 0 real : +oo, function behaving like cos(at)
     783             :  - [[+oo], a*I], a < 0 real : +oo, function behaving like sin(at)
     784             :  and similarly at -oo */
     785             : static GEN
     786        1953 : f_getycplx(GEN a, long prec)
     787             : {
     788             :   long s;
     789             :   GEN tmp, a2R, a2I;
     790             : 
     791        1953 :   if (lg(a) == 2 || gequal0(gel(a,2))) return gen_1;
     792        1911 :   a2R = real_i(gel(a,2));
     793        1911 :   a2I = imag_i(gel(a,2));
     794        1911 :   s = gsigne(a2I); if (s < 0) a2I = gneg(a2I);
     795        1911 :   tmp = s ? ginv(a2I) : ginv(a2R);
     796        1911 :   if (gprecision(tmp) < prec) tmp = gprec_w(tmp, prec);
     797        1911 :   return tmp;
     798             : }
     799             : 
     800             : static void
     801          21 : err_code(GEN a, const char *name)
     802             : {
     803          21 :   char *s = stack_sprintf("intnum [incorrect %s]", name);
     804          21 :   pari_err_TYPE(s, a);
     805           0 : }
     806             : 
     807             : /* a = [[+/-oo], alpha]*/
     808             : static long
     809        2576 : code_aux(GEN a, const char *name)
     810             : {
     811        2576 :   GEN re, im, alpha = gel(a,2);
     812             :   long s;
     813        2576 :   if (!isinC(alpha)) err_code(a, name);
     814        2569 :   re = real_i(alpha);
     815        2569 :   im = imag_i(alpha);
     816        2569 :   s = gsigne(im);
     817        2569 :   if (s)
     818             :   {
     819         266 :     if(!gequal0(re))
     820           7 :       pari_warn(warner,"real(z)*imag(z)!=0 in endpoint code, real(z) ignored");
     821         266 :     return s > 0 ? f_YOSCC : f_YOSCS;
     822             :   }
     823        2303 :   if (gequal0(re) || gcmpgs(re, -2)<=0) return f_YSLOW;
     824        2065 :   if (gsigne(re) > 0) return f_YFAST;
     825         231 :   if (gcmpgs(re, -1) >= 0) err_code(a, name);
     826         231 :   return f_YVSLO;
     827             : }
     828             : 
     829             : static long
     830       10948 : transcode(GEN a, const char *name)
     831             : {
     832             :   GEN a1, a2;
     833       10948 :   switch(typ(a))
     834             :   {
     835        2681 :     case t_VEC: break;
     836             :     case t_INFINITY:
     837         140 :       return inf_get_sign(a) == 1 ? f_YSLOW: -f_YSLOW;
     838        8127 :     default: if (!isinC(a)) err_code(a,name);
     839        8120 :       return f_REG;
     840             :   }
     841        2681 :   switch(lg(a))
     842             :   {
     843          14 :     case 2: return gsigne(gel(a,1)) > 0 ? f_YSLOW : -f_YSLOW;
     844        2660 :     case 3: break;
     845           7 :     default: err_code(a,name);
     846             :   }
     847        2660 :   a1 = gel(a,1);
     848        2660 :   a2 = gel(a,2);
     849        2660 :   switch(typ(a1))
     850             :   {
     851             :     case t_VEC:
     852          14 :       if (lg(a1) != 2) err_code(a,name);
     853          14 :       return gsigne(gel(a1,1)) * code_aux(a, name);
     854             :     case t_INFINITY:
     855        2562 :       return inf_get_sign(a1) * code_aux(a, name);
     856             :     default:
     857          84 :       if (!isinC(a1) || !isinR(a2) || gcmpgs(a2, -1) <= 0) err_code(a,name);
     858          84 :       return gsigne(a2) < 0 ? f_SING : f_REG;
     859             :   }
     860             : }
     861             : 
     862             : /* computes the necessary tabs, knowing a, b and m */
     863             : static GEN
     864         357 : homtab(GEN tab, GEN k)
     865             : {
     866             :   GEN z;
     867         357 :   if (gequal0(k) || gequal(k, gen_1)) return tab;
     868         217 :   if (gsigne(k) < 0) k = gneg(k);
     869         217 :   z = cgetg(LGTAB, t_VEC);
     870         217 :   TABx0(z) = gmul(TABx0(tab), k);
     871         217 :   TABw0(z) = gmul(TABw0(tab), k);
     872         217 :   TABxp(z) = gmul(TABxp(tab), k);
     873         217 :   TABwp(z) = gmul(TABwp(tab), k);
     874         217 :   TABxm(z) = gmul(TABxm(tab), k);
     875         217 :   TABwm(z) = gmul(TABwm(tab), k);
     876         217 :   TABh(z) = rcopy(TABh(tab)); return z;
     877             : }
     878             : 
     879             : static GEN
     880         238 : expvec(GEN v, GEN ea, long prec)
     881             : {
     882         238 :   long lv = lg(v), i;
     883         238 :   GEN z = cgetg(lv, t_VEC);
     884         238 :   for (i = 1; i < lv; i++) gel(z,i) = gpow(gel(v,i),ea,prec);
     885         238 :   return z;
     886             : }
     887             : 
     888             : static GEN
     889      128643 : expscalpr(GEN vnew, GEN xold, GEN wold, GEN ea)
     890             : {
     891      128643 :   pari_sp av = avma;
     892      128643 :   return gerepileupto(av, gdiv(gmul(gmul(vnew, wold), ea), xold));
     893             : }
     894             : static GEN
     895         238 : expvecpr(GEN vnew, GEN xold, GEN wold, GEN ea)
     896             : {
     897         238 :   long lv = lg(vnew), i;
     898         238 :   GEN z = cgetg(lv, t_VEC);
     899      128762 :   for (i = 1; i < lv; i++)
     900      128524 :     gel(z,i) = expscalpr(gel(vnew,i), gel(xold,i), gel(wold,i), ea);
     901         238 :   return z;
     902             : }
     903             : 
     904             : /* here k < -1 */
     905             : static GEN
     906         119 : exptab(GEN tab, GEN k, long prec)
     907             : {
     908             :   GEN v, ea;
     909             : 
     910         119 :   if (gcmpgs(k, -2) <= 0) return tab;
     911         119 :   ea = ginv(gsubsg(-1, k));
     912         119 :   v = cgetg(LGTAB, t_VEC);
     913         119 :   TABx0(v) = gpow(TABx0(tab), ea, prec);
     914         119 :   TABw0(v) = expscalpr(TABx0(v), TABx0(tab), TABw0(tab), ea);
     915         119 :   TABxp(v) = expvec(TABxp(tab), ea, prec);
     916         119 :   TABwp(v) = expvecpr(TABxp(v), TABxp(tab), TABwp(tab), ea);
     917         119 :   TABxm(v) = expvec(TABxm(tab), ea, prec);
     918         119 :   TABwm(v) = expvecpr(TABxm(v), TABxm(tab), TABwm(tab), ea);
     919         119 :   TABh(v) = rcopy(TABh(tab));
     920         119 :   return v;
     921             : }
     922             : 
     923             : static GEN
     924         350 : init_fin(GEN b, long codeb, long m, long l, long prec)
     925             : {
     926         350 :   switch(labs(codeb))
     927             :   {
     928             :     case f_REG:
     929          77 :     case f_SING:  return inittanhsinh(m,l);
     930         112 :     case f_YSLOW: return initexpsinh(m,l);
     931          70 :     case f_YVSLO: return exptab(initexpsinh(m,l), gel(b,2), prec);
     932          42 :     case f_YFAST: return homtab(initexpexp(m,l), f_getycplx(b,l));
     933             :     /* f_YOSCS, f_YOSCC */
     934          49 :     default: return homtab(initnumsine(m,l),f_getycplx(b,l));
     935             :   }
     936             : }
     937             : 
     938             : static GEN
     939         637 : intnuminit_i(GEN a, GEN b, long m, long prec)
     940             : {
     941             :   long codea, codeb, l;
     942             :   GEN T, kma, kmb, tmp;
     943             : 
     944         637 :   if (m > 30) pari_err_OVERFLOW("intnuminit [m]");
     945         637 :   if (m < 0) pari_err_DOMAIN("intnuminit", "m", "<", gen_0, stoi(m));
     946         630 :   l = prec+EXTRAPREC;
     947         630 :   codea = transcode(a, "a");
     948         623 :   codeb = transcode(b, "b");
     949         609 :   if (labs(codea) > labs(codeb)) { swap(a, b); lswap(codea, codeb); }
     950         609 :   if (codea == f_REG)
     951             :   {
     952         329 :     T = init_fin(b, codeb, m,l,prec);
     953         329 :     switch(labs(codeb))
     954             :     {
     955          42 :       case f_YOSCS: if (gequal0(a)) break;
     956           7 :       case f_YOSCC: T = mkvec2(inittanhsinh(m,l), T);
     957             :     }
     958         329 :     return T;
     959             :   }
     960         280 :   if (codea == f_SING)
     961             :   {
     962          21 :     T = init_fin(b,codeb, m,l,prec);
     963          21 :     T = mkvec2(inittanhsinh(m,l), T);
     964          21 :     return T;
     965             :   }
     966             :   /* now a and b are infinite */
     967         259 :   if (codea * codeb > 0) return gen_0;
     968         245 :   kma = f_getycplx(a,l); codea = labs(codea);
     969         245 :   kmb = f_getycplx(b,l); codeb = labs(codeb);
     970         245 :   if (codea == f_YSLOW && codeb == f_YSLOW) return initsinhsinh(m, l);
     971         231 :   if (codea == f_YFAST && codeb == f_YFAST && gequal(kma, kmb))
     972         126 :     return homtab(initsinh(m,l), kmb);
     973         105 :   T = cgetg(3, t_VEC);
     974         105 :   switch (codea)
     975             :   {
     976             :     case f_YSLOW:
     977             :     case f_YVSLO:
     978          56 :       tmp = initexpsinh(m,l);
     979          56 :       gel(T,1) = codea == f_YSLOW? tmp: exptab(tmp, gel(a,2), prec);
     980          56 :       switch (codeb)
     981             :       {
     982          14 :         case f_YVSLO: gel(T,2) = exptab(tmp, gel(b,2), prec); return T;
     983          21 :         case f_YFAST: gel(T,2) = homtab(initexpexp(m,l), kmb); return T;
     984             :         /* YOSC[CS] */
     985          21 :         default: gel(T,2) = homtab(initnumsine(m,l), kmb); return T;
     986             :       }
     987             :       break;
     988             :     case f_YFAST:
     989          21 :       tmp = initexpexp(m, l);
     990          21 :       gel(T,1) = homtab(tmp, kma);
     991          21 :       switch (codeb)
     992             :       {
     993           7 :         case f_YFAST: gel(T,2) = homtab(tmp, kmb); return T;
     994             :         /* YOSC[CS] */
     995          14 :         default: gel(T,2) = homtab(initnumsine(m, l), kmb); return T;
     996             :       }
     997             :     default: /* YOSC[CS] */
     998          28 :       tmp = initnumsine(m, l);
     999          28 :       gel(T,1) = homtab(tmp,kma);
    1000          28 :       if (codea == f_YOSCC && codeb == f_YOSCC && !gequal(kma, kmb))
    1001          14 :         gel(T,2) = mkvec2(inittanhsinh(m,l), homtab(tmp,kmb));
    1002             :       else
    1003          14 :         gel(T,2) = homtab(tmp,kmb);
    1004          28 :       return T;
    1005             :   }
    1006             : }
    1007             : GEN
    1008         497 : intnuminit(GEN a, GEN b, long m, long prec)
    1009             : {
    1010         497 :   pari_sp av = avma;
    1011         497 :   return gerepilecopy(av, intnuminit_i(a,b,m,prec));
    1012             : }
    1013             : 
    1014             : static GEN
    1015        4718 : intnuminit0(GEN a, GEN b, GEN tab, long prec)
    1016             : {
    1017             :   long m;
    1018        4718 :   if (!tab) m = 0;
    1019        4452 :   else if (typ(tab) != t_INT)
    1020             :   {
    1021        4445 :     if (!checktab(tab)) pari_err_TYPE("intnuminit0",tab);
    1022        4445 :     return tab;
    1023             :   }
    1024             :   else
    1025           7 :     m = itos(tab);
    1026         273 :   return intnuminit(a, b, m, prec);
    1027             : }
    1028             : 
    1029             : /* Assigns the values of the function weighted by w[k] at quadrature points x[k]
    1030             :  * [replacing the weights]. Return the index of the last non-zero coeff */
    1031             : static long
    1032         252 : weight(void *E, GEN (*eval)(void *, GEN), GEN x, GEN w)
    1033             : {
    1034         252 :   long k, l = lg(x);
    1035         252 :   for (k = 1; k < l; k++) gel(w,k) = gmul(gel(w,k), eval(E, gel(x,k)));
    1036         252 :   k--; while (k >= 1) if (!gequal0(gel(w,k--))) break;
    1037         252 :   return k;
    1038             : }
    1039             : /* compute the necessary tabs, weights multiplied by f(t) */
    1040             : static GEN
    1041         126 : intfuncinit_i(void *E, GEN (*eval)(void*, GEN), GEN tab)
    1042             : {
    1043         126 :   GEN tabxp = TABxp(tab), tabwp = TABwp(tab);
    1044         126 :   GEN tabxm = TABxm(tab), tabwm = TABwm(tab);
    1045         126 :   long L, L0 = lg(tabxp);
    1046             : 
    1047         126 :   TABw0(tab) = gmul(TABw0(tab), eval(E, TABx0(tab)));
    1048         126 :   if (lg(tabxm) == 1)
    1049             :   {
    1050         126 :     TABxm(tab) = tabxm = gneg(tabxp);
    1051         126 :     TABwm(tab) = tabwm = leafcopy(tabwp);
    1052             :   }
    1053             :   /* update wp and wm in place */
    1054         126 :   L = minss(weight(E, eval, tabxp, tabwp), weight(E, eval, tabxm, tabwm));
    1055         126 :   if (L < L0)
    1056             :   { /* catch up functions whose growth at oo was not adequately described */
    1057         126 :     setlg(tabxp, L+1);
    1058         126 :     setlg(tabwp, L+1);
    1059         126 :     if (lg(tabxm) > 1) { setlg(tabxm, L+1); setlg(tabwm, L+1); }
    1060             :   }
    1061         126 :   return tab;
    1062             : }
    1063             : 
    1064             : GEN
    1065         140 : intfuncinit(void *E, GEN (*eval)(void*, GEN), GEN a, GEN b, long m, long prec)
    1066             : {
    1067         140 :   pari_sp av = avma;
    1068         140 :   GEN tab = intnuminit_i(a, b, m, prec);
    1069             : 
    1070         140 :   if (lg(tab) == 3)
    1071           7 :     pari_err_IMPL("intfuncinit with hard endpoint behavior");
    1072         259 :   if (is_fin_f(transcode(a,"intfuncinit")) ||
    1073         126 :       is_fin_f(transcode(b,"intfuncinit")))
    1074           7 :     pari_err_IMPL("intfuncinit with finite endpoints");
    1075         126 :   return gerepilecopy(av, intfuncinit_i(E, eval, tab));
    1076             : }
    1077             : 
    1078             : static GEN
    1079        4718 : intnum_i(void *E, GEN (*eval)(void*, GEN), GEN a, GEN b, GEN tab, long prec)
    1080             : {
    1081        4718 :   GEN S = gen_0, res1, res2, kma, kmb;
    1082        4718 :   long sb, sgns = 1, codea = transcode(a, "a"), codeb = transcode(b, "b");
    1083             : 
    1084        4718 :   if (codea == f_REG && typ(a) == t_VEC) a = gel(a,1);
    1085        4718 :   if (codeb == f_REG && typ(b) == t_VEC) b = gel(b,1);
    1086        4718 :   if (codea == f_REG && codeb == f_REG) return intn(E, eval, a, b, tab);
    1087         994 :   if (labs(codea) > labs(codeb)) { swap(a,b); lswap(codea,codeb); sgns = -1; }
    1088             :   /* now labs(codea) <= labs(codeb) */
    1089         994 :   if (codeb == f_SING)
    1090             :   {
    1091          21 :     if (codea == f_REG)
    1092           7 :       S = intnsing(E, eval, b, a, tab, prec), sgns = -sgns;
    1093             :     else
    1094             :     {
    1095          14 :       GEN c = gmul2n(gadd(gel(a,1), gel(b,1)), -1);
    1096          14 :       res1 = intnsing(E, eval, a, c, gel(tab,1), prec);
    1097          14 :       res2 = intnsing(E, eval, b, c, gel(tab,2), prec);
    1098          14 :       S = gsub(res1, res2);
    1099             :     }
    1100          21 :     return (sgns < 0) ? gneg(S) : S;
    1101             :   }
    1102             :   /* now b is infinite */
    1103         973 :   sb = codeb > 0 ? 1 : -1;
    1104         973 :   codeb = labs(codeb);
    1105         973 :   if (codea == f_REG && codeb != f_YOSCC
    1106         266 :       && (codeb != f_YOSCS || gequal0(a)))
    1107             :   {
    1108         266 :     S = intninfpm(E, eval, a, sb*codeb, tab);
    1109         266 :     return sgns*sb < 0 ? gneg(S) : S;
    1110             :   }
    1111         707 :   if (is_fin_f(codea))
    1112             :   { /* either codea == f_SING  or codea == f_REG and codeb = f_YOSCC
    1113             :      * or (codeb == f_YOSCS and !gequal0(a)) */
    1114             :     GEN c;
    1115          14 :     GEN pi2p = gmul(Pi2n(1,prec), f_getycplx(b, prec));
    1116          14 :     GEN pis2p = gmul2n(pi2p, -2);
    1117          14 :     c = real_i(codea == f_SING ? gel(a,1) : a);
    1118          14 :     switch(codeb)
    1119             :     {
    1120             :       case f_YOSCC: case f_YOSCS:
    1121           7 :         if (codeb == f_YOSCC) c = gadd(c, pis2p);
    1122           7 :         c = gdiv(c, pi2p);
    1123           7 :         if (sb > 0)
    1124           7 :           c = addui(1, gceil(c));
    1125             :         else
    1126           0 :           c = subiu(gfloor(c), 1);
    1127           7 :         c = gmul(pi2p, c);
    1128           7 :         if (codeb == f_YOSCC) c = gsub(c, pis2p);
    1129           7 :         break;
    1130           7 :       default: c = addui(1, gceil(c));
    1131           7 :         break;
    1132             :     }
    1133          21 :     res1 = codea==f_SING? intnsing(E, eval, a, c, gel(tab,1), prec)
    1134          21 :                         : intn    (E, eval, a, c, gel(tab,1));
    1135          14 :     res2 = intninfpm(E, eval, c, sb*codeb,gel(tab,2));
    1136          14 :     if (sb < 0) res2 = gneg(res2);
    1137          14 :     res1 = gadd(res1, res2);
    1138          14 :     return sgns < 0 ? gneg(res1) : res1;
    1139             :   }
    1140             :   /* now a and b are infinite */
    1141         693 :   if (codea * sb > 0)
    1142             :   {
    1143          14 :     if (codea > 0) pari_warn(warner, "integral from oo to oo");
    1144          14 :     if (codea < 0) pari_warn(warner, "integral from -oo to -oo");
    1145          14 :     return gen_0;
    1146             :   }
    1147         679 :   if (sb < 0) sgns = -sgns;
    1148         679 :   codea = labs(codea);
    1149         679 :   kma = f_getycplx(a, prec);
    1150         679 :   kmb = f_getycplx(b, prec);
    1151         679 :   if ((codea == f_YSLOW && codeb == f_YSLOW)
    1152         672 :    || (codea == f_YFAST && codeb == f_YFAST && gequal(kma, kmb)))
    1153         581 :     S = intninfinf(E, eval, tab);
    1154             :   else
    1155             :   {
    1156          98 :     GEN pis2 = Pi2n(-1, prec);
    1157          98 :     GEN ca = (codea == f_YOSCC)? gmul(pis2, kma): gen_0;
    1158          98 :     GEN cb = (codeb == f_YOSCC)? gmul(pis2, kmb): gen_0;
    1159          98 :     GEN c = codea == f_YOSCC ? ca : cb;
    1160          98 :     GEN SP, SN = intninfpm(E, eval, c, -sb*codea, gel(tab,1)); /*signe(a)=-sb*/
    1161          98 :     if (codea != f_YOSCC)
    1162          84 :       SP = intninfpm(E, eval, cb, sb*codeb, gel(tab,2));
    1163             :     /* codea = codeb = f_YOSCC */
    1164          14 :     else if (gequal(kma, kmb))
    1165           0 :       SP = intninfpm(E, eval, cb, sb*codeb, gel(tab,2));
    1166             :     else
    1167             :     {
    1168          14 :       tab = gel(tab,2);
    1169          14 :       SP = intninfpm(E, eval, cb, sb*codeb, gel(tab,2));
    1170          14 :       SP = gadd(SP, intn(E, eval, ca, cb, gel(tab,1)));
    1171             :     }
    1172          98 :     S = gadd(SN, SP);
    1173             :   }
    1174         679 :   if (sgns < 0) S = gneg(S);
    1175         679 :   return S;
    1176             : }
    1177             : 
    1178             : GEN
    1179        4718 : intnum(void *E, GEN (*eval)(void*, GEN), GEN a, GEN b, GEN tab, long prec)
    1180             : {
    1181        4718 :   pari_sp ltop = avma;
    1182        4718 :   long l = prec+EXTRAPREC;
    1183             :   GEN S;
    1184             : 
    1185        4718 :   tab = intnuminit0(a, b, tab, prec);
    1186        4718 :   S = intnum_i(E, eval, gprec_w(a, l), gprec_w(b, l), tab, prec);
    1187        4718 :   return gerepilecopy(ltop, gprec_wtrunc(S, prec));
    1188             : }
    1189             : 
    1190             : typedef struct auxint_s {
    1191             :   GEN a, R, pi;
    1192             :   GEN (*f)(void*, GEN);
    1193             :   GEN (*w)(GEN, long);
    1194             :   long prec;
    1195             :   void *E;
    1196             : } auxint_t;
    1197             : 
    1198             : static GEN
    1199        3675 : auxcirc(void *E, GEN t)
    1200             : {
    1201        3675 :   auxint_t *D = (auxint_t*) E;
    1202             :   GEN s, c, z;
    1203        3675 :   mpsincos(mulrr(t, D->pi), &s, &c); z = mkcomplex(c,s);
    1204        3675 :   return gmul(z, D->f(D->E, gadd(D->a, gmul(D->R, z))));
    1205             : }
    1206             : 
    1207             : GEN
    1208           7 : intcirc(void *E, GEN (*eval)(void*, GEN), GEN a, GEN R, GEN tab, long prec)
    1209             : {
    1210             :   auxint_t D;
    1211             :   GEN z;
    1212             : 
    1213           7 :   D.a = a;
    1214           7 :   D.R = R;
    1215           7 :   D.pi = mppi(prec);
    1216           7 :   D.f = eval;
    1217           7 :   D.E = E;
    1218           7 :   z = intnum(&D, &auxcirc, real_m1(prec), real_1(prec), tab, prec);
    1219           7 :   return gmul2n(gmul(R, z), -1);
    1220             : }
    1221             : 
    1222             : GEN
    1223        4536 : intnum0(GEN a, GEN b, GEN code, GEN tab, long prec)
    1224        4536 : { EXPR_WRAP(code, intnum(EXPR_ARG, a, b, tab, prec)); }
    1225             : GEN
    1226           7 : intcirc0(GEN a, GEN R, GEN code, GEN tab, long prec)
    1227           7 : { EXPR_WRAP(code, intcirc(EXPR_ARG, a, R, tab, prec)); }
    1228             : GEN
    1229         140 : intfuncinit0(GEN a, GEN b, GEN code, long m, long prec)
    1230         140 : { EXPR_WRAP(code, intfuncinit(EXPR_ARG, a, b, m, prec)); }
    1231             : 
    1232             : #if 0
    1233             : /* Two variable integration */
    1234             : 
    1235             : typedef struct auxf_s {
    1236             :   GEN x;
    1237             :   GEN (*f)(void *, GEN, GEN);
    1238             :   void *E;
    1239             : } auxf_t;
    1240             : 
    1241             : typedef struct indi_s {
    1242             :   GEN (*c)(void*, GEN);
    1243             :   GEN (*d)(void*, GEN);
    1244             :   GEN (*f)(void *, GEN, GEN);
    1245             :   void *Ec;
    1246             :   void *Ed;
    1247             :   void *Ef;
    1248             :   GEN tabintern;
    1249             :   long prec;
    1250             : } indi_t;
    1251             : 
    1252             : static GEN
    1253             : auxf(GEN y, void *E)
    1254             : {
    1255             :   auxf_t *D = (auxf_t*) E;
    1256             :   return D->f(D->E, D->x, y);
    1257             : }
    1258             : 
    1259             : static GEN
    1260             : intnumdoubintern(GEN x, void *E)
    1261             : {
    1262             :   indi_t *D = (indi_t*) E;
    1263             :   GEN c = D->c(x, D->Ec), d = D->d(x, D->Ed);
    1264             :   auxf_t A;
    1265             : 
    1266             :   A.x = x;
    1267             :   A.f = D->f;
    1268             :   A.E = D->Ef;
    1269             :   return intnum(&A, &auxf, c, d, D->tabintern, D->prec);
    1270             : }
    1271             : 
    1272             : GEN
    1273             : intnumdoub(void *Ef, GEN (*evalf)(void *, GEN, GEN), void *Ec, GEN (*evalc)(void*, GEN), void *Ed, GEN (*evald)(void*, GEN), GEN a, GEN b, GEN tabext, GEN tabint, long prec)
    1274             : {
    1275             :   indi_t E;
    1276             : 
    1277             :   E.c = evalc;
    1278             :   E.d = evald;
    1279             :   E.f = evalf;
    1280             :   E.Ec = Ec;
    1281             :   E.Ed = Ed;
    1282             :   E.Ef = Ef;
    1283             :   E.prec = prec;
    1284             :   if (typ(tabint) == t_INT)
    1285             :   {
    1286             :     GEN C = evalc(a, Ec), D = evald(a, Ed);
    1287             :     if (typ(C) != t_VEC && typ(D) != t_VEC) { C = gen_0; D = gen_1; }
    1288             :     E.tabintern = intnuminit0(C, D, tabint, prec);
    1289             :   }
    1290             :   else E.tabintern = tabint;
    1291             :   return intnum(&E, &intnumdoubintern, a, b, tabext, prec);
    1292             : }
    1293             : 
    1294             : GEN
    1295             : intnumdoub0(GEN a, GEN b, int nc, int nd, int nf, GEN tabext, GEN tabint, long prec)
    1296             : {
    1297             :   GEN z;
    1298             :   push_lex(NULL);
    1299             :   push_lex(NULL);
    1300             :   z = intnumdoub(chf, &gp_eval2, chc, &gp_eval, chd, &gp_eval, a, b, tabext, tabint, prec);
    1301             :   pop_lex(1); pop_lex(1); return z;
    1302             : }
    1303             : #endif
    1304             : 
    1305             : 
    1306             : /* The quotient-difference algorithm. Given a vector M, convert the series
    1307             :  * S = \sum_{n >= 0} M[n+1]z^n into a continued fraction.
    1308             :  * Compute the c[n] such that
    1309             :  * S = c[1] / (1 + c[2]z / (1+c[3]z/(1+...c[lim]z))),
    1310             :  * Compute A[n] and B[n] such that
    1311             :  * S = M[1]/ (1+A[1]*z+B[1]*z^2 / (1+A[2]*z+B[2]*z^2/ (1+...1/(1+A[lim\2]*z)))),
    1312             :  * Assume lim <= #M.
    1313             :  * Does not work for certain M. */
    1314             : 
    1315             : /* Given a continued fraction CF output by the quodif program,
    1316             : convert it into an Euler continued fraction A(n), B(n), where
    1317             : $1/(1+c[2]z/(1+c[3]z/(1+..c[lim]z)))
    1318             : =1/(1+A[1]*z+B[1]*z^2/(1+A[2]*z+B[2]*z^2/(1+...1/(1+A[lim\2]*z)))). */
    1319             : static GEN
    1320        3073 : contfrac_Euler(GEN CF)
    1321             : {
    1322        3073 :   long lima, limb, i, lim = lg(CF)-1;
    1323             :   GEN A, B;
    1324        3073 :   lima = lim/2;
    1325        3073 :   limb = (lim - 1)/2;
    1326        3073 :   A = cgetg(lima+1, t_VEC);
    1327        3073 :   B = cgetg(limb+1, t_VEC);
    1328        3073 :   gel (A, 1) = gel(CF, 2);
    1329        3073 :   for (i=2; i <= lima; ++i) gel(A,i) = gadd(gel(CF, 2*i), gel(CF, 2*i-1));
    1330        3073 :   for (i=1; i <= limb; ++i) gel(B,i) = gneg(gmul(gel(CF, 2*i+1), gel(CF, 2*i)));
    1331        3073 :   return mkvec2(A, B);
    1332             : }
    1333             : 
    1334             : static GEN
    1335        3290 : contfracinit_i(GEN M, long lim)
    1336             : {
    1337             :   pari_sp av;
    1338             :   GEN e, q, c;
    1339             :   long lim2;
    1340             :   long j, k;
    1341        3290 :   e = zerovec(lim);
    1342        3290 :   c = zerovec(lim+1); gel(c, 1) = gel(M, 1);
    1343        3290 :   q = cgetg(lim+1, t_VEC);
    1344        3290 :   for (k = 1; k <= lim; ++k) gel(q, k) = gdiv(gel(M, k+1), gel(M, k));
    1345        3290 :   lim2 = lim/2; av = avma;
    1346      117011 :   for (j = 1; j <= lim2; ++j)
    1347             :   {
    1348      113721 :     long l = lim - 2*j;
    1349      113721 :     gel(c, 2*j) = gneg(gel(q, 1));
    1350    10112366 :     for (k = 0; k <= l; ++k)
    1351     9998645 :       gel(e, k+1) = gsub(gadd(gel(e, k+2), gel(q, k+2)), gel(q, k+1));
    1352     9998645 :     for (k = 0; k < l; ++k)
    1353     9884924 :       gel(q, k+1) = gdiv(gmul(gel(q, k+2), gel(e, k+2)), gel(e, k+1));
    1354      113721 :     gel(c, 2*j+1) = gneg(gel(e, 1));
    1355      113721 :     if (gc_needed(av, 3))
    1356             :     {
    1357          98 :       if (DEBUGMEM>1) pari_warn(warnmem,"contfracinit, %ld/%ld",j,lim2);
    1358          98 :       gerepileall(av, 3, &e, &c, &q);
    1359             :     }
    1360             :   }
    1361        3290 :   if (odd(lim)) gel(c, lim+1) = gneg(gel(q, 1));
    1362        3290 :   return c;
    1363             : }
    1364             : 
    1365             : GEN
    1366        3094 : contfracinit(GEN M, long lim)
    1367             : {
    1368        3094 :   pari_sp ltop = avma;
    1369             :   GEN c;
    1370        3094 :   switch(typ(M))
    1371             :   {
    1372             :     case t_RFRAC:
    1373           0 :       if (lim < 0) pari_err_TYPE("contfracinit",M);
    1374           0 :       M = gadd(M, zeroser(gvar(M), lim + 2)); /*fall through*/
    1375          28 :     case t_SER: M = gtovec(M); break;
    1376           7 :     case t_POL: M = gtovecrev(M); break;
    1377        3052 :     case t_VEC: case t_COL: break;
    1378           7 :     default: pari_err_TYPE("contfracinit", M);
    1379             :   }
    1380        3087 :   if (lim < 0)
    1381          28 :     lim = lg(M)-2;
    1382        3059 :   else if (lg(M)-1 <= lim)
    1383           0 :     pari_err_COMPONENT("contfracinit", "<", stoi(lg(M)-1), stoi(lim));
    1384        3087 :   if (lim < 0) retmkvec2(cgetg(1,t_VEC),cgetg(1,t_VEC));
    1385        3073 :   c = contfracinit_i(M, lim);
    1386        3073 :   return gerepilecopy(ltop, contfrac_Euler(c));
    1387             : }
    1388             : 
    1389             : /* Evaluate at 1/tinv the nlim first terms of the continued fraction output by
    1390             :  * contfracinit. */
    1391             : /* Not stack clean */
    1392             : GEN
    1393     1808045 : contfraceval_inv(GEN CF, GEN tinv, long nlim)
    1394             : {
    1395             :   pari_sp btop;
    1396             :   long j;
    1397     1808045 :   GEN S = gen_0, S1, S2, A, B;
    1398     1808045 :   if (typ(CF) != t_VEC || lg(CF) != 3) pari_err_TYPE("contfraceval", CF);
    1399     1808045 :   A = gel(CF, 1); if (typ(A) != t_VEC) pari_err_TYPE("contfraceval", CF);
    1400     1808045 :   B = gel(CF, 2); if (typ(B) != t_VEC) pari_err_TYPE("contfraceval", CF);
    1401     1808045 :   if (nlim < 0)
    1402          14 :     nlim = lg(A)-1;
    1403     1808031 :   else if (lg(A) <= nlim)
    1404           7 :     pari_err_COMPONENT("contfraceval", ">", stoi(lg(A)-1), stoi(nlim));
    1405     1808038 :   if (lg(B)+1 <= nlim)
    1406           0 :     pari_err_COMPONENT("contfraceval", ">", stoi(lg(B)), stoi(nlim));
    1407     1808038 :   btop = avma;
    1408     1808038 :   if (nlim <= 1) return lg(A)==1? gen_0: gdiv(tinv, gadd(gel(A, 1), tinv));
    1409     1794143 :   switch(nlim % 3)
    1410             :   {
    1411             :     case 2:
    1412      639124 :       S = gdiv(gel(B, nlim-1), gadd(gel(A, nlim), tinv));
    1413      639124 :       nlim--; break;
    1414             : 
    1415             :     case 0:
    1416      606277 :       S1 = gadd(gel(A, nlim), tinv);
    1417      606277 :       S2 = gadd(gmul(gadd(gel(A, nlim-1), tinv), S1), gel(B, nlim-1));
    1418      606277 :       S = gdiv(gmul(gel(B, nlim-2), S1), S2);
    1419      606277 :       nlim -= 2; break;
    1420             :   }
    1421             :   /* nlim = 1 (mod 3) */
    1422     8272436 :   for (j = nlim; j >= 4; j -= 3)
    1423             :   {
    1424             :     GEN S3;
    1425     6478293 :     S1 = gadd(gadd(gel(A, j), tinv), S);
    1426     6478293 :     S2 = gadd(gmul(gadd(gel(A, j-1), tinv), S1), gel(B, j-1));
    1427     6478293 :     S3 = gadd(gmul(gadd(gel(A, j-2), tinv), S2), gmul(gel(B, j-2), S1));
    1428     6478293 :     S = gdiv(gmul(gel(B, j-3), S2), S3);
    1429     6478293 :     if (gc_needed(btop, 3)) S = gerepilecopy(btop, S);
    1430             :   }
    1431     1794143 :   return gdiv(tinv, gadd(gadd(gel(A, 1), tinv), S));
    1432             : }
    1433             : 
    1434             : GEN
    1435          35 : contfraceval(GEN CF, GEN t, long nlim)
    1436             : {
    1437          35 :   pari_sp ltop = avma;
    1438          35 :   return gerepileupto(ltop, contfraceval_inv(CF, ginv(t), nlim));
    1439             : }
    1440             : 
    1441             : /* MONIEN SUMMATION */
    1442             : 
    1443             : /* basic Newton, find x ~ z such that Q(x) = 0 */
    1444             : static GEN
    1445        1561 : monrefine(GEN Q, GEN QP, GEN z, long prec)
    1446             : {
    1447        1561 :   pari_sp av = avma;
    1448        1561 :   GEN pr = poleval(Q, z);
    1449             :   for(;;)
    1450             :   {
    1451             :     GEN prnew;
    1452        7021 :     z = gsub(z, gdiv(pr, poleval(QP, z)));
    1453        7021 :     prnew = poleval(Q, z);
    1454        7021 :     if (gcmp(gabs(prnew, prec), gabs(pr, prec)) >= 0) break;
    1455        5460 :     pr = prnew;
    1456        5460 :   }
    1457        1561 :   z = gprec_w(z, 2*prec-2);
    1458        1561 :   z = gsub(z, gdiv(poleval(Q, z), poleval(QP, z)));
    1459        1561 :   return gerepileupto(av, z);
    1460             : }
    1461             : 
    1462             : static GEN
    1463         217 : RX_realroots(GEN x, long prec)
    1464         217 : { return realroots(gprec_wtrunc(x,prec), NULL, prec); }
    1465             : 
    1466             : /* (real) roots of Q, assuming QP = Q' and that half the roots are close to
    1467             :  * k+1, ..., k+m, m = deg(Q)/2-1. N.B. All roots are real and >= 1 */
    1468             : static GEN
    1469         119 : monroots(GEN Q, GEN QP, long k, long prec)
    1470             : {
    1471         119 :   long j, n = degpol(Q), m = n/2 - 1;
    1472         119 :   GEN v2, v1 = cgetg(m+1, t_VEC);
    1473         119 :   for (j = 1; j <= m; ++j) gel(v1, j) = monrefine(Q, QP, stoi(k+j), prec);
    1474         119 :   Q = gdivent(Q, roots_to_pol(v1, varn(Q)));
    1475         119 :   v2 = RX_realroots(Q, prec); settyp(v2, t_VEC);
    1476         119 :   return shallowconcat(v1, v2);
    1477             : }
    1478             : 
    1479             : static void
    1480         217 : Pade(GEN M, GEN *pP, GEN *pQ)
    1481             : {
    1482         217 :   pari_sp av = avma;
    1483         217 :   long n = lg(M)-2, i;
    1484         217 :   GEN v = contfracinit_i(M, n), P = pol_0(0), Q = pol_1(0);
    1485             :   /* evaluate continued fraction => Pade approximants */
    1486       12803 :   for (i = n-1; i >= 1; i--)
    1487             :   { /* S = P/Q: S -> v[i]*x / (1+S) */
    1488       12586 :     GEN R = RgX_shift_shallow(RgX_Rg_mul(Q,gel(v,i)), 1);
    1489       12586 :     Q = RgX_add(P,Q); P = R;
    1490       12586 :     if (gc_needed(av, 3))
    1491             :     {
    1492           0 :       if (DEBUGMEM>1) pari_warn(warnmem,"Pade, %ld/%ld",i,n-1);
    1493           0 :       gerepileall(av, 3, &P, &Q, &v);
    1494             :     }
    1495             :   }
    1496             :   /* S -> 1+S */
    1497         217 :   *pP = RgX_add(P,Q);
    1498         217 :   *pQ = Q;
    1499         217 : }
    1500             : 
    1501             : static GEN
    1502           7 : veczetaprime(GEN a, GEN b, long N, long prec)
    1503             : {
    1504           7 :   long newprec, fpr = prec2nbits(prec), pr = (long)ceil(fpr * 1.5);
    1505           7 :   long l = nbits2prec(pr), e = fpr / 2;
    1506             :   GEN eps, A, B;
    1507           7 :   newprec = nbits2prec(pr + BITS_IN_LONG);
    1508           7 :   a = gprec_w(a, newprec);
    1509           7 :   b = gprec_w(b, newprec);
    1510           7 :   eps = real2n(-e, l);
    1511           7 :   A = veczeta(a, gsub(b, eps), N, newprec);
    1512           7 :   B = veczeta(a, gadd(b, eps), N, newprec);
    1513           7 :   return gmul2n(RgV_sub(B, A), e-1);
    1514             : }
    1515             : 
    1516             : struct mon_w {
    1517             :   GEN w, a, b;
    1518             :   long n, j, prec;
    1519             : };
    1520             : 
    1521             : /* w(x) / x^(a*(j+k)+b), k >= 1; w a t_CLOSURE or t_INT [encodes log(x)^w] */
    1522             : static GEN
    1523       34384 : wrapmonw(void* E, GEN x)
    1524             : {
    1525       34384 :   struct mon_w *W = (struct mon_w*)E;
    1526       34384 :   long k, j = W->j, n = W->n, prec = W->prec, l = 2*n+4-j;
    1527      103152 :   GEN wx = typ(W->w) == t_CLOSURE? closure_callgen1prec(W->w, x, prec)
    1528       68768 :                                  : powgi(glog(x, prec), W->w);
    1529       34384 :   GEN v = cgetg(l, t_VEC);
    1530       34384 :   GEN xa = gpow(x, gneg(W->a), prec), w = gmul(wx, gpowgs(xa, j));
    1531       34384 :   w = gdiv(w, gpow(x,W->b,prec));
    1532       34384 :   for (k = 1; k < l; k++) { gel(v,k) = w; w = gmul(w, xa); }
    1533       34384 :   return v;
    1534             : }
    1535             : /* w(x) / x^(a*j+b) */
    1536             : static GEN
    1537       18819 : wrapmonw2(void* E, GEN x)
    1538             : {
    1539       18819 :   struct mon_w *W = (struct mon_w*)E;
    1540       18819 :   GEN wnx = closure_callgen1prec(W->w, x, W->prec);
    1541       18819 :   return gdiv(wnx, gpow(x, gadd(gmulgs(W->a, W->j), W->b), W->prec));
    1542             : }
    1543             : static GEN
    1544          35 : M_from_wrapmon(struct mon_w *S, GEN wfast, GEN n0)
    1545             : {
    1546          35 :   long j, N = 2*S->n+2;
    1547          35 :   GEN M = cgetg(N+1, t_VEC), faj = gsub(wfast, S->b);
    1548          42 :   for (j = 1; j <= N; j++)
    1549             :   {
    1550          42 :     faj = gsub(faj, S->a);
    1551          42 :     if (gcmpgs(faj, -2) <= 0)
    1552             :     {
    1553          28 :       S->j = j; setlg(M,j);
    1554          28 :       M = shallowconcat(M, sumnum((void*)S, wrapmonw, n0, NULL, S->prec));
    1555          28 :       break;
    1556             :     }
    1557          14 :     S->j = j;
    1558          14 :     gel(M,j) = sumnum((void*)S, wrapmonw2, mkvec2(n0,faj), NULL, S->prec);
    1559             :   }
    1560          28 :   return M;
    1561             : }
    1562             : 
    1563             : static void
    1564         168 : checkmonroots(GEN vr, long n)
    1565             : {
    1566         168 :   if (lg(vr) != n+1)
    1567           0 :     pari_err_IMPL("recovery when missing roots in sumnummonieninit");
    1568         168 : }
    1569             : 
    1570             : static GEN
    1571         175 : sumnummonieninit_i(GEN a, GEN b, GEN w, GEN n0, long prec)
    1572             : {
    1573         175 :   GEN c, M, P, Q, Qp, vr, vabs, vwt, ga = gadd(a, b);
    1574         175 :   double bit = 2*prec2nbits(prec) / gtodouble(ga), D = bit*LOG2;
    1575         175 :   double da = maxdd(1., gtodouble(a));
    1576         175 :   long n = (long)ceil(D / (da*(log(D)-1)));
    1577         175 :   long j, prec2, prec0 = prec + EXTRAPRECWORD;
    1578         175 :   double bit0 = ceil((2*n+1)*LOG2_10);
    1579         175 :   int neg = 1;
    1580             :   struct mon_w S;
    1581             : 
    1582             :   /* 2.05 is heuristic; with 2.03, sumnummonien(n=1,1/n^2) loses
    1583             :    * 19 decimals at \p1500 */
    1584         175 :   prec = nbits2prec(maxdd(2.05*da*bit, bit0));
    1585         175 :   prec2 = nbits2prec(maxdd(1.3*da*bit, bit0));
    1586         175 :   S.w = w;
    1587         175 :   S.a = a = gprec_w(a, 2*prec-2);
    1588         175 :   S.b = b = gprec_w(b, 2*prec-2);
    1589         175 :   S.n = n;
    1590         175 :   S.j = 1;
    1591         175 :   S.prec = prec;
    1592         175 :   if (typ(w) == t_INT)
    1593             :   { /* f(n) ~ \sum_{i > 0} f_i log(n)^k / n^(a*i + b); a > 0, a+b > 1 */
    1594         140 :     long k = itos(w);
    1595         140 :     if (k == 0)
    1596         133 :       M = veczeta(a, ga, 2*n+2, prec);
    1597           7 :     else if (k == 1)
    1598           7 :       M = veczetaprime(a, ga, 2*n+2, prec);
    1599             :     else
    1600           0 :       M = M_from_wrapmon(&S, gen_0, gen_1);
    1601         140 :     if (odd(k)) neg = 0;
    1602             :   }
    1603             :   else
    1604             :   {
    1605          35 :     GEN wfast = gen_0;
    1606          35 :     if (typ(w) == t_VEC) { wfast = gel(w,2); w = gel(w,1); }
    1607          35 :     M = M_from_wrapmon(&S, wfast, n0);
    1608             :   }
    1609             :   /* M[j] = sum(n >= n0, w(n) / n^(a*j+b) */
    1610         168 :   Pade(M, &P,&Q);
    1611         168 :   Qp = RgX_deriv(Q);
    1612         168 :   if (gequal1(a)) a = NULL;
    1613         168 :   if (!a && typ(w) == t_INT)
    1614             :   {
    1615         119 :     vabs = vr = monroots(Q, Qp, signe(w)? 1: 0, prec2);
    1616         119 :     checkmonroots(vr, n);
    1617         119 :     c = b;
    1618             :   }
    1619             :   else
    1620             :   {
    1621          49 :     vr = RX_realroots(Q, prec2); settyp(vr, t_VEC);
    1622          49 :     checkmonroots(vr, n);
    1623          49 :     if (!a) { vabs = vr; c = b; }
    1624             :     else
    1625             :     {
    1626          35 :       GEN ai = ginv(a);
    1627          35 :       vabs = cgetg(n+1, t_VEC);
    1628          35 :       for (j = 1; j <= n; j++) gel(vabs,j) = gpow(gel(vr,j), ai, prec2);
    1629          35 :       c = gdiv(b,a);
    1630             :     }
    1631             :   }
    1632             : 
    1633         168 :   vwt = cgetg(n+1, t_VEC);
    1634         168 :   c = gsubgs(c,1); if (gequal0(c)) c = NULL;
    1635        5012 :   for (j = 1; j <= n; j++)
    1636             :   {
    1637        4844 :     GEN r = gel(vr,j), t = gdiv(poleval(P,r), poleval(Qp,r));
    1638        4844 :     if (c) t = gmul(t, gpow(r, c, prec));
    1639        4844 :     gel(vwt,j) = neg? gneg(t): t;
    1640             :   }
    1641         168 :   if (typ(w) == t_INT && !equali1(n0))
    1642             :   {
    1643          28 :     GEN h = subiu(n0,1);
    1644          28 :     for (j = 1; j <= n; j++) gel(vabs,j) = gadd(gel(vabs,j), h);
    1645             :   }
    1646         168 :   return mkvec3(gprec_wtrunc(vabs,prec0), gprec_wtrunc(vwt,prec0), n0);
    1647             : }
    1648             : 
    1649             : GEN
    1650         105 : sumnummonieninit(GEN asymp, GEN w, GEN n0, long prec)
    1651             : {
    1652         105 :   pari_sp av = avma;
    1653         105 :   const char *fun = "sumnummonieninit";
    1654             :   GEN a, b;
    1655         105 :   if (!n0) n0 = gen_1; else if (typ(n0) != t_INT) pari_err_TYPE(fun, n0);
    1656         105 :   if (!asymp) a = b = gen_1;
    1657             :   else
    1658             :   {
    1659          77 :     if (typ(asymp) == t_VEC)
    1660             :     {
    1661          70 :       if (lg(asymp) != 3) pari_err_TYPE(fun, asymp);
    1662          70 :       a = gel(asymp,1);
    1663          70 :       b = gel(asymp,2);
    1664             :     }
    1665             :     else
    1666             :     {
    1667           7 :       a = gen_1;
    1668           7 :       b = asymp;
    1669             :     }
    1670          77 :     if (gsigne(a) <= 0) pari_err_DOMAIN(fun, "a", "<=", gen_0, a);
    1671          70 :     if (gcmpgs(gadd(a,b), 1) <= 0)
    1672           7 :       pari_err_DOMAIN(fun, "a+b", "<=", gen_1, mkvec2(a,b));
    1673             :   }
    1674          91 :   if (!w) w = gen_0;
    1675          42 :   else switch(typ(w))
    1676             :   {
    1677             :     case t_INT:
    1678           7 :       if (signe(w) < 0) pari_err_IMPL("log power < 0 in sumnummonieninit");
    1679          35 :     case t_CLOSURE: break;
    1680             :     case t_VEC:
    1681           7 :       if (lg(w) == 3 && typ(gel(w,1)) == t_CLOSURE) break;
    1682           0 :     default: pari_err_TYPE(fun, w);
    1683             :   }
    1684          91 :   return gerepilecopy(av, sumnummonieninit_i(a, b, w, n0, prec));
    1685             : }
    1686             : 
    1687             : GEN
    1688         175 : sumnummonien(void *E, GEN (*eval)(void*,GEN), GEN n0, GEN tab, long prec)
    1689             : {
    1690         175 :   pari_sp av = avma;
    1691             :   GEN vabs, vwt, S;
    1692             :   long l, i;
    1693         175 :   if (typ(n0) != t_INT) pari_err_TYPE("sumnummonien", n0);
    1694         175 :   if (!tab)
    1695          84 :     tab = sumnummonieninit_i(gen_1, gen_1, gen_0, n0, prec);
    1696             :   else
    1697             :   {
    1698          91 :     if (lg(tab) != 4 || typ(tab) != t_VEC) pari_err_TYPE("sumnummonien", tab);
    1699          91 :     if (!equalii(n0, gel(tab,3)))
    1700           7 :       pari_err(e_MISC, "incompatible initial value %Ps != %Ps", gel(tab,3),n0);
    1701             :   }
    1702         168 :   vabs= gel(tab,1); l = lg(vabs);
    1703         168 :   vwt = gel(tab,2);
    1704         168 :   if (typ(vabs) != t_VEC || typ(vwt) != t_VEC || lg(vwt) != l)
    1705           0 :     pari_err_TYPE("sumnummonien", tab);
    1706         168 :   S = gen_0;
    1707         168 :   for (i = 1; i < l; i++) S = gadd(S, gmul(gel(vwt,i), eval(E, gel(vabs,i))));
    1708         168 :   return gerepileupto(av, gprec_w(S, prec));
    1709             : }
    1710             : 
    1711             : static GEN
    1712         196 : get_oo(GEN fast) { return mkvec2(mkoo(), fast); }
    1713             : 
    1714             : GEN
    1715         119 : sumnuminit(GEN fast, long prec)
    1716             : {
    1717             :   pari_sp av;
    1718         119 :   GEN s, v, d, C, res = cgetg(6, t_VEC);
    1719         119 :   long bitprec = prec2nbits(prec), N, k, k2, m;
    1720             :   double w;
    1721             : 
    1722         119 :   gel(res, 1) = d = mkfrac(gen_1, utoipos(4)); /* 1/4 */
    1723         119 :   av = avma;
    1724         119 :   w = gtodouble(glambertW(ginv(d), LOWDEFAULTPREC));
    1725         119 :   N = (long)ceil(LOG2*bitprec/(w*(1+w))+5);
    1726         119 :   k = (long)ceil(N*w); if (k&1) k--;
    1727             : 
    1728         119 :   prec += EXTRAPRECWORD;
    1729         119 :   k2 = k/2;
    1730         119 :   s = RgX_to_ser(monomial(real_1(prec),1,0), k+3);
    1731         119 :   s = gmul2n(gasinh(s, prec), 2); /* asinh(x)/d, d = 1/4 */
    1732         119 :   gel(s, 2) = utoipos(4);
    1733         119 :   s = gsub(ser_inv(gexpm1(s,prec)), ser_inv(s));
    1734         119 :   C = matpascal(k-1);
    1735         119 :   v = cgetg(k2+1, t_VEC);
    1736        8449 :   for (m = 1; m <= k2; m++)
    1737             :   {
    1738        8330 :     pari_sp av = avma;
    1739        8330 :     GEN S = real_0(prec);
    1740             :     long j;
    1741      484169 :     for (j = m; j <= k2; j++)
    1742             :     { /* s[X^(2j-1)] * binomial(2*j-1, j-m) */
    1743      475839 :       GEN t = gmul(gel(s,2*j+1), gcoeff(C, 2*j,j-m+1));
    1744      475839 :       t = gmul2n(t, 1-2*j);
    1745      475839 :       S = odd(j)? gsub(S,t): gadd(S,t);
    1746             :     }
    1747        8330 :     if (odd(m)) S = gneg(S);
    1748        8330 :     gel(v,m) = gerepileupto(av, S);
    1749             :   }
    1750         119 :   v = RgC_gtofp(v,prec); settyp(v, t_VEC);
    1751         119 :   gel(res, 4) = gerepileupto(av, v);
    1752         119 :   gel(res, 2) = utoi(N);
    1753         119 :   gel(res, 3) = utoi(k);
    1754         119 :   if (!fast) fast = get_oo(gen_0);
    1755         119 :   gel(res, 5) = intnuminit(gel(res,2), fast, 0, prec - EXTRAPRECWORD);
    1756         119 :   return res;
    1757             : }
    1758             : 
    1759             : static int
    1760          28 : checksumtab(GEN T)
    1761             : {
    1762          28 :   if (typ(T) != t_VEC || lg(T) != 6) return 0;
    1763          21 :   return typ(gel(T,2))==t_INT && typ(gel(T,3))==t_INT && typ(gel(T,4))==t_VEC;
    1764             : }
    1765             : GEN
    1766         133 : sumnum(void *E, GEN (*eval)(void*, GEN), GEN a, GEN tab, long prec)
    1767             : {
    1768         133 :   pari_sp av = avma, av2;
    1769             :   GEN v, tabint, S, d, fast;
    1770             :   long as, N, k, m, prec2;
    1771         133 :   if (!a) { a = gen_1; fast = get_oo(gen_0); }
    1772         133 :   else switch(typ(a))
    1773             :   {
    1774             :   case t_VEC:
    1775          49 :     if (lg(a) != 3) pari_err_TYPE("sumnum", a);
    1776          49 :     fast = get_oo(gel(a,2));
    1777          49 :     a = gel(a,1); break;
    1778             :   default:
    1779          84 :     fast = get_oo(gen_0);
    1780             :   }
    1781         133 :   if (typ(a) != t_INT) pari_err_TYPE("sumnum", a);
    1782         133 :   if (!tab) tab = sumnuminit(fast, prec);
    1783          28 :   else if (!checksumtab(tab)) pari_err_TYPE("sumnum",tab);
    1784         126 :   as = itos(a);
    1785         126 :   d = gel(tab,1);
    1786         126 :   N = maxss(as, itos(gel(tab,2)));
    1787         126 :   k = itos(gel(tab,3));
    1788         126 :   v = gel(tab,4);
    1789         126 :   tabint = gel(tab,5);
    1790         126 :   prec2 = prec+EXTRAPRECWORD;
    1791         126 :   av2 = avma;
    1792         126 :   S = gmul(eval(E, stoi(N)), real2n(-1,prec2));
    1793       15765 :   for (m = as; m < N; m++)
    1794             :   {
    1795       15639 :     S = gadd(S, eval(E, stoi(m)));
    1796       15639 :     if (gc_needed(av, 2))
    1797             :     {
    1798           0 :       if (DEBUGMEM>1) pari_warn(warnmem,"sumnum [1], %ld/%ld",m,N);
    1799           0 :       S = gerepileupto(av2, S);
    1800             :     }
    1801             :   }
    1802        9611 :   for (m = 1; m <= k/2; m++)
    1803             :   {
    1804        9485 :     GEN t = gmulsg(2*m-1, d);
    1805        9485 :     GEN s = gsub(eval(E, gsubsg(N,t)), eval(E, gaddsg(N,t)));
    1806        9485 :     S = gadd(S, gmul(gel(v,m), s));
    1807        9485 :     if (gc_needed(av2, 2))
    1808             :     {
    1809           0 :       if (DEBUGMEM>1) pari_warn(warnmem,"sumnum [2], %ld/%ld",m,k/2);
    1810           0 :       S = gerepileupto(av2, S);
    1811             :     }
    1812             :   }
    1813         126 :   S = gadd(S, intnum(E, eval,stoi(N), fast, tabint, prec2));
    1814         126 :   return gerepilecopy(av, gprec_w(S, prec));
    1815             : }
    1816             : 
    1817             : GEN
    1818         175 : sumnummonien0(GEN a, GEN code, GEN tab, long prec)
    1819         175 : { EXPR_WRAP(code, sumnummonien(EXPR_ARG, a, tab, prec)); }
    1820             : GEN
    1821          91 : sumnum0(GEN a, GEN code, GEN tab, long prec)
    1822          91 : { EXPR_WRAP(code, sumnum(EXPR_ARG, a, tab, prec)); }
    1823             : 
    1824             : /* Abel-Plana summation */
    1825             : 
    1826             : static GEN
    1827          49 : intnumgauexpinit(long prec)
    1828             : {
    1829          49 :   pari_sp ltop = avma;
    1830             :   GEN V, N, E, P, Q, R, vabs, vwt;
    1831          49 :   long l, n, k, j, prec2, prec0 = prec + EXTRAPRECWORD, bit = prec2nbits(prec);
    1832             : 
    1833          49 :   n = (long)ceil(bit*0.226);
    1834          49 :   n |= 1; /* make n odd */
    1835          49 :   prec = nbits2prec(1.5*bit + 32);
    1836          49 :   prec2 = maxss(prec0, nbits2prec(1.15*bit + 32));
    1837          49 :   V = cgetg(n + 4, t_VEC);
    1838        3045 :   for (k = 1; k <= n + 3; ++k)
    1839        2996 :     gel(V, k) = gtofp(gdivgs(bernfrac(2*k), odd(k)? 2*k: -2*k), prec);
    1840          49 :   Pade(V, &P, &Q);
    1841          49 :   N = RgX_recip(gsub(P, Q));
    1842          49 :   E = RgX_recip(Q);
    1843          49 :   R = gdivgs(gdiv(N, RgX_deriv(E)), 2);
    1844          49 :   vabs = RX_realroots(E,prec2);
    1845          49 :   l = lg(vabs); settyp(vabs, t_VEC);
    1846          49 :   vwt = cgetg(l, t_VEC);
    1847        1498 :   for (j = 1; j < l; ++j)
    1848             :   {
    1849        1449 :     GEN a = gel(vabs,j);
    1850        1449 :     gel(vwt, j) = gprec_wtrunc(poleval(R, a), prec0);
    1851        1449 :     gel(vabs, j) = gprec_wtrunc(sqrtr_abs(a), prec0);
    1852             :   }
    1853          49 :   return gerepilecopy(ltop, mkvec2(vabs, vwt));
    1854             : }
    1855             : 
    1856             : /* Compute $\int_{-\infty}^\infty w(x)f(x)\,dx$, where $w(x)=x/(exp(2\pi x)-1)$
    1857             :  * for $x>0$ and $w(-x)=w(x)$. For Abel-Plana (sumnumap). */
    1858             : static GEN
    1859          49 : intnumgauexp(void *E, GEN (*eval)(void*,GEN), GEN gN, GEN tab, long prec)
    1860             : {
    1861          49 :   pari_sp av = avma;
    1862          49 :   GEN U = mkcomplex(gN, gen_0), V = mkcomplex(gel(U,1), gen_0), S = gen_0;
    1863          49 :   GEN vabs = gel(tab, 1), vwt = gel(tab, 2);
    1864          49 :   long l = lg(vabs), i;
    1865          49 :   if (lg(vwt) != l || typ(vabs) != t_VEC || typ(vwt) != t_VEC)
    1866           0 :     pari_err_TYPE("sumnumap", tab);
    1867        1498 :   for (i = 1; i < l; i++)
    1868             :   { /* I * (w_i/a_i) * (f(N + I*a_i) - f(N - I*a_i)) */
    1869        1449 :     GEN x = gel(vabs,i), w = gel(vwt,i), t;
    1870        1449 :     gel(U,2) = x;
    1871        1449 :     gel(V,2) = gneg(x);
    1872        1449 :     t = mulcxI(gsub(eval(E,U), eval(E,V)));
    1873        1449 :     if (typ(t) == t_COMPLEX && gequal0(gel(t,2))) t = gel(t,1);
    1874        1449 :     S = gadd(S, gmul(gdiv(w,x), t));
    1875             :   }
    1876          49 :   return gerepileupto(av, gprec_w(S, prec));
    1877             : }
    1878             : 
    1879             : GEN
    1880          49 : sumnumapinit(GEN fast, long prec)
    1881             : {
    1882          49 :   if (!fast) fast = mkoo();
    1883          49 :   retmkvec2(intnumgauexpinit(prec), intnuminit(gen_1, fast, 0, prec));
    1884             : }
    1885             : 
    1886             : typedef struct {
    1887             :   GEN (*f)(void *E, GEN);
    1888             :   void *E;
    1889             :   long N;
    1890             : } expfn;
    1891             : 
    1892             : /* f(Nx) */
    1893             : static GEN
    1894       31157 : _exfn(void *E, GEN x)
    1895             : {
    1896       31157 :   expfn *S = (expfn*)E;
    1897       31157 :   return S->f(S->E, gmulsg(S->N, x));
    1898             : }
    1899             : 
    1900             : GEN
    1901          56 : sumnumap(void *E, GEN (*eval)(void*,GEN), GEN a, GEN tab, long prec)
    1902             : {
    1903          56 :   pari_sp av = avma;
    1904             :   expfn T;
    1905             :   GEN S, fast, gN;
    1906             :   long as, m, N;
    1907          56 :   if (!a) { a = gen_1; fast = get_oo(gen_0); }
    1908          56 :   else switch(typ(a))
    1909             :   {
    1910             :     case t_VEC:
    1911          21 :       if (lg(a) != 3) pari_err_TYPE("sumnumap", a);
    1912          21 :       fast = get_oo(gel(a,2));
    1913          21 :       a = gel(a,1); break;
    1914             :     default:
    1915          35 :       fast = get_oo(gen_0);
    1916             :   }
    1917          56 :   if (typ(a) != t_INT) pari_err_TYPE("sumnumap", a);
    1918          56 :   if (!tab) tab = sumnumapinit(fast, prec);
    1919          14 :   else if (typ(tab) != t_VEC || lg(tab) != 3) pari_err_TYPE("sumnumap",tab);
    1920          49 :   as = itos(a);
    1921          49 :   T.N = N = maxss(as + 1, (long)ceil(prec2nbits(prec)*0.327));
    1922          49 :   T.E = E;
    1923          49 :   T.f = eval;
    1924          49 :   gN = stoi(N);
    1925          49 :   S = gtofp(gmul2n(eval(E, gN), -1), prec);
    1926          49 :   for (m = as; m < N; ++m) S = gadd(S, eval(E, stoi(m)));
    1927          49 :   S = gadd(S, gmulsg(N, intnum(&T, &_exfn, gen_1, fast, gel(tab, 2), prec)));
    1928          49 :   S = gadd(S, intnumgauexp(E, eval, gN, gel(tab, 1), prec));
    1929          49 :   return gerepileupto(av, S);
    1930             : }
    1931             : 
    1932             : GEN
    1933          56 : sumnumap0(GEN a, GEN code, GEN tab, long prec)
    1934          56 : { EXPR_WRAP(code, sumnumap(EXPR_ARG, a, tab, prec)); }
    1935             : 
    1936             : 
    1937             : /* max (1, |zeros|), P a t_POL or scalar */
    1938             : static GEN
    1939         168 : polmax(GEN P)
    1940             : {
    1941             :   GEN r;
    1942         168 :   if (typ(P) != t_POL || degpol(P) <= 0) return gen_1;
    1943         168 :   r = polrootsbound(P);
    1944         168 :   if (gcmpgs(r, 1) <= 0) return gen_1;
    1945         119 :   return r;
    1946             : }
    1947             : 
    1948             : /* max (1, |poles|), F a t_POL or t_RFRAC or scalar */
    1949             : static GEN
    1950          84 : ratpolemax(GEN F)
    1951             : {
    1952          84 :   if (typ(F) == t_POL) return gen_1;
    1953          84 :   return polmax(gel(F, 2));
    1954             : }
    1955             : /* max (1, |poles|, |zeros|), sets *p = max(1, |poles|)) */
    1956             : static GEN
    1957          42 : ratpolemax2(GEN F, GEN *p)
    1958             : {
    1959             :   GEN t;
    1960          42 :   if (typ(F) == t_POL) { if (p) *p = gen_1; return polmax(F); }
    1961          42 :   t = polmax(gel(F,2)); if (p) *p = t;
    1962          42 :   return gmax(t, polmax(gel(F,1)));
    1963             : }
    1964             : 
    1965             : static GEN
    1966       19593 : sercoeff(GEN x, long n)
    1967             : {
    1968       19593 :   long N = n - valp(x);
    1969       19593 :   return (N < 0)? gen_0: gel(x,N+2);
    1970             : }
    1971             : 
    1972             : /* Compute the integral from N to infinity of a rational function F, deg(F) < -1
    1973             :  * We must have N > 2 * r, r = max(1, |poles F|). */
    1974             : static GEN
    1975          70 : intnumainfrat(GEN F, long N, double r, long prec)
    1976             : {
    1977          70 :   long B = prec2nbits(prec), v, k, lim;
    1978             :   GEN S, ser;
    1979          70 :   pari_sp av = avma;
    1980             : 
    1981          70 :   lim = (long)ceil(B/log2(N/r)) + 1;
    1982          70 :   ser = gmul(F, real_1(prec + EXTRAPREC));
    1983          70 :   ser = rfracrecip_to_ser_absolute(ser, lim);
    1984          70 :   v = valp(ser);
    1985          70 :   S = gdivgs(sercoeff(ser,lim+1), lim*N);
    1986             :   /* goes down to 2, but coeffs are 0 in degree < v */
    1987        4193 :   for (k = lim; k >= v; k--) /* S <- (S + coeff(ser,k)/(k-1)) / N */
    1988        4123 :     S = gdivgs(gadd(S, gdivgs(sercoeff(ser,k), k-1)), N);
    1989          70 :   if (v-2) S = gdiv(S, powuu(N, v-2));
    1990          70 :   return gerepileupto(av, gprec_w(S, prec));
    1991             : }
    1992             : 
    1993             : static GEN
    1994         112 : rfrac_eval0(GEN R)
    1995             : {
    1996         112 :   GEN N, n, D = gel(R,2), d = constant_coeff(D);
    1997         112 :   if (gcmp0(d)) return NULL;
    1998          84 :   N = gel(R,1);
    1999          84 :   n = typ(N)==t_POL? constant_coeff(N): N;
    2000          84 :   return gdiv(n, d);
    2001             : }
    2002             : static GEN
    2003        6356 : rfrac_eval(GEN R, GEN a)
    2004             : {
    2005        6356 :   GEN D = gel(R,2), d = poleval(D,a);
    2006        6356 :   return gcmp0(d)? NULL: gdiv(poleval(gel(R,1),a), d);
    2007             : }
    2008             : /* R = \sum_i vR[i], eval at a omitting poles */
    2009             : static GEN
    2010        6286 : RFRAC_eval(GEN R, GEN vR, GEN a)
    2011             : {
    2012        6286 :   GEN S = rfrac_eval(R,a);
    2013        6286 :   if (!S && vR)
    2014             :   {
    2015          14 :     long i, l = lg(vR);
    2016          56 :     for (i = 1; i < l; i++)
    2017             :     {
    2018          42 :       GEN z = rfrac_eval(gel(vR,i), a);
    2019          42 :       if (z) S = S? gadd(S,z): z;
    2020             :     }
    2021             :   }
    2022        6286 :   return S;
    2023             : }
    2024             : static GEN
    2025        6286 : add_RFRAC_eval(GEN S, GEN R, GEN vR, GEN a)
    2026             : {
    2027        6286 :   GEN z = RFRAC_eval(R, vR, a);
    2028        6286 :   return z? gadd(S, z): S;
    2029             : }
    2030             : static GEN
    2031          63 : add_sumrfrac(GEN S, GEN R, GEN vR, long b)
    2032             : {
    2033             :   long m;
    2034          63 :   for (m = b; m >= 1; m--) S = add_RFRAC_eval(S,R,vR,utoipos(m));
    2035          63 :   return S;
    2036             : }
    2037             : static void
    2038          70 : get_kN(long r, long B, long *pk, long *pN)
    2039             : {
    2040          70 :   long k = maxss(50, (long)ceil(0.241*B));
    2041             :   GEN z;
    2042          70 :   if (k&1L) k++;
    2043          70 :   *pk = k;
    2044          70 :   z = sqrtnr_abs(gmul2n(gtofp(bernfrac(k), LOWDEFAULTPREC), B), k);
    2045          70 :   *pN = maxss(2*r, r + 1 + itos(gceil(z)));
    2046          70 : }
    2047             : /* F a t_RFRAC, F0 = F(0) or NULL [pole], vF a vector of t_RFRAC summing to F
    2048             :  * or NULL [F atomic] */
    2049             : static GEN
    2050          70 : sumnumrat_i(GEN F, GEN F0, GEN vF, long prec)
    2051             : {
    2052          70 :   long B = prec2nbits(prec), vx, j, k, N;
    2053             :   GEN S, S1, S2, intf;
    2054             :   double r;
    2055          70 :   if (poldegree(F, -1) > -2) pari_err(e_MISC, "sum diverges in sumnumrat");
    2056          63 :   vx = varn(gel(F,2));
    2057          63 :   r = gtodouble(ratpolemax(F));
    2058          63 :   get_kN((long)ceil(r), B, &k,&N);
    2059          63 :   intf = intnumainfrat(F, N, r, prec);
    2060             :   /* N > ratpolemax(F) is not a pole */
    2061          63 :   S1 = gmul2n(gtofp(poleval(F, utoipos(N)), prec), -1);
    2062          63 :   S1 = add_sumrfrac(S1, F, vF, N-1);
    2063          63 :   if (F0) S1 = gadd(S1, F0);
    2064          63 :   S = gmul(real_1(prec), gsubst(F, vx, gaddgs(pol_x(vx), N)));
    2065          63 :   S = rfrac_to_ser(S, k + 2);
    2066          63 :   S2 = gen_0;
    2067        3024 :   for (j = 2; j <= k; j += 2)
    2068        2961 :     S2 = gadd(S2, gmul(gdivgs(bernfrac(j),j), sercoeff(S, j-1)));
    2069          63 :   return gadd(intf, gsub(S1, S2));
    2070             : }
    2071             : /* sum_{n >= a} F(n) */
    2072             : GEN
    2073          56 : sumnumrat(GEN F, GEN a, long prec)
    2074             : {
    2075          56 :   pari_sp av = avma;
    2076             :   long vx;
    2077             :   GEN vF, F0;
    2078             : 
    2079          56 :   switch(typ(F))
    2080             :   {
    2081          42 :     case t_RFRAC: break;
    2082             :     case t_INT: case t_REAL: case t_COMPLEX: case t_POL:
    2083          14 :       if (gequal0(F)) return real_0(prec);
    2084           7 :     default: pari_err_TYPE("sumnumrat",F);
    2085             :   }
    2086          42 :   vx = varn(gel(F,2));
    2087          42 :   switch(typ(a))
    2088             :   {
    2089             :     case t_INT:
    2090          21 :       if (signe(a)) F = gsubst(F, vx, deg1pol_shallow(gen_1,a,vx));
    2091          21 :       F0 = rfrac_eval0(F);
    2092          21 :       vF = NULL;
    2093          21 :       break;
    2094             :     case t_INFINITY:
    2095          21 :       if (inf_get_sign(a) == -1)
    2096             :       { /* F(x) + F(-x). Could divide degree by 2, as G(x^2): pb with poles */
    2097          14 :         GEN F2 = gsubst(F, vx, RgX_neg(pol_x(vx)));
    2098          14 :         vF = mkvec2(F,F2);
    2099          14 :         F = gadd(F, F2);
    2100          14 :         if (gequal0(F)) { avma = av; return real_0(prec); }
    2101           7 :         F0 = rfrac_eval0(gel(vF,1));
    2102           7 :         break;
    2103             :       }
    2104             :     default:
    2105           7 :       pari_err_TYPE("sumnumrat",a);
    2106             :       return NULL; /* LCOV_EXCL_LINE */
    2107             :   }
    2108          28 :   return gerepileupto(av, sumnumrat_i(F, F0, vF, prec));
    2109             : }
    2110             : static GEN
    2111          28 : _add(GEN x, GEN y)
    2112             : {
    2113          28 :   if (!x) return y;
    2114          28 :   return y? gadd(x, y): x;
    2115             : }
    2116             : static GEN
    2117          28 : _sub(GEN x, GEN y)
    2118             : {
    2119          28 :   if (!x) return y? gneg(y): NULL;
    2120          14 :   return y? gsub(x, y): x;
    2121             : }
    2122             : static GEN
    2123     1844164 : _mpmul(GEN x, GEN y)
    2124             : {
    2125     1844164 :   if (!x) return y;
    2126     1839306 :   return y? mpmul(x, y): x;
    2127             : }
    2128             : 
    2129             : GEN
    2130          77 : sumaltrat(GEN F, GEN ga, long prec)
    2131             : {
    2132          77 :   pari_sp av = avma;
    2133          77 :   GEN G, res, vG = NULL, G0, F1, F2, X1, X2;
    2134             :   long a, vx;
    2135             : 
    2136          77 :   switch(typ(F))
    2137             :   {
    2138          63 :     case t_RFRAC: break;
    2139             :     case t_INT: case t_REAL: case t_COMPLEX: case t_POL:
    2140          14 :       if (gequal0(F)) return real_0(prec);
    2141           7 :     default: pari_err_TYPE("sumaltrat",F);
    2142             :   }
    2143          63 :   vx = varn(gel(F,2));
    2144          63 :   X2 = deg1pol_shallow(gen_2, gen_0, vx);
    2145          63 :   X1 = deg1pol_shallow(gen_2, gen_1, vx);
    2146          63 :   switch(typ(ga))
    2147             :   {
    2148             :     case t_INT:
    2149          28 :       a = itos(ga);
    2150          28 :       if (signe(ga)) F = gsubst(F, vx, deg1pol_shallow(gen_1,ga,vx));
    2151          28 :       G0 = NULL;
    2152          28 :       break;
    2153             :     case t_INFINITY:
    2154          35 :       if (inf_get_sign(ga) == -1)
    2155             :       {
    2156          28 :         GEN Fm = gsubst(F, vx, RgX_neg(pol_x(vx)));
    2157          28 :         G0 = _sub(rfrac_eval0(F), rfrac_eval(F, gen_1));
    2158          28 :         vG = mkvec2(F,Fm);
    2159          28 :         F = gadd(F, Fm);
    2160          28 :         if (gequal0(F)) { avma = av; return real_0(prec); }
    2161          21 :         a = 0;
    2162          21 :         break;
    2163             :       }
    2164             :     default:
    2165           7 :       pari_err_TYPE("sumnumaltrat",ga);
    2166             :       return NULL; /* LCOV_EXCL_LINE */
    2167             :   }
    2168          49 :   F2 = gsubst(F, vx, X2); /* F(2n) */
    2169          49 :   F1 = gneg(gsubst(F, vx, X1)); /* - F(2n+1) */
    2170          49 :   if (!G0) G0 = _add(rfrac_eval0(F2), rfrac_eval0(F1));
    2171          49 :   G = gadd(F2, F1);
    2172          49 :   if (gequal0(G)) { avma = av; return real_0(prec); }
    2173          42 :   if (vG)
    2174          14 :     vG = shallowconcat(gsubst(vG, vx, X2), gsubst(vG, vx, X1));
    2175             :   else
    2176          28 :     vG = mkvec2(F2, F1);
    2177          42 :   res = sumnumrat_i(G, G0, vG, prec);
    2178          42 :   if (odd(a)) res = gneg(res);
    2179          42 :   return gerepileupto(av, res);
    2180             : }
    2181             : 
    2182             : /* prod_{n >= a} F(n) */
    2183             : GEN
    2184          28 : prodnumrat(GEN F, long a, long prec)
    2185             : {
    2186          28 :   pari_sp ltop = avma;
    2187          28 :   long B = prec2nbits(prec), j, k, m, N, vx;
    2188          28 :   GEN S, S1, S2, intf, G, F1 = gsubgs(F,1);
    2189             :   double r;
    2190             : 
    2191          28 :   switch(typ(F1))
    2192             :   {
    2193          14 :     case t_RFRAC: break;
    2194             :     case t_INT: case t_REAL: case t_COMPLEX: case t_POL:
    2195          14 :       if (gequal0(F1)) return real_1(prec);
    2196           7 :     default: pari_err_TYPE("prodnumrat",F);
    2197             :   }
    2198          14 :   if (poldegree(F1,-1) > -2) pari_err(e_MISC, "product diverges in prodnumrat");
    2199           7 :   vx = varn(gel(F,2));
    2200           7 :   if (a) F = gsubst(F, vx, gaddgs(pol_x(vx), a));
    2201           7 :   r = gtodouble(ratpolemax2(F, NULL));
    2202           7 :   get_kN((long)ceil(r), B, &k,&N);
    2203           7 :   G = gdiv(deriv(F, vx), F);
    2204           7 :   intf = intnumainfrat(gmul(pol_x(vx),G), N, r, prec);
    2205           7 :   intf = gneg(gadd(intf, gmulsg(N, glog(gsubst(F, vx, stoi(N)), prec))));
    2206           7 :   S = gmul(real_1(prec), gsubst(G, vx, gaddgs(pol_x(vx), N)));
    2207           7 :   S = rfrac_to_ser(S, k + 2);
    2208           7 :   S1 = gsqrt(gsubst(F, vx, utoipos(N)), prec);
    2209           7 :   for (m = 0; m < N; m++) S1 = gmul(S1, gsubst(F, vx, utoi(m)));
    2210           7 :   S2 = gen_0;
    2211         336 :   for (j = 2; j <= k; j += 2)
    2212         329 :     S2 = gadd(S2, gmul(gdivgs(bernfrac(j),j*(j-1)), sercoeff(S, j-2)));
    2213           7 :   return gerepileupto(ltop, gmul(S1, gexp(gsub(intf, S2), prec)));
    2214             : }
    2215             : 
    2216             : static GEN
    2217        3374 : sdmob(GEN ser, long n)
    2218             : {
    2219        3374 :   GEN D = divisorsu(n), S = gen_0;
    2220        3374 :   long i, l = lg(D);
    2221       19481 :   for (i = 1; i < l; ++i)
    2222             :   {
    2223       16107 :     long d = D[i], mob = moebiusu(d);
    2224       16107 :     if (mob) S = gadd(S, gdivgs(sercoeff(ser, n/d), mob*d));
    2225             :   }
    2226        3374 :   return S;
    2227             : }
    2228             : static GEN
    2229        1498 : logzetan(GEN s, GEN P, long prec)
    2230             : {
    2231        1498 :   GEN negs = gneg(s), Z = gzeta(gprec_w(s,prec), prec);
    2232        1498 :   long i, l = lg(P);
    2233        1498 :   for (i = 1; i < l; i++) Z = gmul(Z, gsubsg(1, gpow(gel(P,i), negs, prec)));
    2234        1498 :   return glog(Z, prec);
    2235             : }
    2236             : static GEN
    2237          42 : sumlogzeta(GEN ser, GEN s, long N, long vF, long lim, long prec)
    2238             : {
    2239          42 :   GEN z = gen_0, P = primes_interval(gen_2, utoipos(N));
    2240             :   long n;
    2241        3416 :   for (n = lim; n >= vF; n--)
    2242             :   {
    2243        3374 :     GEN t = sdmob(ser, n);
    2244        3374 :     if (!gequal0(t))
    2245             :     {
    2246        1498 :       long e = gexpo(t), prec2 = e <= 0? prec: prec + nbits2prec(e);
    2247        1498 :       z = gadd(z, gmul(logzetan(gmulsg(n,s), P, prec2), t));
    2248        1498 :       z = gprec_w(z, prec);
    2249             :     }
    2250             :   }
    2251          42 :   return z;
    2252             : }
    2253             : 
    2254             : /* sum_{p prime, p >= a} F(p^s), F rational function */
    2255             : GEN
    2256          35 : sumeulerrat(GEN F, GEN s, long a, long prec)
    2257             : {
    2258          35 :   pari_sp av = avma;
    2259             :   forprime_t T;
    2260             :   GEN r, ser, res, rsg;
    2261             :   double rs, RS;
    2262          35 :   long B = prec2nbits(prec), vx, vF, p, N, lim;
    2263             : 
    2264          35 :   switch(typ(F))
    2265             :   {
    2266          21 :     case t_RFRAC: break;
    2267             :     case t_INT: case t_REAL: case t_COMPLEX: case t_POL:
    2268          14 :       if (gequal0(F)) return real_0(prec);
    2269           7 :     default: pari_err_TYPE("sumeulerrat",F);
    2270             :   }
    2271          21 :   if (!s) s = gen_1;
    2272          21 :   if (a < 2) a = 2;
    2273          21 :   vx = varn(gel(F,2));
    2274          21 :   vF = -poldegree(F, -1);
    2275          21 :   rsg = real_i(s);
    2276          21 :   rs = gtodouble(rsg);
    2277          21 :   r = ratpolemax(F);
    2278          21 :   RS = maxdd(1./vF, dbllog2(r) / log2((double)a));
    2279          21 :   if (rs <= RS)
    2280           7 :     pari_err_DOMAIN("sumeulerrat", "real(s)", "<=",  dbltor(RS), dbltor(rs));
    2281          14 :   N = maxss(maxss(30, a), (long)ceil(2*gtodouble(r)));
    2282          14 :   lim = (long)ceil(B / dbllog2(gdiv(gpow(stoi(N), rsg, LOWDEFAULTPREC), r)))+1;
    2283          14 :   ser = gmul(real_1(prec + EXTRAPREC), F);
    2284          14 :   ser = rfracrecip_to_ser_absolute(ser, lim);
    2285          14 :   res = sumlogzeta(ser, s, N, vF, lim, prec);
    2286          14 :   u_forprime_init(&T, a, N);
    2287         168 :   while ( (p = u_forprime_next(&T)) )
    2288         140 :     res = gadd(res, gsubst(F, vx, gpow(utoipos(p), s, prec)));
    2289          14 :   return gerepileupto(av, gprec_w(res, prec));
    2290             : }
    2291             : 
    2292             : /* prod_{p prime, p >= a} F(p^s), F rational function */
    2293             : GEN
    2294          49 : prodeulerrat(GEN F, GEN s, long a, long prec)
    2295             : {
    2296          49 :   pari_sp ltop = avma;
    2297             :   forprime_t T;
    2298             :   GEN F1, r, r1, ser, res, rsg;
    2299             :   double rs, RS;
    2300          49 :   long B = prec2nbits(prec), vx = gvar(F), vF, p, N, lim;
    2301             : 
    2302          49 :   F1 = gsubgs(F, 1);
    2303          49 :   switch(typ(F))
    2304             :   {
    2305          35 :     case t_RFRAC: break;
    2306             :     case t_INT: case t_REAL: case t_COMPLEX: case t_POL:
    2307          14 :       if (gequal0(F1)) return real_1(prec);
    2308           7 :     default: pari_err_TYPE("prodeulerrat",F);
    2309             :   }
    2310          35 :   if (!s) s = gen_1;
    2311          35 :   vF = -poldegree(F1, -1);
    2312          35 :   rsg = real_i(s);
    2313          35 :   rs = gtodouble(rsg);
    2314          35 :   r = ratpolemax2(F, &r1);
    2315          35 :   RS = maxdd(1./vF, dbllog2(r1) / log2((double)a));
    2316          35 :   if (rs <= RS)
    2317           7 :     pari_err_DOMAIN("prodeulerrat", "real(s)", "<=",  dbltor(RS), dbltor(rs));
    2318          28 :   N = maxss(maxss(30, a), (long)ceil(2*gtodouble(r)));
    2319          28 :   lim = (long)ceil(B / dbllog2(gdiv(gpow(stoi(N), rsg, LOWDEFAULTPREC), r)))+1;
    2320          28 :   ser = gmul(real_1(prec + EXTRAPREC), F1);
    2321          28 :   ser = gaddsg(1, rfracrecip_to_ser_absolute(ser, lim));
    2322          28 :   ser = glog(ser, 0);
    2323          28 :   res = sumlogzeta(ser, s, N, vF, lim, prec);
    2324          28 :   res = gexp(res, prec);
    2325          28 :   u_forprime_init(&T, a, N);
    2326         336 :   while ( (p = u_forprime_next(&T)) )
    2327         280 :     res = gmul(res, gsubst(F, vx, gpow(utoipos(p), s, prec)));
    2328          28 :   return gerepileupto(ltop, gprec_w(res, prec));
    2329             : }
    2330             : 
    2331             : /* Compute $\sum_{n\ge a}c(n)$ using Lagrange extrapolation.
    2332             : Assume that the $N$-th remainder term of the series has a
    2333             : regular asymptotic expansion in integral powers of $1/N$. */
    2334             : static GEN
    2335          35 : sumnumlagrange1init(GEN c1, long flag, long prec)
    2336             : {
    2337          35 :   pari_sp av = avma;
    2338             :   GEN V, W, T;
    2339             :   double c1d;
    2340          35 :   long B = prec2nbits(prec), prec2;
    2341             :   ulong n, N;
    2342          35 :   c1d = c1 ? gtodouble(c1) : 0.332;
    2343          35 :   N = (ulong)ceil(c1d*B); if ((N&1L) == 0) N++;
    2344          35 :   prec2 = nbits2prec(B+(long)ceil(1.8444*N) + 16);
    2345          35 :   W = vecbinome(N);
    2346          35 :   T = vecpowuu(N, N);
    2347          35 :   V = cgetg(N+1, t_VEC); gel(V,N) = gel(T,N);
    2348        3773 :   for (n = N-1; n >= 1; n--)
    2349             :   {
    2350        3738 :     pari_sp av = avma;
    2351        3738 :     GEN t = mulii(gel(W, n+1), gel(T,n));
    2352        3738 :     if (!odd(n)) togglesign_safe(&t);
    2353        3738 :     if (flag) t = addii(gel(V, n+1), t);
    2354        3738 :     gel(V, n) = gerepileuptoint(av, t);
    2355             :   }
    2356          35 :   V = gdiv(RgV_gtofp(V, prec2), mpfact(N));
    2357          35 :   return gerepilecopy(av, mkvec4(gen_1, stoi(prec2), gen_1, V));
    2358             : }
    2359             : 
    2360             : static GEN
    2361           7 : sumnumlagrange2init(GEN c1, long flag, long prec)
    2362             : {
    2363           7 :   pari_sp av = avma;
    2364             :   GEN V, W, T, told;
    2365           7 :   double c1d = c1 ? gtodouble(c1) : 0.228;
    2366           7 :   long B = prec2nbits(prec), prec2;
    2367             :   ulong n, N;
    2368             : 
    2369           7 :   N = (ulong)ceil(c1d*B); if ((N&1L) == 0) N++;
    2370           7 :   prec2 = nbits2prec(B+(long)ceil(1.18696*N) + 16);
    2371           7 :   W = vecbinome(2*N);
    2372           7 :   T = vecpowuu(N, 2*N);
    2373           7 :   V = cgetg(N+1, t_VEC); gel(V, N) = told = gel(T,N);
    2374         623 :   for (n = N-1; n >= 1; n--)
    2375             :   {
    2376         616 :     GEN tnew = mulii(gel(W, N-n+1), gel(T,n));
    2377         616 :     if (!odd(n)) togglesign_safe(&tnew);
    2378         616 :     told = addii(told, tnew);
    2379         616 :     if (flag) told = addii(gel(V, n+1), told);
    2380         616 :     gel(V, n) = told; told = tnew;
    2381             :   }
    2382           7 :   V = gdiv(RgV_gtofp(V, prec2), mpfact(2*N));
    2383           7 :   return gerepilecopy(av, mkvec4(gen_2, stoi(prec2), gen_1, V));
    2384             : }
    2385             : 
    2386             : /* Used only for al = 2, 1, 1/2, 1/3, 1/4. */
    2387             : static GEN
    2388          49 : sumnumlagrangeinit_i(GEN al, GEN c1, long flag, long prec)
    2389             : {
    2390          49 :   pari_sp av = avma;
    2391             :   GEN V, W;
    2392          49 :   double c1d = 0.0, c2;
    2393          49 :   long B = prec2nbits(prec), B1, prec2, dal;
    2394             :   ulong j, n, N;
    2395             : 
    2396          49 :   if (typ(al) == t_INT)
    2397             :   {
    2398          28 :     switch(itos_or_0(al))
    2399             :     {
    2400          21 :       case 1: return sumnumlagrange1init(c1, flag, prec);
    2401           7 :       case 2: return sumnumlagrange2init(c1, flag, prec);
    2402           0 :       default: pari_err_IMPL("sumnumlagrange for this alpha");
    2403             :     }
    2404             :   }
    2405          21 :   if (typ(al) != t_FRAC) pari_err_TYPE("sumnumlagrangeinit",al);
    2406          21 :   dal = itos_or_0(gel(al,2));
    2407          21 :   if (dal > 4 || !equali1(gel(al,1)))
    2408           7 :     pari_err_IMPL("sumnumlagrange for this alpha");
    2409          14 :   switch(dal)
    2410             :   {
    2411           7 :     case 2: c2 = 2.6441; c1d = 0.62; break;
    2412           7 :     case 3: c2 = 3.1578; c1d = 1.18; break;
    2413           0 :     case 4: c2 = 3.5364; c1d = 3.00; break;
    2414             :     default: return NULL; /* LCOV_EXCL_LINE */
    2415             :   }
    2416          14 :   if (c1)
    2417             :   {
    2418           0 :     c1d = gtodouble(c1);
    2419           0 :     if (c1d <= 0)
    2420           0 :       pari_err_DOMAIN("sumnumlagrangeinit", "c1", "<=", gen_0, c1);
    2421             :   }
    2422          14 :   N = (ulong)ceil(c1d*B); if ((N&1L) == 0) N++;
    2423          14 :   B1 = B + (long)ceil(c2*N) + 16;
    2424          14 :   prec2 = nbits2prec(B1);
    2425          14 :   V = vecpowug(N, al, prec2);
    2426          14 :   W = cgetg(N+1, t_VEC);
    2427        4872 :   for (n = 1; n <= N; ++n)
    2428             :   {
    2429        4858 :     pari_sp av2 = avma;
    2430        4858 :     GEN t = NULL, vn = gel(V, n);
    2431     1853880 :     for (j = 1; j <= N; j++)
    2432     1849022 :       if (j != n) t = _mpmul(t, mpsub(vn, gel(V, j)));
    2433        4858 :     gel(W, n) = gerepileuptoleaf(av2, mpdiv(gpowgs(vn, N-1), t));
    2434             :   }
    2435          14 :   if (flag)
    2436          14 :     for (n = N-1; n >= 1; n--) gel(W, n) = gadd(gel(W, n+1), gel(W, n));
    2437          14 :   return gerepilecopy(av, mkvec4(al, stoi(prec2), gen_1, W));
    2438             : }
    2439             : 
    2440             : GEN
    2441          63 : sumnumlagrangeinit(GEN al, GEN c1, long prec)
    2442             : {
    2443          63 :   pari_sp ltop = avma;
    2444             :   GEN V, W, S, be;
    2445             :   long n, prec2, fl, N;
    2446             : 
    2447          63 :   if (!al) return sumnumlagrange1init(c1, 1, prec);
    2448          49 :   if (typ(al) != t_VEC) al = mkvec2(gen_1, al);
    2449          35 :   else if (lg(al) != 3) pari_err_TYPE("sumnumlagrangeinit",al);
    2450          49 :   be = gel(al,2);
    2451          49 :   al = gel(al,1);
    2452          49 :   if (gequal0(be)) return sumnumlagrangeinit_i(al, c1, 1, prec);
    2453          14 :   V = sumnumlagrangeinit_i(al, c1, 0, prec);
    2454          14 :   switch(typ(be))
    2455             :   {
    2456           0 :     case t_CLOSURE: fl = 1; break;
    2457           7 :     case t_INT: case t_FRAC: case t_REAL: fl = 0; break;
    2458           7 :     default: pari_err_TYPE("sumnumlagrangeinit", be);
    2459             :              return NULL; /* LCOV_EXCL_LINE */
    2460             :   }
    2461           7 :   prec2 = itos(gel(V, 2));
    2462           7 :   W = gel(V, 4);
    2463           7 :   N = lg(W) - 1;
    2464           7 :   S = gen_0; V = cgetg(N+1, t_VEC);
    2465         910 :   for (n = N; n >= 1; n--)
    2466             :   {
    2467         903 :     GEN tmp, gn = stoi(n);
    2468         903 :     tmp = fl ? closure_callgen1prec(be, gn, prec2) : gpow(gn, gneg(be), prec2);
    2469         903 :     tmp = gdiv(gel(W, n), tmp);
    2470         903 :     S = gadd(S, tmp);
    2471         903 :     gel(V, n) = (n == N)? tmp: gadd(gel(V, n+1), tmp);
    2472             :   }
    2473           7 :   return gerepilecopy(ltop, mkvec4(al, stoi(prec2), S, V));
    2474             : }
    2475             : 
    2476             : /* - sum_{n=1}^{as-1} f(n) */
    2477             : static GEN
    2478          14 : sumaux(void *E, GEN (*eval)(void*,GEN,long), long as, long prec)
    2479             : {
    2480          14 :   GEN S = gen_0;
    2481             :   long n;
    2482          14 :   if (as > 1)
    2483             :   {
    2484           7 :     for (n = 1; n < as; ++n) S = gadd(S, eval(E, stoi(n), prec));
    2485           7 :     S = gneg(S);
    2486             :   }
    2487             :   else
    2488           7 :     for (n = as; n <= 0; ++n) S = gadd(S, eval(E, stoi(n), prec));
    2489          14 :   return S;
    2490             : }
    2491             : 
    2492             : GEN
    2493          84 : sumnumlagrange(void *E, GEN (*eval)(void*,GEN,long), GEN a, GEN tab, long prec)
    2494             : {
    2495          84 :   pari_sp av = avma;
    2496             :   GEN s, S, al, V;
    2497             :   long as, prec2;
    2498             :   ulong n, l;
    2499             : 
    2500          84 :   if (typ(a) != t_INT) pari_err_TYPE("sumnumlagrange", a);
    2501          84 :   if (!tab) tab = sumnumlagrangeinit(NULL, tab, prec);
    2502          70 :   else if (lg(tab) != 5 || typ(gel(tab,2)) != t_INT || typ(gel(tab,4)) != t_VEC)
    2503           0 :     pari_err_TYPE("sumnumlagrange", tab);
    2504             : 
    2505          84 :   as = itos(a);
    2506          84 :   al = gel(tab, 1);
    2507          84 :   prec2 = itos(gel(tab, 2));
    2508          84 :   S = gel(tab, 3);
    2509          84 :   V = gel(tab, 4);
    2510          84 :   l = lg(V);
    2511          84 :   if (gequal(al, gen_2))
    2512             :   {
    2513          14 :     s = sumaux(E, eval, as, prec2);
    2514          14 :     as = 1;
    2515             :   }
    2516             :   else
    2517          70 :     s = gen_0;
    2518       16464 :   for (n = 1; n < l; n++)
    2519       16380 :     s = gadd(s, gmul(gel(V, n), eval(E, stoi(n+as-1), prec2)));
    2520          84 :   if (!gequal1(S)) s = gdiv(s,S);
    2521          84 :   return gerepileupto(av, gprec_w(s, prec));
    2522             : }
    2523             : 
    2524             : GEN
    2525          84 : sumnumlagrange0(GEN a, GEN code, GEN tab, long prec)
    2526          84 : { EXPR_WRAP(code, sumnumlagrange(EXPR_ARGPREC, a, tab, prec)); }

Generated by: LCOV version 1.11