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 - sumiter.c (source / functions) Hit Total Coverage
Test: PARI/GP v2.10.0 lcov report (development 22307-7f6745a) Lines: 1036 1079 96.0 %
Date: 2018-04-22 06:16:17 Functions: 87 87 100.0 %
Legend: Lines: hit not hit

          Line data    Source code
       1             : /* Copyright (C) 2000  The PARI group.
       2             : 
       3             : This file is part of the PARI/GP package.
       4             : 
       5             : PARI/GP is free software; you can redistribute it and/or modify it under the
       6             : terms of the GNU General Public License as published by the Free Software
       7             : Foundation. It is distributed in the hope that it will be useful, but WITHOUT
       8             : ANY WARRANTY WHATSOEVER.
       9             : 
      10             : Check the License for details. You should have received a copy of it, along
      11             : with the package; see the file 'COPYING'. If not, write to the Free Software
      12             : Foundation, Inc., 51 Franklin Street, Fifth Floor, Boston, MA 02110-1301 USA. */
      13             : 
      14             : #include "pari.h"
      15             : #include "paripriv.h"
      16             : 
      17             : GEN
      18     1271739 : iferrpari(GEN a, GEN b, GEN c)
      19             : {
      20             :   GEN res;
      21             :   struct pari_evalstate state;
      22     1271739 :   evalstate_save(&state);
      23     1271739 :   pari_CATCH(CATCH_ALL)
      24             :   {
      25             :     GEN E;
      26       68862 :     if (!b&&!c) return gnil;
      27       34438 :     E = evalstate_restore_err(&state);
      28       34438 :     if (c)
      29             :     {
      30         292 :       push_lex(E,c);
      31         292 :       res = closure_evalnobrk(c);
      32         285 :       pop_lex(1);
      33         285 :       if (gequal0(res))
      34           7 :         pari_err(0, E);
      35             :     }
      36       34424 :     if (!b) return gnil;
      37       34424 :     push_lex(E,b);
      38       34424 :     res = closure_evalgen(b);
      39       34424 :     pop_lex(1);
      40       34424 :     return res;
      41             :   } pari_TRY {
      42     1271739 :     res = closure_evalgen(a);
      43     1237301 :   } pari_ENDCATCH;
      44     1237301 :   return res;
      45             : }
      46             : 
      47             : /********************************************************************/
      48             : /**                                                                **/
      49             : /**                        ITERATIONS                              **/
      50             : /**                                                                **/
      51             : /********************************************************************/
      52             : 
      53             : static void
      54     4926706 : forparii(GEN a, GEN b, GEN code)
      55             : {
      56     4926706 :   pari_sp av, av0 = avma;
      57             :   GEN aa;
      58     9853354 :   if (gcmp(b,a) < 0) return;
      59     4866472 :   if (typ(b) != t_INFINITY) b = gfloor(b);
      60     4866468 :   aa = a = setloop(a);
      61     4866473 :   av=avma;
      62     4866473 :   push_lex(a,code);
      63    92927852 :   while (gcmp(a,b) <= 0)
      64             :   {
      65    83632370 :     closure_evalvoid(code); if (loop_break()) break;
      66    83267184 :     a = get_lex(-1);
      67    83206459 :     if (a == aa)
      68             :     {
      69    83206431 :       a = incloop(a);
      70    83194882 :       if (a != aa) { set_lex(-1,a); aa = a; }
      71             :     }
      72             :     else
      73             :     { /* 'code' modified a ! Be careful (and slow) from now on */
      74          28 :       a = gaddgs(a,1);
      75          28 :       if (gc_needed(av,1))
      76             :       {
      77           0 :         if (DEBUGMEM>1) pari_warn(warnmem,"forparii");
      78           0 :         a = gerepileupto(av,a);
      79             :       }
      80          28 :       set_lex(-1,a);
      81             :     }
      82             :   }
      83     4806640 :   pop_lex(1);  avma = av0;
      84             : }
      85             : 
      86             : void
      87     4926714 : forpari(GEN a, GEN b, GEN code)
      88             : {
      89     4926714 :   pari_sp ltop=avma, av;
      90     9853368 :   if (typ(a) == t_INT) { forparii(a,b,code); return; }
      91           7 :   b = gcopy(b); /* Kludge to work-around the a+(a=2) bug */
      92           7 :   av=avma;
      93           7 :   push_lex(a,code);
      94          35 :   while (gcmp(a,b) <= 0)
      95             :   {
      96          21 :     closure_evalvoid(code); if (loop_break()) break;
      97          21 :     a = get_lex(-1); a = gaddgs(a,1);
      98          21 :     if (gc_needed(av,1))
      99             :     {
     100           0 :       if (DEBUGMEM>1) pari_warn(warnmem,"forpari");
     101           0 :       a = gerepileupto(av,a);
     102             :     }
     103          21 :     set_lex(-1, a);
     104             :   }
     105           7 :   pop_lex(1); avma = ltop;
     106             : }
     107             : 
     108             : /* 0 < a <= b. Using small consecutive chunks to 1) limit memory use, 2) allow
     109             :  * cheap early abort */
     110             : static int
     111          56 : forfactoredpos(ulong a, ulong b, GEN code)
     112             : {
     113          56 :   const ulong step = 1024;
     114          56 :   pari_sp av = avma;
     115             :   ulong x1;
     116        6881 :   for(x1 = a;; x1 += step, avma = av)
     117             :   { /* beware overflow, fuse last two bins (avoid a tiny remainder) */
     118        6881 :     ulong j, lv, x2 = (b >= 2*step && b - 2*step >= x1)? x1-1 + step: b;
     119        6881 :     GEN v = vecfactoru(x1, x2);
     120        6881 :     lv = lg(v);
     121     7008386 :     for (j = 1; j < lv; j++)
     122             :     {
     123     7001519 :       ulong n = x1-1 + j;
     124     7001519 :       set_lex(-1, mkvec2(utoipos(n), Flm_to_ZM(gel(v,j))));
     125     7001519 :       closure_evalvoid(code);
     126     7001519 :       if (loop_break()) return 1;
     127             :     }
     128        6867 :     if (x2 == b) break;
     129        6825 :     set_lex(-1, gen_0);
     130        6825 :   }
     131          42 :   return 0;
     132             : }
     133             : 
     134             : /* vector of primes to squarefree factorization */
     135             : static GEN
     136     4255531 : zv_to_ZM(GEN v)
     137     4255531 : { return mkmat2(zc_to_ZC(v), const_col(lg(v)-1,gen_1)); }
     138             : /* vector of primes to negative squarefree factorization */
     139             : static GEN
     140     4255531 : zv_to_mZM(GEN v)
     141             : {
     142     4255531 :   long i, l = lg(v);
     143     4255531 :   GEN w = cgetg(l+1, t_COL);
     144     4255531 :   gel(w,1) = gen_m1; for (i = 1; i < l; i++) gel(w,i+1) = utoipos(v[i]);
     145     4255531 :   return mkmat2(w, const_col(l,gen_1));
     146             : }
     147             : /* 0 <= a <= b. Using small consecutive chunks to 1) limit memory use, 2) allow
     148             :  * cheap early abort */
     149             : static void
     150          14 : forsquarefreepos(ulong a, ulong b, GEN code)
     151             : {
     152          14 :   const ulong step = 1024;
     153          14 :   pari_sp av = avma;
     154             :   ulong x1;
     155        6839 :   for(x1 = a;; x1 += step, avma = av)
     156             :   { /* beware overflow, fuse last two bins (avoid a tiny remainder) */
     157        6839 :     ulong j, lv, x2 = (b >= 2*step && b - 2*step >= x1)? x1-1 + step: b;
     158        6839 :     GEN v = vecfactorsquarefreeu(x1, x2);
     159        6839 :     lv = lg(v);
     160     7006916 :     for (j = 1; j < lv; j++) if (gel(v,j))
     161             :     {
     162     4255531 :       ulong n = x1-1 + j;
     163     4255531 :       set_lex(-1, mkvec2(utoipos(n), zv_to_ZM(gel(v,j))));
     164     4255545 :       closure_evalvoid(code); if (loop_break()) return;
     165             :     }
     166        6839 :     if (x2 == b) break;
     167        6825 :     set_lex(-1, gen_0);
     168        6825 :   }
     169             : }
     170             : /* 0 <= a <= b. Loop from -b, ... -a through squarefree integers */
     171             : static void
     172          14 : forsquarefreeneg(ulong a, ulong b, GEN code)
     173             : {
     174          14 :   const ulong step = 1024;
     175          14 :   pari_sp av = avma;
     176             :   ulong x2;
     177        6839 :   for(x2 = b;; x2 -= step, avma = av)
     178             :   { /* beware overflow, fuse last two bins (avoid a tiny remainder) */
     179        6839 :     ulong j, x1 = (x2 >= 2*step && x2-2*step >= a)? x2+1 - step: a;
     180        6839 :     GEN v = vecfactorsquarefreeu(x1, x2);
     181     7006916 :     for (j = lg(v)-1; j > 0; j--) if (gel(v,j))
     182             :     {
     183     4255531 :       ulong n = x1-1 + j;
     184     4255531 :       set_lex(-1, mkvec2(utoineg(n), zv_to_mZM(gel(v,j))));
     185     4255545 :       closure_evalvoid(code); if (loop_break()) return;
     186             :     }
     187        6839 :     if (x1 == a) break;
     188        6825 :     set_lex(-1, gen_0);
     189        6825 :   }
     190             : }
     191             : void
     192          28 : forsquarefree(GEN a, GEN b, GEN code)
     193             : {
     194          28 :   pari_sp av = avma;
     195             :   long s;
     196          28 :   if (typ(a) != t_INT) pari_err_TYPE("forsquarefree", a);
     197          28 :   if (typ(b) != t_INT) pari_err_TYPE("forsquarefree", b);
     198          56 :   if (cmpii(a,b) > 0) return;
     199          28 :   s = signe(a);
     200          28 :   if (s * signe(b) < 0) pari_err_TYPE("forsquarefree [!= signs]", mkvec2(a,b));
     201          28 :   push_lex(NULL,code);
     202          28 :   if (s < 0) forsquarefreeneg(itou(b), itou(a), code);
     203          14 :   else       forsquarefreepos(itou(a), itou(b), code);
     204          28 :   pop_lex(1); avma = av;
     205             : }
     206             : 
     207             : /* convert factoru(n) to factor(-n); M pre-allocated factorization matrix
     208             :  * with (-1)^1 already set */
     209             : static void
     210     7001582 : Flm2negfact(GEN v, GEN M)
     211             : {
     212     7001582 :   GEN p = gel(v,1), e = gel(v,2), P = gel(M,1), E = gel(M,2);
     213     7001582 :   long i, l = lg(p);
     214    26980058 :   for (i = 1; i < l; i++)
     215             :   {
     216    19978476 :     gel(P,i+1) = utoipos(p[i]);
     217    19978476 :     gel(E,i+1) = utoipos(e[i]);
     218             :   }
     219     7001582 :   setlg(P,l+1);
     220     7001582 :   setlg(E,l+1);
     221     7001582 : }
     222             : /* 0 < a <= b, from -b to -a */
     223             : static int
     224          84 : forfactoredneg(ulong a, ulong b, GEN code)
     225             : {
     226          84 :   const ulong step = 1024;
     227             :   GEN P, E, M;
     228             :   pari_sp av;
     229             :   ulong x2;
     230             : 
     231          84 :   P = cgetg(18, t_COL); gel(P,1) = gen_m1;
     232          84 :   E = cgetg(18, t_COL); gel(E,1) = gen_1;
     233          84 :   M = mkmat2(P,E);
     234          84 :   av = avma;
     235        6909 :   for(x2 = b;; x2 -= step, avma = av)
     236             :   { /* beware overflow, fuse last two bins (avoid a tiny remainder) */
     237        6909 :     ulong j, x1 = (x2 >= 2*step && x2-2*step >= a)? x2+1 - step: a;
     238        6909 :     GEN v = vecfactoru(x1, x2);
     239     7008470 :     for (j = lg(v)-1; j; j--)
     240             :     { /* run backward: from factor(x1..x2) to factor(-x2..-x1) */
     241     7001582 :       ulong n = x1-1 + j;
     242     7001582 :       Flm2negfact(gel(v,j), M);
     243     7001582 :       set_lex(-1, mkvec2(utoineg(n), M));
     244     7001582 :       closure_evalvoid(code); if (loop_break()) return 1;
     245             :     }
     246        6888 :     if (x1 == a) break;
     247        6825 :     set_lex(-1, gen_0);
     248        6825 :   }
     249          63 :   return 0;
     250             : }
     251             : static int
     252          70 : eval0(GEN code)
     253             : {
     254          70 :   pari_sp av = avma;
     255          70 :   set_lex(-1, mkvec2(gen_0, mkmat2(mkcol(gen_0),mkcol(gen_1))));
     256          70 :   closure_evalvoid(code); avma = av;
     257          70 :   return loop_break();
     258             : }
     259             : void
     260         133 : forfactored(GEN a, GEN b, GEN code)
     261             : {
     262         133 :   pari_sp av = avma;
     263         133 :   long sa, sb, stop = 0;
     264         133 :   if (typ(a) != t_INT) pari_err_TYPE("forfactored", a);
     265         133 :   if (typ(b) != t_INT) pari_err_TYPE("forfactored", b);
     266         266 :   if (cmpii(a,b) > 0) return;
     267         126 :   push_lex(NULL,code);
     268         126 :   sa = signe(a);
     269         126 :   sb = signe(b);
     270         126 :   if (sa < 0)
     271             :   {
     272          84 :     stop = forfactoredneg((sb < 0)? b[2]: 1UL, itou(a), code);
     273          84 :     if (!stop && sb >= 0) stop = eval0(code);
     274          84 :     if (!stop && sb > 0) forfactoredpos(1UL, b[2], code);
     275             :   }
     276             :   else
     277             :   {
     278          42 :     if (!sa) stop = eval0(code);
     279          42 :     if (!stop && sb) forfactoredpos(sa? a[2]: 1UL, itou(b), code);
     280             :   }
     281         126 :   pop_lex(1); avma = av;
     282             : }
     283             : void
     284     1793205 : whilepari(GEN a, GEN b)
     285             : {
     286     1793205 :   pari_sp av = avma;
     287             :   for(;;)
     288             :   {
     289    18690607 :     GEN res = closure_evalnobrk(a);
     290    18690607 :     if (gequal0(res)) break;
     291    16897402 :     avma = av;
     292    16897402 :     closure_evalvoid(b); if (loop_break()) break;
     293    16897402 :   }
     294     1793205 :   avma = av;
     295     1793205 : }
     296             : 
     297             : void
     298      222074 : untilpari(GEN a, GEN b)
     299             : {
     300      222074 :   pari_sp av = avma;
     301             :   for(;;)
     302             :   {
     303             :     GEN res;
     304     2194143 :     closure_evalvoid(b); if (loop_break()) break;
     305     2194143 :     res = closure_evalnobrk(a);
     306     2194143 :     if (!gequal0(res)) break;
     307     1972069 :     avma = av;
     308     1972069 :   }
     309      222074 :   avma = av;
     310      222074 : }
     311             : 
     312          28 : static int negcmp(GEN x, GEN y) { return gcmp(y,x); }
     313             : 
     314             : void
     315        1484 : forstep(GEN a, GEN b, GEN s, GEN code)
     316             : {
     317             :   long ss, i;
     318        1484 :   pari_sp av, av0 = avma;
     319        1484 :   GEN v = NULL;
     320             :   int (*cmp)(GEN,GEN);
     321             : 
     322        1484 :   b = gcopy(b);
     323        1484 :   s = gcopy(s); av = avma;
     324        1484 :   switch(typ(s))
     325             :   {
     326           7 :     case t_VEC: case t_COL: ss = gsigne(vecsum(s)); v = s; break;
     327           0 :     case t_INTMOD: a = gadd(a, gmod(gsub(gel(s,2),a), gel(s,1)));
     328           0 :                    s = gel(s,1);
     329        1477 :     default: ss = gsigne(s);
     330             :   }
     331        1484 :   if (!ss) pari_err_DOMAIN("forstep","step","=",gen_0,s);
     332        1477 :   cmp = (ss > 0)? &gcmp: &negcmp;
     333        1477 :   i = 0;
     334        1477 :   push_lex(a,code);
     335       50484 :   while (cmp(a,b) <= 0)
     336             :   {
     337       47530 :     closure_evalvoid(code); if (loop_break()) break;
     338       47530 :     if (v)
     339             :     {
     340          42 :       if (++i >= lg(v)) i = 1;
     341          42 :       s = gel(v,i);
     342             :     }
     343       47530 :     a = get_lex(-1); a = gadd(a,s);
     344             : 
     345       47530 :     if (gc_needed(av,1))
     346             :     {
     347           0 :       if (DEBUGMEM>1) pari_warn(warnmem,"forstep");
     348           0 :       a = gerepileupto(av,a);
     349             :     }
     350       47530 :     set_lex(-1,a);
     351             :   }
     352        1477 :   pop_lex(1); avma = av0;
     353        1477 : }
     354             : 
     355             : static void
     356          14 : _fordiv(GEN a, GEN code, GEN (*D)(GEN))
     357             : {
     358             :   long i, l;
     359          14 :   pari_sp av2, av = avma;
     360          14 :   GEN t = D(a);
     361          14 :   push_lex(gen_0,code); l=lg(t); av2 = avma;
     362         105 :   for (i=1; i<l; i++)
     363             :   {
     364          91 :     set_lex(-1,gel(t,i));
     365          91 :     closure_evalvoid(code); if (loop_break()) break;
     366          91 :     avma = av2;
     367             :   }
     368          14 :   pop_lex(1); avma=av;
     369          14 : }
     370             : void
     371           7 : fordiv(GEN a, GEN code) { return _fordiv(a, code, &divisors); }
     372             : void
     373           7 : fordivfactored(GEN a, GEN code) { return _fordiv(a, code, &divisors_factored); }
     374             : 
     375             : /* Embedded for loops:
     376             :  *   fl = 0: execute ch (a), where a = (ai) runs through all n-uplets in
     377             :  *     [m1,M1] x ... x [mn,Mn]
     378             :  *   fl = 1: impose a1 <= ... <= an
     379             :  *   fl = 2:        a1 <  ... <  an
     380             :  */
     381             : /* increment and return d->a [over integers]*/
     382             : static GEN
     383      181472 : _next_i(forvec_t *d)
     384             : {
     385      181472 :   long i = d->n;
     386      181472 :   if (d->first) { d->first = 0; return (GEN)d->a; }
     387             :   for (;;) {
     388      233116 :     if (cmpii(d->a[i], d->M[i]) < 0) {
     389      181236 :       d->a[i] = incloop(d->a[i]);
     390      181236 :       return (GEN)d->a;
     391             :     }
     392       51880 :     d->a[i] = resetloop(d->a[i], d->m[i]);
     393       51880 :     if (--i <= 0) return NULL;
     394       51762 :   }
     395             : }
     396             : /* increment and return d->a [generic]*/
     397             : static GEN
     398          63 : _next(forvec_t *d)
     399             : {
     400          63 :   long i = d->n;
     401          63 :   if (d->first) { d->first = 0; return (GEN)d->a; }
     402             :   for (;;) {
     403          98 :     d->a[i] = gaddgs(d->a[i], 1);
     404          98 :     if (gcmp(d->a[i], d->M[i]) <= 0) return (GEN)d->a;
     405          49 :     d->a[i] = d->m[i];
     406          49 :     if (--i <= 0) return NULL;
     407          42 :   }
     408             : }
     409             : 
     410             : /* non-decreasing order [over integers] */
     411             : static GEN
     412         157 : _next_le_i(forvec_t *d)
     413             : {
     414         157 :   long i = d->n;
     415         157 :   if (d->first) { d->first = 0; return (GEN)d->a; }
     416             :   for (;;) {
     417         231 :     if (cmpii(d->a[i], d->M[i]) < 0)
     418             :     {
     419         117 :       d->a[i] = incloop(d->a[i]);
     420             :       /* m[i] < a[i] <= M[i] <= M[i+1] */
     421         301 :       while (i < d->n)
     422             :       {
     423             :         GEN t;
     424          67 :         i++;
     425          67 :         if (cmpii(d->a[i-1], d->a[i]) <= 0) continue;
     426             :         /* a[i] < a[i-1] <= M[i-1] <= M[i] */
     427          67 :         t = d->a[i-1]; if (cmpii(t, d->m[i]) < 0) t = d->m[i];
     428          67 :         d->a[i] = resetloop(d->a[i], t);/*a[i]:=max(a[i-1],m[i])*/
     429             :       }
     430         117 :       return (GEN)d->a;
     431             :     }
     432         114 :     d->a[i] = resetloop(d->a[i], d->m[i]);
     433         114 :     if (--i <= 0) return NULL;
     434          94 :   }
     435             : }
     436             : /* non-decreasing order [generic] */
     437             : static GEN
     438         154 : _next_le(forvec_t *d)
     439             : {
     440         154 :   long i = d->n;
     441         154 :   if (d->first) { d->first = 0; return (GEN)d->a; }
     442             :   for (;;) {
     443         266 :     d->a[i] = gaddgs(d->a[i], 1);
     444         266 :     if (gcmp(d->a[i], d->M[i]) <= 0)
     445             :     {
     446         350 :       while (i < d->n)
     447             :       {
     448             :         GEN c;
     449          98 :         i++;
     450          98 :         if (gcmp(d->a[i-1], d->a[i]) <= 0) continue;
     451             :         /* M[i] >= M[i-1] >= a[i-1] > a[i] */
     452          98 :         c = gceil(gsub(d->a[i-1], d->a[i]));
     453          98 :         d->a[i] = gadd(d->a[i], c);
     454             :         /* a[i-1] <= a[i] < M[i-1] + 1 => a[i] < M[i]+1 => a[i] <= M[i] */
     455             :       }
     456         126 :       return (GEN)d->a;
     457             :     }
     458         140 :     d->a[i] = d->m[i];
     459         140 :     if (--i <= 0) return NULL;
     460         126 :   }
     461             : }
     462             : /* strictly increasing order [over integers] */
     463             : static GEN
     464     1173546 : _next_lt_i(forvec_t *d)
     465             : {
     466     1173546 :   long i = d->n;
     467     1173546 :   if (d->first) { d->first = 0; return (GEN)d->a; }
     468             :   for (;;) {
     469     1290065 :     if (cmpii(d->a[i], d->M[i]) < 0)
     470             :     {
     471     1159940 :       d->a[i] = incloop(d->a[i]);
     472             :       /* m[i] < a[i] <= M[i] < M[i+1] */
     473     2436385 :       while (i < d->n)
     474             :       {
     475             :         pari_sp av;
     476             :         GEN t;
     477      116505 :         i++;
     478      116505 :         if (cmpii(d->a[i-1], d->a[i]) < 0) continue;
     479      116505 :         av = avma;
     480             :         /* M[i] > M[i-1] >= a[i-1] */
     481      116505 :         t = addiu(d->a[i-1],1); if (cmpii(t, d->m[i]) < 0) t = d->m[i];
     482      116505 :         d->a[i] = resetloop(d->a[i], t);/*a[i]:=max(a[i-1]+1,m[i]) <= M[i]*/
     483      116505 :         avma = av;
     484             :       }
     485     1159940 :       return (GEN)d->a;
     486             :     }
     487      130125 :     d->a[i] = resetloop(d->a[i], d->m[i]);
     488      130125 :     if (--i <= 0) return NULL;
     489      123322 :   }
     490             : }
     491             : /* strictly increasing order [generic] */
     492             : static GEN
     493          84 : _next_lt(forvec_t *d)
     494             : {
     495          84 :   long i = d->n;
     496          84 :   if (d->first) { d->first = 0; return (GEN)d->a; }
     497             :   for (;;) {
     498         133 :     d->a[i] = gaddgs(d->a[i], 1);
     499         133 :     if (gcmp(d->a[i], d->M[i]) <= 0)
     500             :     {
     501         147 :       while (i < d->n)
     502             :       {
     503             :         GEN c;
     504          35 :         i++;
     505          35 :         if (gcmp(d->a[i-1], d->a[i]) < 0) continue;
     506             :         /* M[i] > M[i-1] >= a[i-1] >= a[i] */
     507          35 :         c = addiu(gfloor(gsub(d->a[i-1], d->a[i])), 1); /* > a[i-1] - a[i] */
     508          35 :         d->a[i] = gadd(d->a[i], c);
     509             :         /* a[i-1] < a[i] <= M[i-1] + 1 => a[i] < M[i]+1 => a[i] <= M[i] */
     510             :       }
     511          56 :       return (GEN)d->a;
     512             :     }
     513          77 :     d->a[i] = d->m[i];
     514          77 :     if (--i <= 0) return NULL;
     515          63 :   }
     516             : }
     517             : /* for forvec(v=[],) */
     518             : static GEN
     519          14 : _next_void(forvec_t *d)
     520             : {
     521          14 :   if (d->first) { d->first = 0; return (GEN)d->a; }
     522           7 :   return NULL;
     523             : }
     524             : 
     525             : /* Initialize minima (m) and maxima (M); guarantee M[i] - m[i] integer and
     526             :  *   if flag = 1: m[i-1] <= m[i] <= M[i] <= M[i+1]
     527             :  *   if flag = 2: m[i-1] <  m[i] <= M[i] <  M[i+1],
     528             :  * for all i */
     529             : int
     530        7018 : forvec_init(forvec_t *d, GEN x, long flag)
     531             : {
     532        7018 :   long i, tx = typ(x), l = lg(x), t = t_INT;
     533        7018 :   if (!is_vec_t(tx)) pari_err_TYPE("forvec [not a vector]", x);
     534        7018 :   d->first = 1;
     535        7018 :   d->n = l - 1;
     536        7018 :   d->a = (GEN*)cgetg(l,tx);
     537        7018 :   d->m = (GEN*)cgetg(l,tx);
     538        7018 :   d->M = (GEN*)cgetg(l,tx);
     539        7018 :   if (l == 1) { d->next = &_next_void; return 1; }
     540       21285 :   for (i = 1; i < l; i++)
     541             :   {
     542       14302 :     GEN a, e = gel(x,i), m = gel(e,1), M = gel(e,2);
     543       14302 :     tx = typ(e);
     544       14302 :     if (! is_vec_t(tx) || lg(e)!=3)
     545          14 :       pari_err_TYPE("forvec [expected vector not of type [min,MAX]]",e);
     546       14288 :     if (typ(m) != t_INT) t = t_REAL;
     547       14288 :     if (i > 1) switch(flag)
     548             :     {
     549             :       case 1: /* a >= m[i-1] - m */
     550          55 :         a = gceil(gsub(d->m[i-1], m));
     551          55 :         if (typ(a) != t_INT) pari_err_TYPE("forvec",a);
     552          55 :         if (signe(a) > 0) m = gadd(m, a); else m = gcopy(m);
     553          55 :         break;
     554             :       case 2: /* a > m[i-1] - m */
     555        6852 :         a = gfloor(gsub(d->m[i-1], m));
     556        6852 :         if (typ(a) != t_INT) pari_err_TYPE("forvec",a);
     557        6852 :         a = addiu(a, 1);
     558        6852 :         if (signe(a) > 0) m = gadd(m, a); else m = gcopy(m);
     559        6852 :         break;
     560         384 :       default: m = gcopy(m);
     561         384 :         break;
     562             :     }
     563       14288 :     M = gadd(m, gfloor(gsub(M,m))); /* ensure M-m is an integer */
     564       14281 :     if (gcmp(m,M) > 0) { d->a = NULL; d->next = &_next; return 0; }
     565       14274 :     d->m[i] = m;
     566       14274 :     d->M[i] = M;
     567             :   }
     568        7038 :   if (flag == 1) for (i = l-2; i >= 1; i--)
     569             :   {
     570          55 :     GEN M = d->M[i], a = gfloor(gsub(d->M[i+1], M));
     571          55 :     if (typ(a) != t_INT) pari_err_TYPE("forvec",a);
     572             :     /* M[i]+a <= M[i+1] */
     573          55 :     if (signe(a) < 0) d->M[i] = gadd(M, a);
     574             :   }
     575       13794 :   else if (flag == 2) for (i = l-2; i >= 1; i--)
     576             :   {
     577        6845 :     GEN M = d->M[i], a = gceil(gsub(d->M[i+1], M));
     578        6845 :     if (typ(a) != t_INT) pari_err_TYPE("forvec",a);
     579        6845 :     a = subiu(a, 1);
     580             :     /* M[i]+a < M[i+1] */
     581        6845 :     if (signe(a) < 0) d->M[i] = gadd(M, a);
     582             :   }
     583        6983 :   if (t == t_INT) {
     584       21110 :     for (i = 1; i < l; i++) {
     585       14162 :       d->a[i] = setloop(d->m[i]);
     586       14162 :       if (typ(d->M[i]) != t_INT) d->M[i] = gfloor(d->M[i]);
     587             :     }
     588             :   } else {
     589          35 :     for (i = 1; i < l; i++) d->a[i] = d->m[i];
     590             :   }
     591        6983 :   switch(flag)
     592             :   {
     593         125 :     case 0: d->next = t==t_INT? &_next_i:    &_next; break;
     594          34 :     case 1: d->next = t==t_INT? &_next_le_i: &_next_le; break;
     595        6817 :     case 2: d->next = t==t_INT? &_next_lt_i: &_next_lt; break;
     596           7 :     default: pari_err_FLAG("forvec");
     597             :   }
     598        6976 :   return 1;
     599             : }
     600             : GEN
     601     1355490 : forvec_next(forvec_t *d) { return d->next(d); }
     602             : 
     603             : void
     604        7000 : forvec(GEN x, GEN code, long flag)
     605             : {
     606        7000 :   pari_sp av = avma;
     607             :   forvec_t T;
     608             :   GEN v;
     609       13972 :   if (!forvec_init(&T, x, flag)) { avma = av; return; }
     610        6965 :   push_lex((GEN)T.a, code);
     611        6965 :   while ((v = forvec_next(&T)))
     612             :   {
     613     1348333 :     closure_evalvoid(code);
     614     1348333 :     if (loop_break()) break;
     615             :   }
     616        6965 :   pop_lex(1); avma = av;
     617             : }
     618             : 
     619             : /********************************************************************/
     620             : /**                                                                **/
     621             : /**                              SUMS                              **/
     622             : /**                                                                **/
     623             : /********************************************************************/
     624             : 
     625             : GEN
     626       63406 : somme(GEN a, GEN b, GEN code, GEN x)
     627             : {
     628       63406 :   pari_sp av, av0 = avma;
     629             :   GEN p1;
     630             : 
     631       63406 :   if (typ(a) != t_INT) pari_err_TYPE("sum",a);
     632       63406 :   if (!x) x = gen_0;
     633       63406 :   if (gcmp(b,a) < 0) return gcopy(x);
     634             : 
     635       63406 :   b = gfloor(b);
     636       63406 :   a = setloop(a);
     637       63406 :   av=avma;
     638       63406 :   push_lex(a,code);
     639             :   for(;;)
     640             :   {
     641     1716421 :     p1 = closure_evalnobrk(code);
     642     1716421 :     x=gadd(x,p1); if (cmpii(a,b) >= 0) break;
     643     1653015 :     a = incloop(a);
     644     1653015 :     if (gc_needed(av,1))
     645             :     {
     646           0 :       if (DEBUGMEM>1) pari_warn(warnmem,"sum");
     647           0 :       x = gerepileupto(av,x);
     648             :     }
     649     1653015 :     set_lex(-1,a);
     650     1653015 :   }
     651       63406 :   pop_lex(1); return gerepileupto(av0,x);
     652             : }
     653             : 
     654             : static GEN
     655          21 : sum_init(GEN x0, GEN t)
     656             : {
     657          21 :   long tp = typ(t);
     658             :   GEN x;
     659          21 :   if (is_vec_t(tp))
     660             :   {
     661           7 :     x = const_vec(lg(t)-1, x0);
     662           7 :     settyp(x, tp);
     663             :   }
     664             :   else
     665          14 :     x = x0;
     666          21 :   return x;
     667             : }
     668             : 
     669             : GEN
     670          21 : suminf(void *E, GEN (*eval)(void *, GEN), GEN a, long prec)
     671             : {
     672          21 :   long fl = 0, G = prec2nbits(prec) + 5;
     673          21 :   pari_sp av0 = avma, av;
     674          21 :   GEN x = NULL, _1;
     675             : 
     676          21 :   if (typ(a) != t_INT) pari_err_TYPE("suminf",a);
     677          21 :   a = setloop(a);
     678          21 :   av = avma;
     679             :   for(;;)
     680             :   {
     681       15379 :     GEN t = eval(E, a);
     682       15379 :     if (!x) _1 = x = sum_init(real_1(prec), t);
     683             : 
     684       15379 :     x = gadd(x,t);
     685       15379 :     if (!gequal0(t) && gexpo(t) > gexpo(x)-G)
     686       15316 :       fl = 0;
     687          63 :     else if (++fl == 3)
     688          21 :       break;
     689       15358 :     a = incloop(a);
     690       15358 :     if (gc_needed(av,1))
     691             :     {
     692           0 :       if (DEBUGMEM>1) pari_warn(warnmem,"suminf");
     693           0 :       gerepileall(av,2, &x, &_1);
     694             :     }
     695       15358 :   }
     696          21 :   return gerepileupto(av0, gsub(x, _1));
     697             : }
     698             : GEN
     699          21 : suminf0(GEN a, GEN code, long prec)
     700          21 : { EXPR_WRAP(code, suminf(EXPR_ARG, a, prec)); }
     701             : 
     702             : GEN
     703          49 : sumdivexpr(GEN num, GEN code)
     704             : {
     705          49 :   pari_sp av = avma;
     706          49 :   GEN y = gen_0, t = divisors(num);
     707          49 :   long i, l = lg(t);
     708             : 
     709          49 :   push_lex(gen_0, code);
     710        9289 :   for (i=1; i<l; i++)
     711             :   {
     712        9240 :     set_lex(-1,gel(t,i));
     713        9240 :     y = gadd(y, closure_evalnobrk(code));
     714             :   }
     715          49 :   pop_lex(1); return gerepileupto(av,y);
     716             : }
     717             : GEN
     718          42 : sumdivmultexpr(GEN num, GEN code)
     719             : {
     720          42 :   pari_sp av = avma;
     721          42 :   GEN y = gen_1, P,E;
     722          42 :   int isint = divisors_init(num, &P,&E);
     723          42 :   long i, l = lg(P);
     724             :   GEN (*mul)(GEN,GEN);
     725             : 
     726          42 :   if (l == 1) { avma = av; return gen_1; }
     727          42 :   push_lex(gen_0, code);
     728          42 :   mul = isint? mulii: gmul;
     729         196 :   for (i=1; i<l; i++)
     730             :   {
     731         154 :     GEN p = gel(P,i), q = p, z = gen_1;
     732         154 :     long j, e = E[i];
     733         560 :     for (j = 1; j <= e; j++, q = mul(q, p))
     734             :     {
     735         560 :       set_lex(-1, q);
     736         560 :       z = gadd(z, closure_evalnobrk(code));
     737         560 :       if (j == e) break;
     738             :     }
     739         154 :     y = gmul(y, z);
     740             :   }
     741          42 :   pop_lex(1); return gerepileupto(av,y);
     742             : }
     743             : 
     744             : /********************************************************************/
     745             : /**                                                                **/
     746             : /**                           PRODUCTS                             **/
     747             : /**                                                                **/
     748             : /********************************************************************/
     749             : 
     750             : GEN
     751      120134 : produit(GEN a, GEN b, GEN code, GEN x)
     752             : {
     753      120134 :   pari_sp av, av0 = avma;
     754             :   GEN p1;
     755             : 
     756      120134 :   if (typ(a) != t_INT) pari_err_TYPE("prod",a);
     757      120134 :   if (!x) x = gen_1;
     758      120134 :   if (gcmp(b,a) < 0) return gcopy(x);
     759             : 
     760      115206 :   b = gfloor(b);
     761      115206 :   a = setloop(a);
     762      115206 :   av=avma;
     763      115206 :   push_lex(a,code);
     764             :   for(;;)
     765             :   {
     766      348768 :     p1 = closure_evalnobrk(code);
     767      348768 :     x = gmul(x,p1); if (cmpii(a,b) >= 0) break;
     768      233562 :     a = incloop(a);
     769      233562 :     if (gc_needed(av,1))
     770             :     {
     771           0 :       if (DEBUGMEM>1) pari_warn(warnmem,"prod");
     772           0 :       x = gerepileupto(av,x);
     773             :     }
     774      233562 :     set_lex(-1,a);
     775      233562 :   }
     776      115206 :   pop_lex(1); return gerepileupto(av0,x);
     777             : }
     778             : 
     779             : GEN
     780           7 : prodinf(void *E, GEN (*eval)(void *, GEN), GEN a, long prec)
     781             : {
     782           7 :   pari_sp av0 = avma, av;
     783             :   long fl,G;
     784           7 :   GEN p1,x = real_1(prec);
     785             : 
     786           7 :   if (typ(a) != t_INT) pari_err_TYPE("prodinf",a);
     787           7 :   a = setloop(a);
     788           7 :   av = avma;
     789           7 :   fl=0; G = -prec2nbits(prec)-5;
     790             :   for(;;)
     791             :   {
     792         952 :     p1 = eval(E, a); if (gequal0(p1)) { x = p1; break; }
     793         952 :     x = gmul(x,p1); a = incloop(a);
     794         952 :     p1 = gsubgs(p1, 1);
     795         952 :     if (gequal0(p1) || gexpo(p1) <= G) { if (++fl==3) break; } else fl=0;
     796         945 :     if (gc_needed(av,1))
     797             :     {
     798           0 :       if (DEBUGMEM>1) pari_warn(warnmem,"prodinf");
     799           0 :       x = gerepileupto(av,x);
     800             :     }
     801         945 :   }
     802           7 :   return gerepilecopy(av0,x);
     803             : }
     804             : GEN
     805           7 : prodinf1(void *E, GEN (*eval)(void *, GEN), GEN a, long prec)
     806             : {
     807           7 :   pari_sp av0 = avma, av;
     808             :   long fl,G;
     809           7 :   GEN p1,p2,x = real_1(prec);
     810             : 
     811           7 :   if (typ(a) != t_INT) pari_err_TYPE("prodinf1",a);
     812           7 :   a = setloop(a);
     813           7 :   av = avma;
     814           7 :   fl=0; G = -prec2nbits(prec)-5;
     815             :   for(;;)
     816             :   {
     817         952 :     p2 = eval(E, a); p1 = gaddgs(p2,1);
     818         952 :     if (gequal0(p1)) { x = p1; break; }
     819         952 :     x = gmul(x,p1); a = incloop(a);
     820         952 :     if (gequal0(p2) || gexpo(p2) <= G) { if (++fl==3) break; } else fl=0;
     821         945 :     if (gc_needed(av,1))
     822             :     {
     823           0 :       if (DEBUGMEM>1) pari_warn(warnmem,"prodinf1");
     824           0 :       x = gerepileupto(av,x);
     825             :     }
     826         945 :   }
     827           7 :   return gerepilecopy(av0,x);
     828             : }
     829             : GEN
     830          14 : prodinf0(GEN a, GEN code, long flag, long prec)
     831             : {
     832          14 :   switch(flag)
     833             :   {
     834           7 :     case 0: EXPR_WRAP(code, prodinf (EXPR_ARG, a, prec));
     835           7 :     case 1: EXPR_WRAP(code, prodinf1(EXPR_ARG, a, prec));
     836             :   }
     837           0 :   pari_err_FLAG("prodinf");
     838             :   return NULL; /* LCOV_EXCL_LINE */
     839             : }
     840             : 
     841             : GEN
     842           7 : prodeuler(void *E, GEN (*eval)(void *, GEN), GEN a, GEN b, long prec)
     843             : {
     844           7 :   pari_sp av, av0 = avma;
     845           7 :   GEN x = real_1(prec), prime;
     846             :   forprime_t T;
     847             : 
     848           7 :   av = avma;
     849           7 :   if (!forprime_init(&T, a,b)) { avma = av; return x; }
     850             : 
     851           7 :   av = avma;
     852        8617 :   while ( (prime = forprime_next(&T)) )
     853             :   {
     854        8603 :     x = gmul(x, eval(E, prime));
     855        8603 :     if (gc_needed(av,1))
     856             :     {
     857           0 :       if (DEBUGMEM>1) pari_warn(warnmem,"prodeuler");
     858           0 :       x = gerepilecopy(av, x);
     859             :     }
     860             :   }
     861           7 :   return gerepilecopy(av0,x);
     862             : }
     863             : GEN
     864           7 : prodeuler0(GEN a, GEN b, GEN code, long prec)
     865           7 : { EXPR_WRAP(code, prodeuler(EXPR_ARG, a, b, prec)); }
     866             : GEN
     867         126 : direuler0(GEN a, GEN b, GEN code, GEN c)
     868         126 : { EXPR_WRAP(code, direuler(EXPR_ARG, a, b, c)); }
     869             : 
     870             : /********************************************************************/
     871             : /**                                                                **/
     872             : /**                       VECTORS & MATRICES                       **/
     873             : /**                                                                **/
     874             : /********************************************************************/
     875             : 
     876             : INLINE GEN
     877     3461283 : copyupto(GEN z, GEN t)
     878             : {
     879     3461283 :   if (is_universal_constant(z) || (z>(GEN)pari_mainstack->bot && z<=t))
     880     3461276 :     return z;
     881             :   else
     882           7 :     return gcopy(z);
     883             : }
     884             : 
     885             : GEN
     886      114996 : vecexpr0(GEN vec, GEN code, GEN pred)
     887             : {
     888      114996 :   switch(typ(vec))
     889             :   {
     890             :     case t_LIST:
     891             :     {
     892          21 :       if (list_typ(vec)==t_LIST_MAP)
     893           7 :         vec = mapdomain_shallow(vec);
     894             :       else
     895          14 :         vec = list_data(vec);
     896          21 :       if (!vec) return cgetg(1, t_VEC);
     897          14 :       break;
     898             :     }
     899      114975 :     case t_VEC: case t_COL: case t_MAT: break;
     900           0 :     default: pari_err_TYPE("[_|_<-_,_]",vec);
     901             :   }
     902      114989 :   if (pred && code)
     903         273 :     EXPR_WRAP(code,vecselapply((void*)pred,&gp_evalbool,EXPR_ARGUPTO,vec))
     904      114716 :   else if (code)
     905      114716 :     EXPR_WRAP(code,vecapply(EXPR_ARGUPTO,vec))
     906             :   else
     907           0 :     EXPR_WRAP(pred,vecselect(EXPR_ARGBOOL,vec))
     908             : }
     909             : 
     910             : GEN
     911         175 : vecexpr1(GEN vec, GEN code, GEN pred)
     912             : {
     913         175 :   GEN v = vecexpr0(vec, code, pred);
     914         175 :   return lg(v) == 1? v: shallowconcat1(v);
     915             : }
     916             : 
     917             : GEN
     918     2354433 : vecteur(GEN nmax, GEN code)
     919             : {
     920             :   GEN y, c;
     921     2354433 :   long i, m = gtos(nmax);
     922             : 
     923     2354433 :   if (m < 0)  pari_err_DOMAIN("vector", "dimension", "<", gen_0, stoi(m));
     924     2354419 :   if (!code) return zerovec(m);
     925        8418 :   c = cgetipos(3); /* left on stack */
     926        8418 :   y = cgetg(m+1,t_VEC); push_lex(c, code);
     927     1682293 :   for (i=1; i<=m; i++)
     928             :   {
     929     1673889 :     c[2] = i;
     930     1673889 :     gel(y,i) = copyupto(closure_evalnobrk(code), y);
     931     1673875 :     set_lex(-1,c);
     932             :   }
     933        8404 :   pop_lex(1); return y;
     934             : }
     935             : 
     936             : GEN
     937         763 : vecteursmall(GEN nmax, GEN code)
     938             : {
     939             :   pari_sp av;
     940             :   GEN y, c;
     941         763 :   long i, m = gtos(nmax);
     942             : 
     943         763 :   if (m < 0)  pari_err_DOMAIN("vectorsmall", "dimension", "<", gen_0, stoi(m));
     944         756 :   if (!code) return zero_zv(m);
     945         735 :   c = cgetipos(3); /* left on stack */
     946         735 :   y = cgetg(m+1,t_VECSMALL); push_lex(c,code);
     947         735 :   av = avma;
     948       15974 :   for (i = 1; i <= m; i++)
     949             :   {
     950       15246 :     c[2] = i;
     951       15246 :     y[i] = gtos(closure_evalnobrk(code));
     952       15239 :     avma = av;
     953       15239 :     set_lex(-1,c);
     954             :   }
     955         728 :   pop_lex(1); return y;
     956             : }
     957             : 
     958             : GEN
     959         490 : vvecteur(GEN nmax, GEN n)
     960             : {
     961         490 :   GEN y = vecteur(nmax,n);
     962         483 :   settyp(y,t_COL); return y;
     963             : }
     964             : 
     965             : GEN
     966      138215 : matrice(GEN nlig, GEN ncol, GEN code)
     967             : {
     968             :   GEN c1, c2, y;
     969             :   long i, m, n;
     970             : 
     971      138215 :   n = gtos(nlig);
     972      138215 :   m = ncol? gtos(ncol): n;
     973      138215 :   if (m < 0)  pari_err_DOMAIN("matrix", "nbcols", "<", gen_0, stoi(m));
     974      138208 :   if (n < 0)  pari_err_DOMAIN("matrix", "nbrows", "<", gen_0, stoi(n));
     975      138201 :   if (!m) return cgetg(1,t_MAT);
     976      138159 :   if (!code || !n) return zeromatcopy(n, m);
     977      136003 :   c1 = cgetipos(3); push_lex(c1,code);
     978      136003 :   c2 = cgetipos(3); push_lex(c2,NULL); /* c1,c2 left on stack */
     979      136003 :   y = cgetg(m+1,t_MAT);
     980      526400 :   for (i = 1; i <= m; i++)
     981             :   {
     982      390397 :     GEN z = cgetg(n+1,t_COL);
     983             :     long j;
     984      390397 :     c2[2] = i; gel(y,i) = z;
     985     2177798 :     for (j = 1; j <= n; j++)
     986             :     {
     987     1787401 :       c1[2] = j;
     988     1787401 :       gel(z,j) = copyupto(closure_evalnobrk(code), y);
     989     1787401 :       set_lex(-2,c1);
     990     1787401 :       set_lex(-1,c2);
     991             :     }
     992             :   }
     993      136003 :   pop_lex(2); return y;
     994             : }
     995             : 
     996             : /********************************************************************/
     997             : /**                                                                **/
     998             : /**                         SUMMING SERIES                         **/
     999             : /**                                                                **/
    1000             : /********************************************************************/
    1001             : /* h = (2+2x)g'- g; g has t_INT coeffs */
    1002             : static GEN
    1003        1246 : delt(GEN g, long n)
    1004             : {
    1005        1246 :   GEN h = cgetg(n+3,t_POL);
    1006             :   long k;
    1007        1246 :   h[1] = g[1];
    1008        1246 :   gel(h,2) = gel(g,2);
    1009      359597 :   for (k=1; k<n; k++)
    1010      358351 :     gel(h,k+2) = addii(mului(k+k+1,gel(g,k+2)), mului(k<<1,gel(g,k+1)));
    1011        1246 :   gel(h,n+2) = mului(n<<1, gel(g,n+1)); return h;
    1012             : }
    1013             : 
    1014             : #ifdef _MSC_VER /* Bill Daly: work around a MSVC bug */
    1015             : #pragma optimize("g",off)
    1016             : #endif
    1017             : /* P = polzagier(n,m)(-X), unnormalized (P(0) != 1) */
    1018             : static GEN
    1019          49 : polzag1(long n, long m)
    1020             : {
    1021          49 :   const long d = n - m, d2 = d<<1, r = (m+1)>>1, D = (d+1)>>1;
    1022             :   long i, k;
    1023          49 :   pari_sp av = avma;
    1024             :   GEN g, T;
    1025             : 
    1026          49 :   if (d <= 0 || m < 0) return pol_0(0);
    1027          49 :   g = cgetg(d+2, t_POL);
    1028          49 :   g[1] = evalsigne(1)|evalvarn(0);
    1029          49 :   T = cgetg(d+1,t_VEC);
    1030             :   /* T[k+1] = binomial(2d,2k+1), 0 <= k < d */
    1031          49 :   gel(T,1) = utoipos(d2);
    1032        1267 :   for (k = 1; k < D; k++)
    1033             :   {
    1034        1218 :     long k2 = k<<1;
    1035        2436 :     gel(T,k+1) = diviiexact(mulii(gel(T,k), muluu(d2-k2+1, d2-k2)),
    1036        1218 :                             muluu(k2,k2+1));
    1037             :   }
    1038          49 :   for (; k < d; k++) gel(T,k+1) = gel(T,d-k);
    1039          49 :   gel(g,2) = gel(T,d); /* binomial(2d, 2(d-1)+1) */
    1040        2499 :   for (i = 1; i < d; i++)
    1041             :   {
    1042        2450 :     pari_sp av2 = avma;
    1043        2450 :     GEN s, t = gel(T,d-i); /* binomial(2d, 2(d-1-i)+1) */
    1044        2450 :     s = t;
    1045      180082 :     for (k = d-i; k < d; k++)
    1046             :     {
    1047      177632 :       long k2 = k<<1;
    1048      177632 :       t = diviiexact(mulii(t, muluu(d2-k2+1, d-k)), muluu(k2+1,k-(d-i)+1));
    1049      177632 :       s = addii(s, t);
    1050             :     }
    1051             :     /* g_i = sum_{d-1-i <= k < d}, binomial(2*d, 2*k+1)*binomial(k,d-1-i) */
    1052        2450 :     gel(g,i+2) = gerepileuptoint(av2, s);
    1053             :   }
    1054             :   /* sum_{0 <= i < d} g_i x^i * (x+x^2)^r */
    1055          49 :   g = RgX_mulXn(gmul(g, gpowgs(deg1pol(gen_1,gen_1,0),r)), r);
    1056          49 :   if (!odd(m)) g = delt(g, n);
    1057        1274 :   for (i=1; i<=r; i++)
    1058             :   {
    1059        1225 :     g = delt(ZX_deriv(g), n);
    1060        1225 :     if (gc_needed(av,4))
    1061             :     {
    1062           0 :       if (DEBUGMEM>1) pari_warn(warnmem,"polzag, i = %ld/%ld", i,r);
    1063           0 :       g = gerepilecopy(av, g);
    1064             :     }
    1065             :   }
    1066          49 :   return g;
    1067             : }
    1068             : GEN
    1069          28 : polzag(long n, long m)
    1070             : {
    1071          28 :   pari_sp av = avma;
    1072          28 :   GEN g = ZX_z_unscale(polzag1(n,m), -1);
    1073          28 :   return gerepileupto(av, RgX_Rg_div(g,gel(g,2)));
    1074             : }
    1075             : 
    1076             : GEN
    1077           7 : sumalt(void *E, GEN (*eval)(void *, GEN), GEN a, long prec)
    1078             : {
    1079             :   ulong k, N;
    1080           7 :   pari_sp av = avma, av2;
    1081             :   GEN s, az, c, d;
    1082             : 
    1083           7 :   if (typ(a) != t_INT) pari_err_TYPE("sumalt",a);
    1084           7 :   N = (ulong)(0.39322*(prec2nbits(prec) + 7)); /*0.39322 > 1/log_2(3+sqrt(8))*/
    1085           7 :   d = powru(addsr(3, sqrtr(stor(8,prec))), N);
    1086           7 :   d = shiftr(addrr(d, invr(d)),-1);
    1087           7 :   a = setloop(a);
    1088           7 :   az = gen_m1; c = d;
    1089           7 :   s = gen_0;
    1090           7 :   av2 = avma;
    1091         371 :   for (k=0; ; k++) /* k < N */
    1092             :   {
    1093         371 :     c = addir(az,c); s = gadd(s, gmul(c, eval(E, a)));
    1094         371 :     if (k==N-1) break;
    1095         364 :     az = diviuuexact(muluui((N-k)<<1,N+k,az), k+1, (k<<1)+1);
    1096         364 :     a = incloop(a); /* in place! */
    1097         364 :     if (gc_needed(av,4))
    1098             :     {
    1099           0 :       if (DEBUGMEM>1) pari_warn(warnmem,"sumalt, k = %ld/%ld", k,N-1);
    1100           0 :       gerepileall(av2, 3, &az,&c,&s);
    1101             :     }
    1102         364 :   }
    1103           7 :   return gerepileupto(av, gdiv(s,d));
    1104             : }
    1105             : 
    1106             : GEN
    1107           7 : sumalt2(void *E, GEN (*eval)(void *, GEN), GEN a, long prec)
    1108             : {
    1109             :   long k, N;
    1110           7 :   pari_sp av = avma, av2;
    1111             :   GEN s, dn, pol;
    1112             : 
    1113           7 :   if (typ(a) != t_INT) pari_err_TYPE("sumalt",a);
    1114           7 :   N = (long)(0.307073*(prec2nbits(prec) + 5)); /*0.307073 > 1/log_2(\beta_B)*/
    1115           7 :   pol = ZX_div_by_X_1(polzag1(N,N>>1), &dn);
    1116           7 :   a = setloop(a);
    1117           7 :   N = degpol(pol);
    1118           7 :   s = gen_0;
    1119           7 :   av2 = avma;
    1120         280 :   for (k=0; k<=N; k++)
    1121             :   {
    1122         280 :     GEN t = itor(gel(pol,k+2), prec+EXTRAPRECWORD);
    1123         280 :     s = gadd(s, gmul(t, eval(E, a)));
    1124         280 :     if (k == N) break;
    1125         273 :     a = incloop(a); /* in place! */
    1126         273 :     if (gc_needed(av,4))
    1127             :     {
    1128           0 :       if (DEBUGMEM>1) pari_warn(warnmem,"sumalt2, k = %ld/%ld", k,N-1);
    1129           0 :       s = gerepileupto(av2, s);
    1130             :     }
    1131             :   }
    1132           7 :   return gerepileupto(av, gdiv(s,dn));
    1133             : }
    1134             : 
    1135             : GEN
    1136          14 : sumalt0(GEN a, GEN code, long flag, long prec)
    1137             : {
    1138          14 :   switch(flag)
    1139             :   {
    1140           7 :     case 0: EXPR_WRAP(code, sumalt (EXPR_ARG,a,prec));
    1141           7 :     case 1: EXPR_WRAP(code, sumalt2(EXPR_ARG,a,prec));
    1142           0 :     default: pari_err_FLAG("sumalt");
    1143             :   }
    1144             :   return NULL; /* LCOV_EXCL_LINE */
    1145             : }
    1146             : 
    1147             : /* For k > 0, set S[k*2^i] <- g(k*2^i), k*2^i <= N = #S.
    1148             :  * Only needed with k odd (but also works for g even). */
    1149             : static void
    1150        7406 : binsum(GEN S, ulong k, void *E, GEN (*f)(void *, GEN), GEN a,
    1151             :         long G, long prec)
    1152             : {
    1153        7406 :   long e, i, N = lg(S)-1, l = expu(N / k); /* k 2^l <= N < k 2^(l+1) */
    1154             :   pari_sp av;
    1155        7406 :   GEN r, t = gen_0;
    1156             : 
    1157        7406 :   gel(S, k << l) = cgetr(prec); av = avma;
    1158        7406 :   G -= l;
    1159        7406 :   r = utoipos(k<<l);
    1160     3964471 :   for(e=0;;e++) /* compute g(k 2^l) with absolute error ~ 2^(G-l) */
    1161             :   {
    1162     3964471 :     GEN u = gtofp(f(E, addii(a,r)), prec);
    1163     3964471 :     if (typ(u) != t_REAL) pari_err_TYPE("sumpos",u);
    1164     3964471 :     if (!signe(u)) break;
    1165     3964282 :     if (!e)
    1166        7217 :       t = u;
    1167             :     else {
    1168     3957065 :       shiftr_inplace(u, e);
    1169     3957065 :       t = addrr(t,u);
    1170     3957065 :       if (expo(u) < G) break;
    1171             :     }
    1172     3957065 :     r = shifti(r,1);
    1173     3957065 :   }
    1174        7406 :   gel(S, k << l) = t = gerepileuptoleaf(av, t);
    1175             :   /* g(j) = 2g(2j) + f(a+j) for all j > 0 */
    1176       14812 :   for(i = l-1; i >= 0; i--)
    1177             :   { /* t ~ g(2 * k*2^i) with error ~ 2^(G-i-1) */
    1178             :     GEN u;
    1179        7406 :     av = avma; u = gtofp(f(E, addiu(a, k << i)), prec);
    1180        7406 :     if (typ(u) != t_REAL) pari_err_TYPE("sumpos",u);
    1181        7406 :     t = addrr(gtofp(u,prec), shiftr(t,1)); /* ~ g(k*2^i) */
    1182        7406 :     gel(S, k << i) = t = gerepileuptoleaf(av, t);
    1183             :   }
    1184        7406 : }
    1185             : /* For k > 0, let g(k) := \sum_{e >= 0} 2^e f(a + k*2^e).
    1186             :  * Return [g(k), 1 <= k <= N] */
    1187             : static GEN
    1188          70 : sumpos_init(void *E, GEN (*f)(void *, GEN), GEN a, long N, long prec)
    1189             : {
    1190          70 :   GEN S = cgetg(N+1,t_VEC);
    1191          70 :   long k, G = -prec2nbits(prec) - 5;
    1192          70 :   for (k=1; k<=N; k+=2) binsum(S,k, E,f, a,G,prec);
    1193          70 :   return S;
    1194             : }
    1195             : 
    1196             : GEN
    1197          56 : sumpos(void *E, GEN (*eval)(void *, GEN), GEN a, long prec)
    1198             : {
    1199             :   ulong k, N;
    1200          56 :   pari_sp av = avma;
    1201             :   GEN s, az, c, d, S;
    1202             : 
    1203          56 :   if (typ(a) != t_INT) pari_err_TYPE("sumpos",a);
    1204          56 :   a = subiu(a,1);
    1205          56 :   N = (ulong)(0.4*(prec2nbits(prec) + 7));
    1206          56 :   if (odd(N)) N++; /* extra precision for free */
    1207          56 :   d = powru(addsr(3, sqrtr(stor(8,prec))), N);
    1208          56 :   d = shiftr(addrr(d, invr(d)),-1);
    1209          56 :   az = gen_m1; c = d;
    1210             : 
    1211          56 :   S = sumpos_init(E, eval, a, N, prec);
    1212          56 :   s = gen_0;
    1213       10360 :   for (k=0; k<N; k++)
    1214             :   {
    1215             :     GEN t;
    1216       10360 :     c = addir(az,c);
    1217       10360 :     t = mulrr(gel(S,k+1), c);
    1218       10360 :     s = odd(k)? mpsub(s, t): mpadd(s, t);
    1219       10360 :     if (k == N-1) break;
    1220       10304 :     az = diviuuexact(muluui((N-k)<<1,N+k,az), k+1, (k<<1)+1);
    1221             :   }
    1222          56 :   return gerepileupto(av, gdiv(s,d));
    1223             : }
    1224             : 
    1225             : GEN
    1226          14 : sumpos2(void *E, GEN (*eval)(void *, GEN), GEN a, long prec)
    1227             : {
    1228             :   ulong k, N;
    1229          14 :   pari_sp av = avma;
    1230             :   GEN s, pol, dn, S;
    1231             : 
    1232          14 :   if (typ(a) != t_INT) pari_err_TYPE("sumpos2",a);
    1233          14 :   a = subiu(a,1);
    1234          14 :   N = (ulong)(0.31*(prec2nbits(prec) + 5));
    1235             : 
    1236          14 :   if (odd(N)) N++; /* extra precision for free */
    1237          14 :   S = sumpos_init(E, eval, a, N, prec);
    1238          14 :   pol = ZX_div_by_X_1(polzag1(N,N>>1), &dn);
    1239          14 :   s = gen_0;
    1240        4466 :   for (k=0; k<N; k++)
    1241             :   {
    1242        4452 :     GEN t = mulri(gel(S,k+1), gel(pol,k+2));
    1243        4452 :     s = odd(k)? mpsub(s,t): mpadd(s,t);
    1244             :   }
    1245          14 :   return gerepileupto(av, gdiv(s,dn));
    1246             : }
    1247             : 
    1248             : GEN
    1249          70 : sumpos0(GEN a, GEN code, long flag, long prec)
    1250             : {
    1251          70 :   switch(flag)
    1252             :   {
    1253          56 :     case 0: EXPR_WRAP(code, sumpos (EXPR_ARG,a,prec));
    1254          14 :     case 1: EXPR_WRAP(code, sumpos2(EXPR_ARG,a,prec));
    1255           0 :     default: pari_err_FLAG("sumpos");
    1256             :   }
    1257             :   return NULL; /* LCOV_EXCL_LINE */
    1258             : }
    1259             : 
    1260             : /********************************************************************/
    1261             : /**                                                                **/
    1262             : /**            SEARCH FOR REAL ZEROS of an expression              **/
    1263             : /**                                                                **/
    1264             : /********************************************************************/
    1265             : /* Brent's method, [a,b] bracketing interval */
    1266             : GEN
    1267         889 : zbrent(void *E, GEN (*eval)(void *, GEN), GEN a, GEN b, long prec)
    1268             : {
    1269             :   long sig, iter, itmax;
    1270         889 :   pari_sp av = avma;
    1271             :   GEN c, d, e, tol, fa, fb, fc;
    1272             : 
    1273         889 :   if (typ(a) != t_REAL || realprec(a) < prec) a = gtofp(a, prec);
    1274         889 :   if (typ(b) != t_REAL || realprec(b) < prec) b = gtofp(b, prec);
    1275         889 :   sig = cmprr(b, a);
    1276         889 :   if (!sig) return gerepileupto(av, a);
    1277         889 :   if (sig < 0) {c = a; a = b; b = c;} else c = b;
    1278         889 :   fa = eval(E, a);
    1279         889 :   fb = eval(E, b);
    1280         889 :   if (gsigne(fa)*gsigne(fb) > 0)
    1281           7 :     pari_err_DOMAIN("solve", "f(a)f(b)", ">", gen_0, mkvec2(fa, fb));
    1282         882 :   itmax = prec2nbits(prec) * 2 + 1;
    1283         882 :   tol = real2n(5-prec2nbits(prec), LOWDEFAULTPREC);
    1284         882 :   fc = fb;
    1285         882 :   e = d = NULL; /* gcc -Wall */
    1286       14442 :   for (iter = 1; iter <= itmax; ++iter)
    1287             :   {
    1288             :     GEN xm, tol1;
    1289       14442 :     if (gsigne(fb)*gsigne(fc) > 0)
    1290             :     {
    1291        7671 :       c = a; fc = fa; e = d = subrr(b, a);
    1292             :     }
    1293       14442 :     if (gcmp(gabs(fc, 0), gabs(fb, 0)) < 0)
    1294             :     {
    1295        3757 :       a = b; b = c; c = a; fa = fb; fb = fc; fc = fa;
    1296             :     }
    1297       14442 :     tol1 = abscmprr(tol, b) > 0? sqrr(tol): mulrr(tol, absr(b));
    1298       14442 :     xm = shiftr(subrr(c, b), -1);
    1299       14442 :     if (abscmprr(xm, tol1) <= 0 || gequal0(fb)) break; /* SUCCESS */
    1300             : 
    1301       13560 :     if (abscmprr(e, tol1) >= 0 && gcmp(gabs(fa, 0), gabs(fb, 0)) > 0)
    1302       12219 :     { /* attempt interpolation */
    1303       12219 :       GEN min1, min2, p, q, s = gdiv(fb, fa);
    1304       12219 :       if (cmprr(a, c) == 0)
    1305             :       {
    1306        7357 :         p = gmul2n(gmul(xm, s), 1);
    1307        7357 :         q = gsubsg(1, s);
    1308             :       }
    1309             :       else
    1310             :       {
    1311        4862 :         GEN r = gdiv(fb, fc);
    1312        4862 :         q = gdiv(fa, fc);
    1313        4862 :         p = gmul2n(gmul(gsub(q, r), gmul(xm, q)), 1);
    1314        4862 :         p = gmul(s, gsub(p, gmul(gsub(b, a), gsubgs(r, 1))));
    1315        4862 :         q = gmul(gmul(gsubgs(q, 1), gsubgs(r, 1)), gsubgs(s, 1));
    1316             :       }
    1317       12219 :       if (gsigne(p) > 0) q = gneg_i(q); else p = gneg_i(p);
    1318       12219 :       min1 = gsub(gmulsg(3, gmul(xm,q)), gabs(gmul(q, tol1), 0));
    1319       12219 :       min2 = gabs(gmul(e, q), 0);
    1320       12219 :       if (gcmp(gmul2n(p, 1), gmin(min1, min2)) < 0)
    1321       10276 :         { e = d; d = gdiv(p, q); } /* interpolation OK */
    1322             :       else
    1323        1943 :         { d = xm; e = d; } /* failed, use bisection */
    1324             :     }
    1325        1341 :     else { d = xm; e = d; } /* bound decreasing too slowly, use bisection */
    1326       13560 :     a = b; fa = fb;
    1327       13560 :     if (gcmp(gabs(d, 0), tol1) > 0) b = gadd(b, d);
    1328         844 :     else if (gsigne(xm) > 0)      b = addrr(b, tol1);
    1329         491 :     else                          b = subrr(b, tol1);
    1330       13560 :     if (realprec(b) < prec) b = rtor(b, prec);
    1331       13560 :     fb = eval(E, b);
    1332             :   }
    1333         882 :   if (iter > itmax) pari_err_IMPL("solve recovery [too many iterations]");
    1334         882 :   return gerepileuptoleaf(av, rcopy(b));
    1335             : }
    1336             : 
    1337             : GEN
    1338          21 : zbrent0(GEN a, GEN b, GEN code, long prec)
    1339          21 : { EXPR_WRAP(code, zbrent(EXPR_ARG, a, b, prec)); }
    1340             : 
    1341             : /* x = solve_start(&D, a, b, prec)
    1342             :  * while (x) {
    1343             :  *   y = ...(x);
    1344             :  *   x = solve_next(&D, y);
    1345             :  * }
    1346             :  * return D.res; */
    1347             : 
    1348             : /* Find zeros of a function in the real interval [a,b] by interval splitting */
    1349             : GEN
    1350          91 : solvestep(void *E, GEN (*f)(void *,GEN), GEN a, GEN b, GEN step, long flag, long prec)
    1351             : {
    1352          91 :   const long ITMAX = 10;
    1353          91 :   pari_sp av = avma;
    1354             :   GEN fa, ainit, binit;
    1355          91 :   long sainit, it, bit = prec2nbits(prec) / 2, ct = 0, s = gcmp(a,b);
    1356             : 
    1357          91 :   if (!s) return gequal0(f(E, a)) ? gcopy(mkvec(a)): cgetg(1,t_VEC);
    1358          91 :   if (s > 0) swap(a, b);
    1359          91 :   if (flag&4)
    1360             :   {
    1361          84 :     if (gcmpgs(step,1)<=0) pari_err_DOMAIN("solvestep","step","<=",gen_1,step);
    1362          84 :     if (gsigne(a) <= 0) pari_err_DOMAIN("solvestep","a","<=",gen_0,a);
    1363             :   }
    1364           7 :   else if (gsigne(step) <= 0)
    1365           0 :     pari_err_DOMAIN("solvestep","step","<=",gen_0,step);
    1366          91 :   ainit = a = gtofp(a, prec); fa = f(E, a);
    1367          91 :   binit = b = gtofp(b, prec); step = gtofp(step, prec);
    1368          91 :   sainit = gsigne(fa);
    1369          91 :   if (gexpo(fa) < -bit) sainit = 0;
    1370          98 :   for (it = 0; it < ITMAX; it++)
    1371             :   {
    1372          98 :     pari_sp av2 = avma;
    1373          98 :     GEN v = cgetg(1, t_VEC);
    1374             :     long sa;
    1375          98 :     a = ainit;
    1376          98 :     b = binit;
    1377          98 :     sa = sainit;
    1378        2457 :     while (gcmp(a,b) < 0)
    1379             :     {
    1380        2261 :       GEN fc, c = (flag&4)? gmul(a, step): gadd(a, step);
    1381             :       long sc;
    1382        2261 :       if (gcmp(c,b) > 0) c = b;
    1383        2261 :       fc = f(E, c);
    1384        2261 :       sc = gsigne(fc);
    1385        2261 :       if (gexpo(fc) < -bit) sc = 0;
    1386        2261 :       if (!sc || sa*sc < 0)
    1387             :       {
    1388             :         long e;
    1389             :         GEN z;
    1390         441 :         z = sc? zbrent(E, f, a, c, prec): c;
    1391         441 :         (void)grndtoi(z, &e);
    1392         441 :         if (e  <= -bit) ct = 1;
    1393         441 :         if ((flag&1) && ((!(flag&8)) || ct)) return gerepileupto(av, z);
    1394         441 :         v = gconcat(v, z);
    1395             :       }
    1396        2261 :       a = c; fa = fc; sa = sc;
    1397             :     }
    1398          98 :     if ((!(flag&2) || lg(v) > 1) && (!(flag&8) || ct))
    1399          91 :       return gerepilecopy(av, v);
    1400           7 :     step = (flag&4)? sqrtr(sqrtr(step)): gmul2n(step, -2);
    1401           7 :     gerepileall(av2, 2, &fa, &step);
    1402             :   }
    1403           0 :   if (it == ITMAX) pari_err_IMPL("solvestep recovery [too many iterations]");
    1404           0 :   return NULL;
    1405             : }
    1406             : 
    1407             : GEN
    1408           7 : solvestep0(GEN a, GEN b, GEN step, GEN code, long flag, long prec)
    1409           7 : { EXPR_WRAP(code, solvestep(EXPR_ARG, a,b, step, flag, prec)); }
    1410             : 
    1411             : /********************************************************************/
    1412             : /**                     Numerical derivation                       **/
    1413             : /********************************************************************/
    1414             : 
    1415             : struct deriv_data
    1416             : {
    1417             :   GEN code;
    1418             :   GEN args;
    1419             : };
    1420             : 
    1421         147 : static GEN deriv_eval(void *E, GEN x, long prec)
    1422             : {
    1423         147 :  struct deriv_data *data=(struct deriv_data *)E;
    1424         147 :  gel(data->args,1)=x;
    1425         147 :  return closure_callgenvecprec(data->code, data->args, prec);
    1426             : }
    1427             : 
    1428             : /* Rationale: (f(2^-e) - f(-2^-e) + O(2^-pr)) / (2 * 2^-e) = f'(0) + O(2^-2e)
    1429             :  * since 2nd derivatives cancel.
    1430             :  *   prec(LHS) = pr - e
    1431             :  *   prec(RHS) = 2e, equal when  pr = 3e = 3/2 fpr (fpr = required final prec)
    1432             :  *
    1433             :  * For f'(x), x far from 0: prec(LHS) = pr - e - expo(x)
    1434             :  * --> pr = 3/2 fpr + expo(x) */
    1435             : GEN
    1436        1029 : derivnum(void *E, GEN (*eval)(void *, GEN, long), GEN x, long prec)
    1437             : {
    1438             :   GEN eps,a,b, y;
    1439             :   long pr, l, e, ex, newprec;
    1440        1029 :   pari_sp av = avma;
    1441        1029 :   long p = precision(x);
    1442        1029 :   long fpr = p ? prec2nbits(p): prec2nbits(prec);
    1443        1029 :   ex = gexpo(x);
    1444        1029 :   if (ex < 0) ex = 0; /* near 0 */
    1445        1029 :   pr = (long)ceil(fpr * 1.5 + ex);
    1446        1029 :   l = nbits2prec(pr);
    1447        1029 :   newprec = nbits2prec(pr + ex + BITS_IN_LONG);
    1448        1029 :   switch(typ(x))
    1449             :   {
    1450             :     case t_REAL:
    1451             :     case t_COMPLEX:
    1452         427 :       x = gprec_w(x, newprec);
    1453             :   }
    1454             : 
    1455        1029 :   e = fpr/2; /* 1/2 required prec (in sig. bits) */
    1456        1029 :   eps = real2n(-e, l);
    1457        1029 :   a = eval(E, gsub(x, eps), newprec);
    1458        1029 :   b = eval(E, gadd(x, eps), newprec);
    1459        1029 :   y = gmul2n(gsub(b,a), e-1);
    1460        1029 :   return gerepileupto(av, gprec_w(y, nbits2prec(fpr)));
    1461             : }
    1462             : 
    1463             : /* Fornberg interpolation algorithm for finite differences coefficients
    1464             : * using N+1 equidistant grid points around 0 [ assume N even >= M ].
    1465             : * Compute \delta[m]_{N,nu} for all derivation orders m = 0..M such that
    1466             : *   h^m * f^{(m)}(0) = \sum_{nu = 0}^n delta[m]_{N,nu}  f(a_nu) + O(h^{N-m+1}),
    1467             : * for step size h.
    1468             : * Return a = [0,-1,1...,-N2,N2] and vector of vectors d: d[m+1][nu+1]
    1469             : * = (M!/m!) * delta[m]_{N,nu}, nu = 0..N */
    1470             : static void
    1471          28 : FD(long M, long N, GEN *pd, GEN *pa)
    1472             : {
    1473             :   GEN d, a, b, W, Wp, t, F, Mfact;
    1474             :   long N2, m, nu, i;
    1475             : 
    1476          28 :   if (odd(N)) N++; /* make it even */
    1477          28 :   N2 = N>>1;
    1478          28 :   F = cgetg(N+2, t_VEC);
    1479          28 :   a = cgetg(N+2, t_VEC);
    1480          28 :   b = cgetg(N2+1, t_VEC);
    1481          28 :   gel(a,1) = gen_0;
    1482         357 :   for (i = 1; i <= N2; i++)
    1483             :   {
    1484         329 :     gel(a,2*i)   = utoineg(i);
    1485         329 :     gel(a,2*i+1) = utoipos(i);
    1486         329 :     gel(b,i) = sqru(i);
    1487             :   }
    1488             :   /* w = \prod (X - a[i]) = x W(x^2) */
    1489          28 :   Mfact = mpfact(M);
    1490          28 :   W = roots_to_pol(b, 0);
    1491          28 :   Wp = ZX_deriv(W);
    1492          28 :   t = gel(W,2); /* w'(0) */
    1493          28 :   t = diviiexact(t, Mfact);
    1494          28 :   gel(F,1) = RgX_Rg_div(RgX_inflate(W,2), t);
    1495         357 :   for (i = 1; i <= N2; i++)
    1496             :   { /* t = w'(a_{2i}) = w'(a_{2i+1}) */
    1497         329 :     GEN r, t = mulii(shifti(gel(b,i),1), poleval(Wp, gel(b,i)));
    1498             :     GEN U, S, T;
    1499         329 :     U = RgX_inflate(RgX_div_by_X_x(W, gel(b,i), &r), 2);
    1500         329 :     U = RgX_shift_shallow(U, 1);
    1501         329 :     U = RgXn_red_shallow(U, M+1); /* higher terms not needed */
    1502         329 :     t = diviiexact(t, Mfact);
    1503         329 :     U = RgX_Rg_div(U, t);
    1504         329 :     S = RgX_shift_shallow(U,1);
    1505         329 :     T = RgX_Rg_mul(U, gel(a,2*i+1));
    1506         329 :     gel(F,2*i)   = RgX_sub(S, T);
    1507         329 :     gel(F,2*i+1) = RgX_add(S, T);
    1508             :   }
    1509             :   /* F[i] = M! w(X) / ((X-a[i])w'(a[i])) + O(X^(M+1)) */
    1510          28 :   d = cgetg(M+2, t_VEC);
    1511         280 :   for (m = 0; m <= M; m++)
    1512             :   {
    1513         252 :     GEN v = cgetg(N+2, t_VEC); /* coeff(F[nu],X^m) */
    1514         252 :     for (nu = 0; nu <= N; nu++) gel(v, nu+1) = gmael(F, nu+1, m+2);
    1515         252 :     gel(d,m+1) = v;
    1516             :   }
    1517          28 :   *pd = d;
    1518          28 :   *pa = a;
    1519          28 : }
    1520             : 
    1521             : static void
    1522         231 : chk_ord(long m)
    1523             : {
    1524         231 :   if (m < 0)
    1525          14 :     pari_err_DOMAIN("derivnumk", "derivation order", "<", gen_0, stoi(m));
    1526         217 : }
    1527             : 
    1528             : GEN
    1529          42 : derivnumk(void *E, GEN (*eval)(void *, GEN, long), GEN x, GEN ind0, long prec)
    1530             : {
    1531             :   GEN A, D, X, F, ind;
    1532             :   long M, fpr, p, i, pr, l, lA, e, ex, eD, newprec;
    1533          42 :   pari_sp av = avma;
    1534          42 :   int allodd = 1;
    1535             : 
    1536          42 :   ind = gtovecsmall(ind0);
    1537          42 :   l = lg(ind);
    1538          42 :   F = cgetg(l, t_VEC);
    1539          42 :   M = vecsmall_max(ind);
    1540          42 :   chk_ord(M);
    1541          42 :   if (!M) /* silly degenerate case */
    1542             :   {
    1543          14 :     X = eval(E, x, prec);
    1544          14 :     for (i = 1; i < l; i++) { chk_ord(ind[i]); gel(F,i) = X; }
    1545           7 :     if (typ(ind0) == t_INT) F = gel(F,1);
    1546           7 :     return gerepilecopy(av, F);
    1547             :   }
    1548          28 :   FD(M, 3*M-1, &D,&A); /* optimal if 'eval' uses quadratic time */
    1549             : 
    1550          28 :   p = precision(x);
    1551          28 :   fpr = p ? prec2nbits(p): prec2nbits(prec);
    1552          28 :   eD = gexpo(gel(D,M));
    1553          28 :   e = (fpr + 3*M*log2((double)M)) / (2*M);
    1554          28 :   ex = gexpo(x);
    1555          28 :   if (ex < 0) ex = 0; /* near 0 */
    1556          28 :   pr = (long)ceil(fpr + e * M); /* ~ 3fpr/2 */
    1557          28 :   newprec = nbits2prec(pr + eD + ex + BITS_IN_LONG);
    1558          28 :   switch(typ(x))
    1559             :   {
    1560             :     case t_REAL:
    1561             :     case t_COMPLEX:
    1562           0 :       x = gprec_w(x, newprec);
    1563             :   }
    1564          28 :   lA = lg(A); X = cgetg(lA, t_VEC);
    1565          42 :   for (i = 1; i < l; i++)
    1566          35 :     if (!odd(ind[i])) { allodd = 0; break; }
    1567             :   /* if only odd derivation orders, the value at 0 (A[1]) is not needed */
    1568          28 :   gel(X, 1) = gen_0;
    1569         707 :   for (i = allodd? 2: 1; i < lA; i++)
    1570             :   {
    1571         679 :     GEN t = eval(E, gadd(x, gmul2n(gel(A,i), -e)), newprec);
    1572         679 :     if (!gprecision(t)) t = gtofp(t, newprec);
    1573         679 :     gel(X, i) = t;
    1574             :   }
    1575             : 
    1576         147 :   for (i = 1; i < l; i++)
    1577             :   {
    1578             :     GEN t;
    1579         126 :     long m = ind[i]; chk_ord(m);
    1580         119 :     t = gmul2n(RgV_dotproduct(gel(D,m+1), X), e*m);
    1581         119 :     if (m < M) t = gdiv(t, mulu_interval(m+1,M));
    1582         119 :     gel(F,i) = t;
    1583             :   }
    1584          21 :   if (typ(ind0) == t_INT) F = gel(F,1);
    1585          21 :   return gerepileupto(av, gprec_w(F, nbits2prec(fpr)));
    1586             : }
    1587             : /* v(t') */
    1588             : static long
    1589          14 : rfrac_val_deriv(GEN t)
    1590             : {
    1591          14 :   long v = varn(gel(t,2));
    1592          14 :   return gvaluation(deriv(t, v), pol_x(v));
    1593             : }
    1594             : 
    1595             : GEN
    1596        1043 : derivfunk(void *E, GEN (*eval)(void *, GEN, long), GEN x, GEN ind0, long prec)
    1597             : {
    1598             :   pari_sp av;
    1599             :   GEN ind, xp, ixp, F, G;
    1600             :   long i, l, vx, M;
    1601        1043 :   if (!ind0) return derivfun(E, eval, x, prec);
    1602          63 :   switch(typ(x))
    1603             :   {
    1604             :   case t_REAL: case t_INT: case t_FRAC: case t_COMPLEX:
    1605          42 :     return derivnumk(E,eval, x, ind0, prec);
    1606             :   case t_POL:
    1607          14 :     ind = gtovecsmall(ind0);
    1608          14 :     M = vecsmall_max(ind);
    1609          14 :     xp = RgX_deriv(x);
    1610          14 :     x = RgX_to_ser(x, precdl+2 + M * (1+RgX_val(xp)));
    1611          14 :     break;
    1612             :   case t_RFRAC:
    1613           7 :     ind = gtovecsmall(ind0);
    1614           7 :     M = vecsmall_max(ind);
    1615           7 :     x = rfrac_to_ser(x, precdl+2 + M * (1+rfrac_val_deriv(x)));
    1616           7 :     xp = derivser(x);
    1617           7 :     break;
    1618             :   case t_SER:
    1619           0 :     ind = gtovecsmall(ind0);
    1620           0 :     M = vecsmall_max(ind);
    1621           0 :     xp = derivser(x);
    1622           0 :     break;
    1623           0 :   default: pari_err_TYPE("numerical derivation",x);
    1624             :     return NULL; /*LCOV_EXCL_LINE*/
    1625             :   }
    1626          21 :   av = avma; chk_ord(M);
    1627          21 :   vx = varn(x);
    1628          21 :   ixp = M? ginv(xp): NULL;
    1629          21 :   F = cgetg(M+2, t_VEC);
    1630          21 :   gel(F,1) = eval(E, x, prec);
    1631          21 :   for (i = 1; i <= M; i++) gel(F,i+1) = gmul(deriv(gel(F,i),vx), ixp);
    1632          21 :   l = lg(ind); G = cgetg(l, t_VEC);
    1633          42 :   for (i = 1; i < l; i++)
    1634             :   {
    1635          21 :     long m = ind[i]; chk_ord(m);
    1636          21 :     gel(G,i) = gel(F,m+1);
    1637             :   }
    1638          21 :   if (typ(ind0) == t_INT) G = gel(G,1);
    1639          21 :   return gerepilecopy(av, G);
    1640             : }
    1641             : 
    1642             : GEN
    1643        1064 : derivfun(void *E, GEN (*eval)(void *, GEN, long), GEN x, long prec)
    1644             : {
    1645        1064 :   pari_sp av = avma;
    1646             :   GEN xp;
    1647             :   long vx;
    1648        1064 :   switch(typ(x))
    1649             :   {
    1650             :   case t_REAL: case t_INT: case t_FRAC: case t_COMPLEX:
    1651        1029 :     return derivnum(E,eval, x, prec);
    1652             :   case t_POL:
    1653          14 :     xp = RgX_deriv(x);
    1654          14 :     x = RgX_to_ser(x, precdl+2+ (1 + RgX_val(xp)));
    1655          14 :     break;
    1656             :   case t_RFRAC:
    1657           7 :     x = rfrac_to_ser(x, precdl+2+ (1 + rfrac_val_deriv(x)));
    1658             :     /* fall through */
    1659             :   case t_SER:
    1660          14 :     xp = derivser(x);
    1661          14 :     break;
    1662           7 :   default: pari_err_TYPE("formal derivation",x);
    1663             :     return NULL; /*LCOV_EXCL_LINE*/
    1664             :   }
    1665          28 :   vx = varn(x);
    1666          28 :   return gerepileupto(av, gdiv(deriv(eval(E, x, prec),vx), xp));
    1667             : }
    1668             : 
    1669             : GEN
    1670          21 : laurentseries(void *E, GEN (*f)(void*,GEN x, long), long M, long v, long prec)
    1671             : {
    1672          21 :   pari_sp av = avma;
    1673             :   long d;
    1674             : 
    1675          21 :   if (v < 0) v = 0;
    1676          21 :   d = maxss(M+1,1);
    1677             :   for (;;)
    1678             :   {
    1679             :     long i, dr, vr;
    1680             :     GEN s;
    1681          35 :     s = cgetg(d+2, t_SER); s[1] = evalsigne(1) | evalvalp(1) | evalvarn(v);
    1682          35 :     gel(s, 2) = gen_1; for (i = 3; i <= d+1; i++) gel(s, i) = gen_0;
    1683          35 :     s = f(E, s, prec);
    1684          35 :     if (typ(s) != t_SER || varn(s) != v) pari_err_TYPE("laurentseries", s);
    1685          35 :     vr = valp(s);
    1686          35 :     if (M < vr) { avma = av; return zeroser(v, M); }
    1687          35 :     dr = lg(s) + vr - 3 - M;
    1688          35 :     if (dr >= 0) return gerepileupto(av, s);
    1689          14 :     avma = av; d -= dr;
    1690          14 :   }
    1691             : }
    1692             : static GEN
    1693          35 : _evalclosprec(void *E, GEN x, long prec)
    1694             : {
    1695             :   GEN s;
    1696          35 :   push_localprec(prec); s = closure_callgen1((GEN)E, x);
    1697          35 :   pop_localprec(); return s;
    1698             : }
    1699             : #define CLOS_ARGPREC __E, &_evalclosprec
    1700             : GEN
    1701          21 : laurentseries0(GEN f, long M, long v, long prec)
    1702             : {
    1703          21 :   if (typ(f) != t_CLOSURE || closure_arity(f) != 1 || closure_is_variadic(f))
    1704           0 :     pari_err_TYPE("laurentseries",f);
    1705          21 :   EXPR_WRAP(f, laurentseries(CLOS_ARGPREC,M,v,prec));
    1706             : }
    1707             : 
    1708             : GEN
    1709        1043 : derivnum0(GEN a, GEN code, GEN ind, long prec)
    1710        1043 : { EXPR_WRAP(code, derivfunk(EXPR_ARGPREC,a,ind,prec)); }
    1711             : 
    1712             : GEN
    1713          84 : derivfun0(GEN code, GEN args, long prec)
    1714             : {
    1715             :   struct deriv_data E;
    1716          84 :   E.code=code; E.args=args;
    1717          84 :   return derivfun((void*)&E, deriv_eval, gel(args,1), prec);
    1718             : }
    1719             : 
    1720             : /********************************************************************/
    1721             : /**                   Numerical extrapolation                      **/
    1722             : /********************************************************************/
    1723             : /* for t_CLOSURE */
    1724             : static double
    1725          84 : fun_getmf(long mul)
    1726             : {
    1727          84 :   const double A[] = {0,
    1728             :     0.331,0.260,0.228,0.210,0.198,
    1729             :     0.188,0.181,0.175,0.170,0.166 };
    1730          84 :   const double B[] = {0,
    1731             :     0.166,0.142,0.132,0.125,0.120,
    1732             :     0.116,0.113,0.111,0.109,0.107 };
    1733          84 :   if (mul <=  10) return A[mul];
    1734          77 :   if (mul <= 109) return B[mul/10]; else return 0.105;
    1735             : }
    1736             : /* for t_VEC */
    1737             : static double
    1738          14 : vec_getmf(long mul)
    1739             : {
    1740          14 :   const double A[] = {0,
    1741             :     0.2062, 0.1760, 0.1613, 0.1519, 0.1459,
    1742             :     0.1401, 0.1360, 0.1327, 0.1299, 0.1280 };
    1743          14 :   const double B[] = {0,
    1744             :     0.1280, 0.1133, 0.1064, 0.1019, 0.0987,
    1745             :     0.0962, 0.0942, 0.0925, 0.0911, 0.0899 };
    1746          14 :   if (mul <=  10) return A[mul];
    1747          14 :   if (mul <= 109) return B[mul/10]; else return 0.0899;
    1748             : }
    1749             : 
    1750             : /* [u(n*muli), u <= N], muli = 1 unless f!=NULL */
    1751             : static GEN
    1752          98 : get_u(void *E, GEN (*f)(void *, GEN, long), long N, long muli, long prec)
    1753             : {
    1754             :   long n;
    1755          98 :   GEN u = cgetg(N+1, t_VEC);
    1756          98 :   if (f)
    1757             :   {
    1758          84 :     for (n = 1; n <= N; n++) gel(u,n) = f(E, stoi(muli*n), prec);
    1759             :   }
    1760             :   else
    1761             :   {
    1762          14 :     GEN v = (GEN)E;
    1763          14 :     long t = lg(v)-1;
    1764          14 :     if (t < N*muli) pari_err_COMPONENT("limitnum","<",stoi(N), stoi(t));
    1765          14 :     for (n = 1; n <= N; n++) gel(u,n) = gel(v, muli*n);
    1766             :   }
    1767        4137 :   for (n = 1; n <= N; n++)
    1768             :   {
    1769        4039 :     GEN un = gel(u,n);
    1770        4039 :     if (is_rational_t(typ(un))) gel(u,n) = gtofp(un, prec);
    1771             :   }
    1772          98 :   return u;
    1773             : }
    1774             : 
    1775             : struct limit
    1776             : {
    1777             :   long prec0; /* target accuracy */
    1778             :   long prec; /* working accuracy */
    1779             :   long N; /* number of terms */
    1780             :   GEN u; /* sequence to extrapolate */
    1781             :   GEN na; /* [n^alpha, n <= N] */
    1782             :   GEN nma; /* [n^-alpha, n <= N] or NULL (alpha = 1) */
    1783             :   GEN coef; /* or NULL (alpha != 1) */
    1784             : };
    1785             : 
    1786             : static void
    1787          98 : limit_init(struct limit *L, void *E, GEN (*f)(void*,GEN,long),
    1788             :            long muli, GEN alpha, long prec)
    1789             : {
    1790          98 :   long bitprec = prec2nbits(prec), n, N;
    1791             :   GEN na;
    1792             : 
    1793          98 :   if (muli <= 0) muli = 20;
    1794          98 :   L->N = N = (long)ceil((f? fun_getmf(muli): vec_getmf(muli)) * bitprec);
    1795          98 :   L->prec = nbits2prec(bitprec + (long)ceil(1.844*N));
    1796          98 :   L->prec0 = prec;
    1797          98 :   L->u = get_u(E, f, N, muli, L->prec);
    1798          98 :   if (alpha && !gequal1(alpha))
    1799          28 :   {
    1800          28 :     long prec2 = gprecision(alpha);
    1801             :     GEN nma;
    1802          28 :     if (!prec2) prec2 = L->prec;
    1803          28 :     na = vecpowug(N, alpha, prec2);
    1804          28 :     L->coef = NULL;
    1805          28 :     L->nma = nma = cgetg(N+1, t_VEC);
    1806          28 :     for (n = 1; n <= N; n++) gel(nma, n) = ginv(gel(na, n));
    1807          28 :     if (muli != 1) na = gmul(na, gpow(utor(muli,prec2), alpha, prec2));
    1808             :   }
    1809             :   else
    1810             :   {
    1811          70 :     GEN coef, C = vecbinomial(N), T = vecpowuu(N, N);
    1812          70 :     na = cgetg(N+1, t_VEC);
    1813          70 :     L->coef = coef = cgetg(N+1, t_VEC);
    1814          70 :     L->nma = NULL;
    1815        2506 :     for (n = 1; n <= N; n++)
    1816             :     {
    1817        2436 :       GEN c = mulii(gel(C,n+1), gel(T,n));
    1818        2436 :       if (odd(N-n)) togglesign_safe(&c);
    1819        2436 :       gel(coef, n) = c;
    1820        2436 :       gel(na, n) = utoipos(n*muli);
    1821             :     }
    1822             :   }
    1823          98 :   L->na = na;
    1824          98 : }
    1825             : 
    1826             : /* Zagier/Lagrange extrapolation */
    1827             : static GEN
    1828         784 : limitnum_i(struct limit *L)
    1829             : {
    1830         784 :   pari_sp av = avma;
    1831             :   GEN S;
    1832         784 :   if (L->nma)
    1833         343 :     S = polint(L->nma, L->u,gen_0,NULL);
    1834             :   else
    1835         441 :     S = gdiv(RgV_dotproduct(L->u,L->coef), mpfact(L->N));
    1836         784 :   return gerepilecopy(av, gprec_w(S, L->prec0));
    1837             : }
    1838             : GEN
    1839          49 : limitnum(void *E, GEN (*f)(void *, GEN, long), long muli, GEN alpha, long prec)
    1840             : {
    1841             :   struct limit L;
    1842          49 :   limit_init(&L, E,f, muli, alpha, prec);
    1843          49 :   return limitnum_i(&L);
    1844             : }
    1845             : GEN
    1846          49 : limitnum0(GEN u, long muli, GEN alpha, long prec)
    1847             : {
    1848          49 :   void *E = (void*)u;
    1849          49 :   GEN (*f)(void*,GEN,long) = NULL;
    1850          49 :   switch(typ(u))
    1851             :   {
    1852             :     case t_COL:
    1853           7 :     case t_VEC: break;
    1854          42 :     case t_CLOSURE: f = gp_callprec; break;
    1855           0 :     default: pari_err_TYPE("limitnum", u);
    1856             :   }
    1857          49 :   return limitnum(E,f, muli,alpha, prec);
    1858             : }
    1859             : 
    1860             : GEN
    1861          49 : asympnum(void *E, GEN (*f)(void *, GEN, long), long muli, GEN alpha, long prec)
    1862             : {
    1863          49 :   const long MAX = 100;
    1864          49 :   pari_sp av = avma;
    1865          49 :   GEN u, vres = vectrunc_init(MAX);
    1866          49 :   long i, B = prec2nbits(prec);
    1867          49 :   double LB = 0.9*expu(B); /* 0.9 and 0.95 below are heuristic */
    1868             :   struct limit L;
    1869          49 :   limit_init(&L, E,f, muli, alpha, prec);
    1870          49 :   if (alpha) LB *= gtodouble(alpha);
    1871          49 :   u = L.u;
    1872         735 :   for(i = 1; i <= MAX; i++)
    1873             :   {
    1874             :     GEN a, s, v, p, q;
    1875             :     long n;
    1876         735 :     s = limitnum_i(&L);
    1877             :     /* NOT bestappr: lindep properly ignores the lower bits */
    1878         735 :     v = lindep_bit(mkvec2(gen_1, s), maxss((long)(0.95*floor(B - i*LB)), 32));
    1879         735 :     if (lg(v) == 1) break;
    1880         728 :     p = negi(gel(v,1));
    1881         728 :     q = gel(v,2);
    1882         728 :     if (!signe(q)) break;
    1883         728 :     a = gdiv(p,q);
    1884         728 :     s = gsub(s, a);
    1885             :     /* |s|q^2 > eps */
    1886         728 :     if (!gequal0(s) && gexpo(s) + 2*expi(q) > -17) break;
    1887         686 :     vectrunc_append(vres, a);
    1888         686 :     for (n = 1; n <= L.N; n++) gel(u,n) = gmul(gsub(gel(u,n), a), gel(L.na,n));
    1889             :   }
    1890          49 :   return gerepilecopy(av, vres);
    1891             : }
    1892             : GEN
    1893          49 : asympnum0(GEN u, long muli, GEN alpha, long prec)
    1894             : {
    1895          49 :   void *E = (void*)u;
    1896          49 :   GEN (*f)(void*,GEN,long) = NULL;
    1897          49 :   switch(typ(u))
    1898             :   {
    1899             :     case t_COL:
    1900           7 :     case t_VEC: break;
    1901          42 :     case t_CLOSURE: f = gp_callprec; break;
    1902           0 :     default: pari_err_TYPE("asympnum", u);
    1903             :   }
    1904          49 :   return asympnum(E,f, muli,alpha, prec);
    1905             : }

Generated by: LCOV version 1.11