Code coverage tests

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

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

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

LCOV - code coverage report
Current view: top level - basemath - bibli2.c (source / functions) Hit Total Coverage
Test: PARI/GP v2.8.0 lcov report (development 19355-c7ae729) Lines: 1033 1069 96.6 %
Date: 2016-08-26 06:12:17 Functions: 92 96 95.8 %
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             : /*******************************************************************/
      18             : /**                                                               **/
      19             : /**                      SPECIAL POLYNOMIALS                      **/
      20             : /**                                                               **/
      21             : /*******************************************************************/
      22             : /* Tchebichev polynomial: T0=1; T1=X; T(n)=2*X*T(n-1)-T(n-2)
      23             :  * T(n) = (n/2) sum_{k=0}^{n/2} a_k x^(n-2k)
      24             :  *   where a_k = (-1)^k 2^(n-2k) (n-k-1)! / k!(n-2k)! is an integer
      25             :  *   and a_0 = 2^(n-1), a_k / a_{k-1} = - (n-2k+2)(n-2k+1) / 4k(n-k) */
      26             : GEN
      27        2156 : polchebyshev1(long n, long v) /* Assume 4*n < LONG_MAX */
      28             : {
      29             :   long k, l;
      30             :   pari_sp av;
      31             :   GEN q,a,r;
      32             : 
      33        2156 :   if (v<0) v = 0;
      34             :   /* polchebyshev(-n,1) = polchebyshev(n,1) */
      35        2156 :   if (n < 0) n = -n;
      36        2156 :   if (n==0) return pol_1(v);
      37        2135 :   if (n==1) return pol_x(v);
      38             : 
      39        2093 :   q = cgetg(n+3, t_POL); r = q + n+2;
      40        2093 :   a = int2n(n-1);
      41        2093 :   gel(r--,0) = a;
      42        2093 :   gel(r--,0) = gen_0;
      43       31955 :   for (k=1,l=n; l>1; k++,l-=2)
      44             :   {
      45       29862 :     av = avma;
      46       29862 :     a = diviuuexact(muluui(l, l-1, a), 4*k, n-k);
      47       29862 :     togglesign(a); a = gerepileuptoint(av, a);
      48       29862 :     gel(r--,0) = a;
      49       29862 :     gel(r--,0) = gen_0;
      50             :   }
      51        2093 :   q[1] = evalsigne(1) | evalvarn(v);
      52        2093 :   return q;
      53             : }
      54             : static void
      55          70 : polchebyshev1_eval_aux(long n, GEN x, GEN *pt1, GEN *pt2)
      56             : {
      57             :   GEN t1, t2, b;
      58          84 :   if (n == 1) { *pt1 = gen_1; *pt2 = x; return; }
      59          56 :   if (n == 0) { *pt1 = x; *pt2 = gen_1; return; }
      60          56 :   polchebyshev1_eval_aux((n+1) >> 1, x, &t1, &t2);
      61          56 :   b = gsub(gmul(gmul2n(t1,1), t2), x);
      62          56 :   if (odd(n)) { *pt1 = gadd(gmul2n(gsqr(t1), 1), gen_m1); *pt2 = b; }
      63          42 :   else        { *pt1 = b; *pt2 = gadd(gmul2n(gsqr(t2), 1), gen_m1); }
      64             : }
      65             : static GEN
      66          14 : polchebyshev1_eval(long n, GEN x)
      67             : {
      68             :   GEN t1, t2;
      69             :   long i, v;
      70             :   pari_sp av;
      71             : 
      72          14 :   if (n < 0) n = -n;
      73          14 :   if (n==0) return gen_1;
      74          14 :   if (n==1) return gcopy(x);
      75          14 :   av = avma;
      76          14 :   v = u_lvalrem(n, 2, (ulong*)&n);
      77          14 :   polchebyshev1_eval_aux((n+1)>>1, x, &t1, &t2);
      78          14 :   if (n != 1) t2 = gsub(gmul(gmul2n(t1,1), t2), x);
      79          14 :   for (i = 1; i <= v; i++) t2 = gadd(gmul2n(gsqr(t2), 1), gen_m1);
      80          14 :   return gerepileupto(av, t2);
      81             : }
      82             : 
      83             : /* Chebychev  polynomial of the second kind U(n,x): the coefficient in front of
      84             :  * x^(n-2*m) is (-1)^m * 2^(n-2m)*(n-m)!/m!/(n-2m)!  for m=0,1,...,n/2 */
      85             : GEN
      86        2135 : polchebyshev2(long n, long v)
      87             : {
      88             :   pari_sp av;
      89             :   GEN q, a, r;
      90             :   long m;
      91        2135 :   int neg = 0;
      92             : 
      93        2135 :   if (v<0) v = 0;
      94             :   /* polchebyshev(-n,2) = -polchebyshev(n-2,2) */
      95        2135 :   if (n < 0) {
      96        1050 :     if (n == -1) return zeropol(v);
      97        1029 :     neg = 1; n = -n-2;
      98             :   }
      99        2114 :   if (n==0) return neg ? scalar_ZX_shallow(gen_m1, v): pol_1(v);
     100             : 
     101        2072 :   q = cgetg(n+3, t_POL); r = q + n+2;
     102        2072 :   a = int2n(n);
     103        2072 :   if (neg) togglesign(a);
     104        2072 :   gel(r--,0) = a;
     105        2072 :   gel(r--,0) = gen_0;
     106       30807 :   for (m=1; 2*m<= n; m++)
     107             :   {
     108       28735 :     av = avma;
     109       28735 :     a = diviuuexact(muluui(n-2*m+2, n-2*m+1, a), 4*m, n-m+1);
     110       28735 :     togglesign(a); a = gerepileuptoint(av, a);
     111       28735 :     gel(r--,0) = a;
     112       28735 :     gel(r--,0) = gen_0;
     113             :   }
     114        2072 :   q[1] = evalsigne(1) | evalvarn(v);
     115        2072 :   return q;
     116             : }
     117             : static void
     118          91 : polchebyshev2_eval_aux(long n, GEN x, GEN *pu1, GEN *pu2)
     119             : {
     120             :   GEN u1, u2, u, mu1;
     121         112 :   if (n == 1) { *pu1 = gen_1; *pu2 = gmul2n(x,1); return; }
     122          70 :   if (n == 0) { *pu1 = gen_0; *pu2 = gen_1; return; }
     123          70 :   polchebyshev2_eval_aux(n >> 1, x, &u1, &u2);
     124          70 :   mu1 = gneg(u1);
     125          70 :   u = gmul(gadd(u2,u1), gadd(u2,mu1));
     126          70 :   if (odd(n)) { *pu1 = u; *pu2 = gmul(gmul2n(u2,1), gadd(gmul(x,u2), mu1)); }
     127          35 :   else        { *pu2 = u; *pu1 = gmul(gmul2n(u1,1), gadd(u2, gmul(x,mu1))); }
     128             : }
     129             : static GEN
     130          35 : polchebyshev2_eval(long n, GEN x)
     131             : {
     132             :   GEN u1, u2, mu1;
     133          35 :   long neg = 0;
     134             :   pari_sp av;
     135             : 
     136          35 :   if (n < 0) {
     137          14 :     if (n == -1) return gen_0;
     138           7 :     neg = 1; n = -n-2;
     139             :   }
     140          28 :   if (n==0) return neg ? gen_m1: gen_1;
     141          21 :   av = avma;
     142          21 :   polchebyshev2_eval_aux(n>>1, x, &u1, &u2);
     143          21 :   mu1 = gneg(u1);
     144          21 :   if (odd(n)) u2 = gmul(gmul2n(u2,1), gadd(gmul(x,u2), mu1));
     145          14 :   else        u2 = gmul(gadd(u2,u1), gadd(u2,mu1));
     146          21 :   if (neg) u2 = gneg(u2);
     147          21 :   return gerepileupto(av, u2);
     148             : }
     149             : 
     150             : GEN
     151        4284 : polchebyshev(long n, long kind, long v)
     152             : {
     153        4284 :   switch (kind)
     154             :   {
     155        2149 :     case 1: return polchebyshev1(n, v);
     156        2135 :     case 2: return polchebyshev2(n, v);
     157           0 :     default: pari_err_FLAG("polchebyshev");
     158             :   }
     159           0 :   return NULL; /* not reached */
     160             : }
     161             : GEN
     162        4333 : polchebyshev_eval(long n, long kind, GEN x)
     163             : {
     164        4333 :   if (!x) return polchebyshev(n, kind, 0);
     165          63 :   if (gequalX(x)) return polchebyshev(n, kind, varn(x));
     166          49 :   switch (kind)
     167             :   {
     168          14 :     case 1: return polchebyshev1_eval(n, x);
     169          35 :     case 2: return polchebyshev2_eval(n, x);
     170           0 :     default: pari_err_FLAG("polchebyshev");
     171             :   }
     172           0 :   return NULL; /* not reached */
     173             : }
     174             : 
     175             : /* Hermite polynomial H(n,x):  H(n+1) = 2x H(n) - 2n H(n-1)
     176             :  * The coefficient in front of x^(n-2*m) is
     177             :  * (-1)^m * n! * 2^(n-2m)/m!/(n-2m)!  for m=0,1,...,n/2.. */
     178             : GEN
     179        1428 : polhermite(long n, long v)
     180             : {
     181             :   long m;
     182             :   pari_sp av;
     183             :   GEN q,a,r;
     184             : 
     185        1428 :   if (v<0) v = 0;
     186        1428 :   if (n < 0) pari_err_DOMAIN("polhermite", "degree", "<", gen_0, stoi(n));
     187        1428 :   if (n==0) return pol_1(v);
     188             : 
     189        1421 :   q = cgetg(n+3, t_POL); r = q + n+2;
     190        1421 :   a = int2n(n);
     191        1421 :   gel(r--,0) = a;
     192        1421 :   gel(r--,0) = gen_0;
     193       40285 :   for (m=1; 2*m<= n; m++)
     194             :   {
     195       38864 :     av = avma;
     196       38864 :     a = diviuexact(muluui(n-2*m+2, n-2*m+1, a), 4*m);
     197       38864 :     togglesign(a);
     198       38864 :     gel(r--,0) = a = gerepileuptoint(av, a);
     199       38864 :     gel(r--,0) = gen_0;
     200             :   }
     201        1421 :   q[1] = evalsigne(1) | evalvarn(v);
     202        1421 :   return q;
     203             : }
     204             : GEN
     205        1442 : polhermite_eval(long n, GEN x)
     206             : {
     207             :   long i;
     208             :   pari_sp av, av2;
     209             :   GEN x2, u, v;
     210             : 
     211        1442 :   if (!x) return polhermite(n, 0);
     212          14 :   if (gequalX(x)) return polhermite(n, varn(x));
     213          14 :   if (n==0) return gen_1;
     214          14 :   if (n==1) return gmul2n(x,1);
     215          14 :   av = avma; x2 = gmul2n(x,1); v = gen_1; u = x2;
     216          14 :   av2= avma;
     217        7035 :   for (i=1; i<n; i++)
     218             :   { /* u = H_i(x), v = H_{i-1}(x), compute t = H_{i+1}(x) */
     219             :     GEN t;
     220        7021 :     if ((i & 0xff) == 0) gerepileall(av2,2,&u, &v);
     221        7021 :     t = gsub(gmul(x2, u), gmulsg(2*i,v));
     222        7021 :     v = u; u = t;
     223             :   }
     224          14 :   return gerepileupto(av, u);
     225             : }
     226             : 
     227             : /* Legendre polynomial
     228             :  * L0=1; L1=X; (n+1)*L(n+1)=(2*n+1)*X*L(n)-n*L(n-1)
     229             :  * L(n) = 2^-n sum_{k=0}^{n/2} a_k x^(n-2k)
     230             :  *   where a_k = (-1)^k (2n-2k)! / k! (n-k)! (n-2k)! is an integer
     231             :  *   and a_0 = binom(2n,n), a_k / a_{k-1} = - (n-2k+1)(n-2k+2) / 2k (2n-2k+1) */
     232             : GEN
     233        2184 : pollegendre(long n, long v)
     234             : {
     235             :   long k, l;
     236             :   pari_sp av;
     237             :   GEN a, r, q;
     238             : 
     239        2184 :   if (v<0) v = 0;
     240             :   /* pollegendre(-n) = pollegendre(n-1) */
     241        2184 :   if (n < 0) n = -n-1;
     242        2184 :   if (n==0) return pol_1(v);
     243        2142 :   if (n==1) return pol_x(v);
     244             : 
     245        2100 :   av = avma;
     246        2100 :   q = cgetg(n+3, t_POL); r = q + n+2;
     247        2100 :   gel(r--,0) = a = binomialuu(n<<1,n);
     248        2100 :   gel(r--,0) = gen_0;
     249       32354 :   for (k=1,l=n; l>1; k++,l-=2)
     250             :   { /* l = n-2*k+2 */
     251       30254 :     av = avma;
     252       30254 :     a = diviuuexact(muluui(l, l-1, a), 2*k, n+l-1);
     253       30254 :     togglesign(a); a = gerepileuptoint(av, a);
     254       30254 :     gel(r--,0) = a;
     255       30254 :     gel(r--,0) = gen_0;
     256             :   }
     257        2100 :   q[1] = evalsigne(1) | evalvarn(v);
     258        2100 :   return gerepileupto(av, gmul2n(q,-n));
     259             : }
     260             : 
     261             : GEN
     262        2163 : pollegendre_eval(long n, GEN x)
     263             : {
     264             :   long i;
     265             :   pari_sp av;
     266             :   GEN u, v;
     267             : 
     268        2163 :   if (!x) return pollegendre(n, 0);
     269          14 :   if (gequalX(x)) return pollegendre(n, varn(x));
     270             :   /* pollegendre(-n) = pollegendre(n-1) */
     271          14 :   if (n < 0) n = -n-1;
     272          14 :   if (n==0) return gen_1;
     273          14 :   if (n==1) return gcopy(x);
     274          14 :   av = avma; v = gen_1; u = x;
     275        7035 :   for (i=1; i<n; i++)
     276             :   { /* u = P_i(x), v = P_{i-1}(x), compute t = P_{i+1}(x) */
     277             :     GEN t;
     278        7021 :     if ((i & 0xff) == 0) gerepileall(av,2,&u, &v);
     279        7021 :     t = gdivgs(gsub(gmul(gmulsg(2*i+1,x), u), gmulsg(i,v)), i+1);
     280        7021 :     v = u; u = t;
     281             :   }
     282          14 :   return gerepileupto(av, u);
     283             : }
     284             : 
     285             : /* polcyclo(p) = X^(p-1) + ... + 1 */
     286             : static GEN
     287       35290 : polcyclo_prime(long p, long v)
     288             : {
     289       35290 :   GEN T = cgetg(p+2, t_POL);
     290             :   long i;
     291       35290 :   T[1] = evalsigne(1) | evalvarn(v);
     292       35290 :   for (i = 2; i < p+2; i++) gel(T,i) = gen_1;
     293       35290 :   return T;
     294             : }
     295             : 
     296             : /* cyclotomic polynomial */
     297             : GEN
     298       40037 : polcyclo(long n, long v)
     299             : {
     300             :   long s, q, i, l;
     301       40037 :   pari_sp av=avma;
     302             :   GEN T, P;
     303             : 
     304       40037 :   if (v<0) v = 0;
     305       40037 :   if (n < 3)
     306        4747 :     switch(n)
     307             :     {
     308         868 :       case 1: return deg1pol_shallow(gen_1, gen_m1, v);
     309        3879 :       case 2: return deg1pol_shallow(gen_1, gen_1, v);
     310           0 :       default: pari_err_DOMAIN("polcyclo", "index", "<=", gen_0, stoi(n));
     311             :     }
     312       35290 :   P = gel(factoru(n), 1); l = lg(P);
     313       35290 :   s = P[1]; T = polcyclo_prime(s, v);
     314       38153 :   for (i = 2; i < l; i++)
     315             :   { /* Phi_{np}(X) = Phi_n(X^p) / Phi_n(X) */
     316        2863 :     s *= P[i];
     317        2863 :     T = RgX_div(RgX_inflate(T, P[i]), T);
     318             :   }
     319             :   /* s = squarefree part of n */
     320       35290 :   q = n / s;
     321       35290 :   if (q == 1) return gerepileupto(av, T);
     322        5125 :   return gerepilecopy(av, RgX_inflate(T,q));
     323             : }
     324             : 
     325             : /* cyclotomic polynomial */
     326             : GEN
     327       14385 : polcyclo_eval(long n, GEN x)
     328             : {
     329       14385 :   pari_sp av= avma;
     330             :   GEN P, md, xd, yneg, ypos;
     331             :   long l, s, i, j, q, tx;
     332       14385 :   long root_of_1 = 0;
     333             : 
     334       14385 :   if (!x) return polcyclo(n, 0);
     335       12943 :   tx = typ(x);
     336       12943 :   if (gequalX(x)) return polcyclo(n, varn(x));
     337       12719 :   if (n <= 0) pari_err_DOMAIN("polcyclo", "index", "<=", gen_0, stoi(n));
     338       12719 :   if (n == 1) return gsubgs(x, 1);
     339       12719 :   if (tx == t_INT && !signe(x)) return gen_1;
     340       12719 :   while ((n & 3) == 0) { n >>= 1; x = gsqr(x); } /* Phi_4n(x) = Phi_2n(x^2) */
     341             :   /* n not divisible by 4 */
     342       12719 :   if (n == 2) return gerepileupto(av, gaddgs(x,1));
     343         819 :   if (!odd(n)) { n >>= 1; x = gneg(x); } /* Phi_2n(x) = Phi_n(-x) for n>1 odd */
     344             :   /* n odd > 2.  s largest squarefree divisor of n */
     345         819 :   P = gel(factoru(n), 1); s = zv_prod(P);
     346             :   /* replace n by largest squarefree divisor */
     347         819 :   q = n/s; if (q != 1) { x = gpowgs(x, q); n = s; }
     348         819 :   l = lg(P)-1;
     349             :   /* n squarefree odd > 2, l distinct prime divisors. Now handle x = 1 or -1 */
     350         819 :   if (tx == t_INT) { /* shortcut */
     351         749 :     if (is_pm1(x))
     352             :     {
     353          56 :       avma = av;
     354          56 :       if (signe(x) > 0 && l == 1) return utoipos(P[1]);
     355          35 :       return gen_1;
     356             :     }
     357             :   } else {
     358          70 :     if (gequal1(x))
     359             :     { /* n is prime, return n; multiply by x to keep the type */
     360          14 :       if (l == 1) return gerepileupto(av, gmulgs(x,n));
     361           7 :       return gerepilecopy(av, x); /* else 1 */
     362             :     }
     363          56 :     if (gequalm1(x)) return gerepileupto(av, gneg(x)); /* -1 */
     364             :   }
     365             :   /* Heuristic: evaluation will probably not improve things */
     366         742 :   if (tx == t_POL || tx == t_MAT || lg(x) > n)
     367          17 :     return gerepileupto(av, poleval(polcyclo(n,0), x));
     368             : 
     369         725 :   xd = cgetg((1L<<l) + 1, t_VEC); /* the x^d, where d | n */
     370         725 :   md = cgetg((1L<<l) + 1, t_VECSMALL); /* the mu(d), where d | n */
     371         725 :   gel(xd, 1) = x;
     372         725 :   md[1] = 1;
     373             :   /* Use Phi_n(x) = Prod_{d|n} (x^d-1)^mu(n/d).
     374             :    * If x has exact order D, n = Dq, then the result is 0 if q = 1. Otherwise
     375             :    * the factors with x^d-1, D|d are omitted and we multiply at the end by
     376             :    *   prod_{d | q} d^mu(q/d) = q if prime, 1 otherwise */
     377             :   /* We store the factors with mu(d)= 1 (resp.-1) in ypos (resp yneg).
     378             :    * At the end we return ypos/yneg if mu(n)=1 and yneg/ypos if mu(n)=-1 */
     379         725 :   ypos = gsubgs(x,1);
     380         725 :   yneg = gen_1;
     381        1520 :   for (i = 1; i <= l; i++)
     382             :   {
     383         795 :     long ti = 1L<<(i-1), p = P[i];
     384        1674 :     for (j = 1; j <= ti; j++) {
     385         879 :       GEN X = gpowgs(gel(xd,j), p), t = gsubgs(X,1);
     386         879 :       gel(xd,ti+j) = X;
     387         879 :       md[ti+j] = -md[j];
     388         879 :       if (gequal0(t))
     389             :       { /* x^d = 1; root_of_1 := the smallest index ti+j such that X == 1
     390             :         * (whose bits code d: bit i-1 is set iff P[i] | d). If no such index
     391             :         * exists, then root_of_1 remains 0. Do not multiply with X-1 if X = 1,
     392             :         * we handle these factors at the end */
     393          28 :         if (!root_of_1) root_of_1 = ti+j;
     394             :       }
     395             :       else
     396             :       {
     397         851 :         if (md[ti+j] == 1) ypos = gmul(ypos, t);
     398         788 :         else               yneg = gmul(yneg, t);
     399             :       }
     400             :     }
     401             :   }
     402         725 :   ypos = odd(l)? gdiv(yneg,ypos): gdiv(ypos,yneg);
     403         725 :   if (root_of_1)
     404             :   {
     405          21 :     GEN X = gel(xd,(1<<l)); /* = x^n = 1 */
     406          21 :     long bitmask_q = (1<<l) - root_of_1;
     407             :     /* bitmask_q encodes q = n/d: bit (i-1) is 1 iff P[i] | q */
     408             : 
     409             :     /* x is a root of unity.  If bitmask_q = 0, then x was a primitive n-th
     410             :      * root of 1 and the result is zero. Return X - 1 to preserve type. */
     411          21 :     if (!bitmask_q) return gerepileupto(av, gsubgs(X, 1));
     412             :     /* x is a primitive d-th root of unity, where d|n and d<n: we
     413             :      * must multiply ypos by if(isprime(n/d), n/d, 1) */
     414           7 :     ypos = gmul(ypos, X); /* multiply by X = 1 to preserve type */
     415             :     /* If bitmask_q = 1<<(i-1) for some i <= l, then q == P[i] and we multiply
     416             :      * by P[i]; otherwise q is composite and nothing more needs to be done */
     417           7 :     if (!(bitmask_q & (bitmask_q-1))) /* detects power of 2, since bitmask!=0 */
     418             :     {
     419           7 :       i = vals(bitmask_q)+1; /* q = P[i] */
     420           7 :       ypos = gmulgs(ypos, P[i]);
     421             :     }
     422             :   }
     423         711 :   return gerepileupto(av, ypos);
     424             : }
     425             : /********************************************************************/
     426             : /**                                                                **/
     427             : /**                  HILBERT & PASCAL MATRICES                     **/
     428             : /**                                                                **/
     429             : /********************************************************************/
     430             : GEN
     431         133 : mathilbert(long n) /* Hilbert matrix of order n */
     432             : {
     433             :   long i,j;
     434             :   GEN p;
     435             : 
     436         133 :   if (n < 0) pari_err_DOMAIN("mathilbert", "dimension", "<", gen_0, stoi(n));
     437         133 :   p = cgetg(n+1,t_MAT);
     438        1120 :   for (j=1; j<=n; j++)
     439             :   {
     440         987 :     gel(p,j) = cgetg(n+1,t_COL);
     441       16583 :     for (i=1+(j==1); i<=n; i++)
     442       15596 :       gcoeff(p,i,j) = mkfrac(gen_1, utoipos(i+j-1));
     443             :   }
     444         133 :   if (n) gcoeff(p,1,1) = gen_1;
     445         133 :   return p;
     446             : }
     447             : 
     448             : /* q-Pascal triangle = (choose(i,j)_q) (ordinary binomial if q = NULL) */
     449             : GEN
     450        1057 : matqpascal(long n, GEN q)
     451             : {
     452             :   long i, j, I;
     453        1057 :   pari_sp av = avma;
     454        1057 :   GEN m, qpow = NULL; /* gcc -Wall */
     455             : 
     456        1057 :   if (n < -1)  pari_err_DOMAIN("matpascal", "n", "<", gen_m1, stoi(n));
     457        1057 :   n++; m = cgetg(n+1,t_MAT);
     458        1057 :   for (j=1; j<=n; j++) gel(m,j) = cgetg(n+1,t_COL);
     459        1057 :   if (q)
     460             :   {
     461          42 :     I = (n+1)/2;
     462          42 :     if (I > 1) { qpow = new_chunk(I+1); gel(qpow,2)=q; }
     463          42 :     for (j=3; j<=I; j++) gel(qpow,j) = gmul(q, gel(qpow,j-1));
     464             :   }
     465       22120 :   for (i=1; i<=n; i++)
     466             :   {
     467       21063 :     I = (i+1)/2; gcoeff(m,i,1)= gen_1;
     468       21063 :     if (q)
     469             :     {
     470         483 :       for (j=2; j<=I; j++)
     471         476 :         gcoeff(m,i,j) = gadd(gmul(gel(qpow,j),gcoeff(m,i-1,j)),
     472         238 :                              gcoeff(m,i-1,j-1));
     473             :     }
     474             :     else
     475             :     {
     476     1023463 :       for (j=2; j<=I; j++)
     477     1002645 :         gcoeff(m,i,j) = addii(gcoeff(m,i-1,j), gcoeff(m,i-1,j-1));
     478             :     }
     479       21063 :     for (   ; j<=i; j++) gcoeff(m,i,j) = gcoeff(m,i,i+1-j);
     480       21063 :     for (   ; j<=n; j++) gcoeff(m,i,j) = gen_0;
     481             :   }
     482        1057 :   return gerepilecopy(av, m);
     483             : }
     484             : 
     485             : /******************************************************************/
     486             : /**                                                              **/
     487             : /**                       PRECISION CHANGES                      **/
     488             : /**                                                              **/
     489             : /******************************************************************/
     490             : 
     491             : GEN
     492         217 : gprec(GEN x, long l)
     493             : {
     494             :   long lx, i;
     495             :   GEN y;
     496             : 
     497         217 :   if (l <= 0) pari_err_DOMAIN("gprec", "precision", "<=", gen_0, stoi(l));
     498         217 :   switch(typ(x))
     499             :   {
     500             :     case t_REAL:
     501          21 :       return rtor(x, ndec2prec(l));
     502             :     case t_COMPLEX:
     503           7 :       y = cgetg(3, t_COMPLEX);
     504           7 :       gel(y,1) = gprec(gel(x,1),l);
     505           7 :       gel(y,2) = gprec(gel(x,2),l);
     506           7 :       break;
     507             :     case t_POL: case t_SER:
     508          42 :       y = cgetg_copy(x, &lx); y[1] = x[1];
     509          42 :       for (i=2; i<lx; i++) gel(y,i) = gprec(gel(x,i),l);
     510          42 :       break;
     511             :     case t_POLMOD: case t_RFRAC: case t_VEC: case t_COL: case t_MAT:
     512          21 :       y = cgetg_copy(x, &lx);
     513          21 :       for (i=1; i<lx; i++) gel(y,i) = gprec(gel(x,i),l);
     514          21 :       break;
     515         126 :     default: y = gcopy(x);
     516             :   }
     517         196 :   return y;
     518             : }
     519             : 
     520             : /* internal: precision given in word length (including codewords) */
     521             : GEN
     522     2100731 : gprec_w(GEN x, long pr)
     523             : {
     524             :   long lx, i;
     525             :   GEN y;
     526             : 
     527     2100731 :   switch(typ(x))
     528             :   {
     529             :     case t_REAL:
     530     1761213 :       if (signe(x)) return rtor(x,pr);
     531       11986 :       i = -prec2nbits(pr);
     532       11986 :       if (i < expo(x)) return real_0_bit(i);
     533       11479 :       y = cgetr(2); y[1] = x[1]; return y;
     534             :     case t_COMPLEX:
     535      167648 :       y = cgetg(3, t_COMPLEX);
     536      167648 :       gel(y,1) = gprec_w(gel(x,1),pr);
     537      167648 :       gel(y,2) = gprec_w(gel(x,2),pr);
     538      167648 :       break;
     539             :    case t_POL: case t_SER:
     540       45723 :       y = cgetg_copy(x, &lx); y[1] = x[1];
     541       45723 :       for (i=2; i<lx; i++) gel(y,i) = gprec_w(gel(x,i),pr);
     542       45723 :       break;
     543             :     case t_POLMOD: case t_RFRAC: case t_VEC: case t_COL: case t_MAT:
     544       26516 :       y = cgetg_copy(x, &lx);
     545       26516 :       for (i=1; i<lx; i++) gel(y,i) = gprec_w(gel(x,i),pr);
     546       26516 :       break;
     547       99631 :     default: return x;
     548             :   }
     549      239887 :   return y;
     550             : }
     551             : 
     552             : /* internal: precision given in word length (including codewords), truncate
     553             :  * mantissa to precision 'pr' but never _increase_ it */
     554             : GEN
     555       22302 : gprec_wtrunc(GEN x, long pr)
     556             : {
     557             :   long lx, i;
     558             :   GEN y;
     559             : 
     560       22302 :   switch(typ(x))
     561             :   {
     562             :     case t_REAL:
     563       17073 :       return (signe(x) && realprec(x) > pr)? rtor(x,pr): x;
     564             :     case t_COMPLEX:
     565        5208 :       y = cgetg(3, t_COMPLEX);
     566        5208 :       gel(y,1) = gprec_wtrunc(gel(x,1),pr);
     567        5208 :       gel(y,2) = gprec_wtrunc(gel(x,2),pr);
     568        5208 :       break;
     569             :     case t_POL:
     570             :     case t_SER:
     571           0 :       y = cgetg_copy(x, &lx); y[1] = x[1];
     572           0 :       for (i=2; i<lx; i++) gel(y,i) = gprec_wtrunc(gel(x,i),pr);
     573           0 :       break;
     574             :     case t_POLMOD: case t_RFRAC: case t_VEC: case t_COL: case t_MAT:
     575           7 :       y = cgetg_copy(x, &lx);
     576           7 :       for (i=1; i<lx; i++) gel(y,i) = gprec_wtrunc(gel(x,i),pr);
     577           7 :       break;
     578          14 :     default: return x;
     579             :   }
     580        5215 :   return y;
     581             : }
     582             : 
     583             : /********************************************************************/
     584             : /**                                                                **/
     585             : /**                      SERIES TRANSFORMS                         **/
     586             : /**                                                                **/
     587             : /********************************************************************/
     588             : /**                  LAPLACE TRANSFORM (OF A SERIES)               **/
     589             : /********************************************************************/
     590             : static GEN
     591           7 : serlaplace(GEN x)
     592             : {
     593           7 :   long i, l = lg(x), e = valp(x);
     594           7 :   GEN t, y = cgetg(l,t_SER);
     595           7 :   if (e < 0) pari_err_DOMAIN("laplace","valuation","<",gen_0,stoi(e));
     596           7 :   t = mpfact(e); y[1] = x[1];
     597         119 :   for (i=2; i<l; i++)
     598             :   {
     599         112 :     gel(y,i) = gmul(t, gel(x,i));
     600         112 :     e++; t = mului(e,t);
     601             :   }
     602           7 :   return y;
     603             : }
     604             : static GEN
     605          14 : pollaplace(GEN x)
     606             : {
     607          14 :   long i, e = 0, l = lg(x);
     608          14 :   GEN t = gen_1, y = cgetg(l,t_POL);
     609          14 :   y[1] = x[1];
     610          63 :   for (i=2; i<l; i++)
     611             :   {
     612          49 :     gel(y,i) = gmul(t, gel(x,i));
     613          49 :     e++; t = mului(e,t);
     614             :   }
     615          14 :   return y;
     616             : }
     617             : GEN
     618          21 : laplace(GEN x)
     619             : {
     620          21 :   pari_sp av = avma;
     621          21 :   switch(typ(x))
     622             :   {
     623          14 :     case t_POL: x = pollaplace(x); break;
     624           7 :     case t_SER: x = serlaplace(x); break;
     625           0 :     default: pari_err_TYPE("laplace",x);
     626             :   }
     627          21 :   return gerepilecopy(av, x);
     628             : }
     629             : 
     630             : /********************************************************************/
     631             : /**              CONVOLUTION PRODUCT (OF TWO SERIES)               **/
     632             : /********************************************************************/
     633             : GEN
     634           7 : convol(GEN x, GEN y)
     635             : {
     636           7 :   long j, lx, ly, ex, ey, vx = varn(x);
     637             :   GEN z;
     638             : 
     639           7 :   if (typ(x) != t_SER) pari_err_TYPE("convol",x);
     640           7 :   if (typ(y) != t_SER) pari_err_TYPE("convol",y);
     641           7 :   if (varn(y) != vx) pari_err_VAR("convol", x,y);
     642           7 :   ex = valp(x);
     643           7 :   ey = valp(y);
     644           7 :   if (ser_isexactzero(x))
     645           0 :     return scalarser(gadd(RgX_get_0(x), RgX_get_0(y)), varn(x), maxss(ex,ey));
     646           7 :   lx = lg(x) + ex; x -= ex;
     647           7 :   ly = lg(y) + ey; y -= ey;
     648             :   /* inputs shifted: x[i] and y[i] now correspond to monomials of same degree */
     649           7 :   if (ly < lx) lx = ly; /* min length */
     650           7 :   if (ex < ey) ex = ey; /* max valuation */
     651           7 :   if (lx - ex < 3) return zeroser(vx, lx-2);
     652             : 
     653           7 :   z = cgetg(lx - ex, t_SER);
     654           7 :   z[1] = evalvalp(ex) | evalvarn(vx);
     655           7 :   for (j = ex+2; j<lx; j++) gel(z,j-ex) = gmul(gel(x,j),gel(y,j));
     656           7 :   return normalize(z);
     657             : }
     658             : 
     659             : /***********************************************************************/
     660             : /*               OPERATIONS ON DIRICHLET SERIES: *, /                  */
     661             : /* (+, -, scalar multiplication are done on the corresponding vectors) */
     662             : /***********************************************************************/
     663             : static long
     664        1302 : dirval(GEN x)
     665             : {
     666        1302 :   long i = 1, lx = lg(x);
     667        1302 :   while (i < lx && gequal0(gel(x,i))) i++;
     668        1302 :   return i;
     669             : }
     670             : 
     671             : GEN
     672         133 : dirmul(GEN x, GEN y)
     673             : {
     674         133 :   pari_sp av = avma, av2;
     675             :   long nx, ny, nz, dx, dy, i, j, k;
     676             :   GEN z;
     677             : 
     678         133 :   if (typ(x)!=t_VEC) pari_err_TYPE("dirmul",x);
     679         133 :   if (typ(y)!=t_VEC) pari_err_TYPE("dirmul",y);
     680         133 :   dx = dirval(x); nx = lg(x)-1;
     681         133 :   dy = dirval(y); ny = lg(y)-1;
     682         133 :   if (ny-dy < nx-dx) { swap(x,y); lswap(nx,ny); lswap(dx,dy); }
     683         133 :   nz = minss(nx*dy,ny*dx);
     684         133 :   y = RgV_kill0(y);
     685         133 :   av2 = avma;
     686         133 :   z = zerovec(nz);
     687        6454 :   for (j=dx; j<=nx; j++)
     688             :   {
     689        6321 :     GEN c = gel(x,j);
     690        6321 :     if (gequal0(c)) continue;
     691        3101 :     if (gequal1(c))
     692             :     {
     693       13027 :       for (k=dy,i=j*dy; i<=nz; i+=j,k++)
     694       11522 :         if (gel(y,k)) gel(z,i) = gadd(gel(z,i),gel(y,k));
     695             :     }
     696        1596 :     else if (gequalm1(c))
     697             :     {
     698        4669 :       for (k=dy,i=j*dy; i<=nz; i+=j,k++)
     699        3563 :         if (gel(y,k)) gel(z,i) = gsub(gel(z,i),gel(y,k));
     700             :     }
     701             :     else
     702             :     {
     703        2233 :       for (k=dy,i=j*dy; i<=nz; i+=j,k++)
     704        1743 :         if (gel(y,k)) gel(z,i) = gadd(gel(z,i),gmul(c,gel(y,k)));
     705             :     }
     706        3101 :     if (gc_needed(av2,3))
     707             :     {
     708           0 :       if (DEBUGMEM>1) pari_warn(warnmem,"dirmul, %ld/%ld",j,nx);
     709           0 :       z = gerepilecopy(av2,z);
     710             :     }
     711             :   }
     712         133 :   return gerepilecopy(av,z);
     713             : }
     714             : 
     715             : GEN
     716         518 : dirdiv(GEN x, GEN y)
     717             : {
     718         518 :   pari_sp av = avma, av2;
     719             :   long nx,ny,nz, dx,dy, i,j,k;
     720             :   GEN p1;
     721             : 
     722         518 :   if (typ(x)!=t_VEC) pari_err_TYPE("dirdiv",x);
     723         518 :   if (typ(y)!=t_VEC) pari_err_TYPE("dirdiv",y);
     724         518 :   dx = dirval(x); nx = lg(x)-1;
     725         518 :   dy = dirval(y); ny = lg(y)-1;
     726         518 :   if (dy != 1 || !ny) pari_err_INV("dirdiv",y);
     727         518 :   nz = minss(nx,ny*dx);
     728         518 :   p1 = gel(y,1);
     729         518 :   if (gequal1(p1)) p1 = NULL; else y = gdiv(y,p1);
     730         518 :   y = RgV_kill0(y);
     731         518 :   av2 = avma;
     732         518 :   x = p1 ? gdiv(x,p1): leafcopy(x);
     733         518 :   for (j=1; j<dx; j++) gel(x,j) = gen_0;
     734         518 :   setlg(x,nz+1);
     735      514486 :   for (j=dx; j<=nz; j++)
     736             :   {
     737      513968 :     GEN c = gel(x,j);
     738      513968 :     if (gequal0(c)) continue;
     739      237489 :     if (gequal1(c))
     740             :     {
     741     1446046 :       for (i=j+j,k=2; i<=nz; i+=j,k++)
     742     1344371 :         if (gel(y,k)) gel(x,i) = gsub(gel(x,i),gel(y,k));
     743             :     }
     744      135814 :     else if (gequalm1(c))
     745             :     {
     746     1321964 :       for (i=j+j,k=2; i<=nz; i+=j,k++)
     747     1217314 :         if (gel(y,k)) gel(x,i) = gadd(gel(x,i),gel(y,k));
     748             :     }
     749             :     else
     750             :     {
     751      132888 :       for (i=j+j,k=2; i<=nz; i+=j,k++)
     752      101724 :         if (gel(y,k)) gel(x,i) = gsub(gel(x,i),gmul(c,gel(y,k)));
     753             :     }
     754      237489 :     if (gc_needed(av2,3))
     755             :     {
     756           0 :       if (DEBUGMEM>1) pari_warn(warnmem,"dirdiv, %ld/%ld",j,nz);
     757           0 :       x = gerepilecopy(av2,x);
     758             :     }
     759             :   }
     760         518 :   return gerepilecopy(av,x);
     761             : }
     762             : 
     763             : /*******************************************************************/
     764             : /**                                                               **/
     765             : /**                       COMBINATORICS                           **/
     766             : /**                                                               **/
     767             : /*******************************************************************/
     768             : /**                      BINOMIAL COEFFICIENTS                    **/
     769             : /*******************************************************************/
     770             : GEN
     771       74032 : binomialuu(ulong n, ulong k)
     772             : {
     773       74032 :   pari_sp ltop = avma;
     774             :   GEN z;
     775       74032 :   if (k > n) return gen_0;
     776       74025 :   k = minuu(k,n-k);
     777       74025 :   if (!k) return gen_1;
     778       65079 :   if (k == 1) return utoipos(n);
     779       58989 :   z = diviiexact(mulu_interval(n-k+1, n), mulu_interval(2UL, k));
     780       58989 :   return gerepileuptoint(ltop,z);
     781             : }
     782             : 
     783             : GEN
     784       99967 : binomial(GEN n, long k)
     785             : {
     786             :   long i, prec;
     787             :   pari_sp av;
     788             :   GEN y;
     789             : 
     790       99967 :   if (k <= 1)
     791             :   {
     792       59339 :     if (is_noncalc_t(typ(n))) pari_err_TYPE("binomial",n);
     793       59339 :     if (k < 0) return gen_0;
     794       59339 :     if (k == 0) return gen_1;
     795       26019 :     return gcopy(n);
     796             :   }
     797       40628 :   av = avma;
     798       40628 :   if (typ(n) == t_INT)
     799             :   {
     800       40614 :     if (signe(n) > 0)
     801             :     {
     802       40614 :       GEN z = subis(n,k);
     803       40614 :       if (cmpis(z,k) < 0)
     804             :       {
     805         889 :         k = itos(z); avma = av;
     806         889 :         if (k <= 1)
     807             :         {
     808         336 :           if (k < 0) return gen_0;
     809         336 :           if (k == 0) return gen_1;
     810         336 :           return icopy(n);
     811             :         }
     812             :       }
     813             :     }
     814             :     /* k > 1 */
     815       40278 :     if (lgefint(n) == 3 && signe(n) > 0)
     816             :     {
     817       40271 :       y = binomialuu(itou(n),(ulong)k);
     818       40271 :       return gerepileupto(av, y);
     819             :     }
     820             :     else
     821             :     {
     822           7 :       y = cgetg(k+1,t_VEC);
     823           7 :       for (i=1; i<=k; i++) gel(y,i) = subis(n,i-1);
     824           7 :       y = ZV_prod(y);
     825             :     }
     826           7 :     y = diviiexact(y, mpfact(k));
     827           7 :     return gerepileuptoint(av, y);
     828             :   }
     829             : 
     830          14 :   prec = precision(n);
     831          14 :   if (prec && k > 200 + 0.8*prec2nbits(prec)) {
     832           7 :     GEN A = mpfactr(k, prec), B = ggamma(gsubgs(n,k-1), prec);
     833           7 :     return gerepileupto(av, gdiv(ggamma(gaddgs(n,1), prec), gmul(A,B)));
     834             :   }
     835             : 
     836           7 :   y = cgetg(k+1,t_VEC);
     837           7 :   for (i=1; i<=k; i++) gel(y,i) = gsubgs(n,i-1);
     838           7 :   return gerepileupto(av, gdiv(RgV_prod(y), mpfact(k)));
     839             : }
     840             : 
     841             : /* Assume n >= 0, return bin, bin[k+1] = binomial(n, k) */
     842             : GEN
     843       13979 : vecbinome(long n)
     844             : {
     845             :   long d, k;
     846             :   GEN C;
     847       13979 :   if (!n) return mkvec(gen_1);
     848       13692 :   C = cgetg(n+2, t_VEC) + 1; /* C[k] = binomial(n, k) */
     849       13692 :   gel(C,0) = gen_1;
     850       13692 :   gel(C,1) = utoipos(n); d = (n + 1) >> 1;
     851      103026 :   for (k=2; k <= d; k++)
     852             :   {
     853       89334 :     pari_sp av = avma;
     854       89334 :     gel(C,k) = gerepileuptoint(av, diviuexact(mului(n-k+1, gel(C,k-1)), k));
     855             :   }
     856       13692 :   for (   ; k <= n; k++) gel(C,k) = gel(C,n-k);
     857       13692 :   return C - 1;
     858             : }
     859             : 
     860             : /********************************************************************/
     861             : /**                  STIRLING NUMBERS                              **/
     862             : /********************************************************************/
     863             : /* Stirling number of the 2nd kind. The number of ways of partitioning
     864             :    a set of n elements into m non-empty subsets. */
     865             : GEN
     866        1694 : stirling2(ulong n, ulong m)
     867             : {
     868        1694 :   pari_sp av = avma;
     869             :   GEN s, bmk;
     870             :   ulong k;
     871        1694 :   if (n==0) return (m == 0)? gen_1: gen_0;
     872        1694 :   if (m > n || m == 0) return gen_0;
     873        1694 :   if (m==n) return gen_1;
     874             :   /* k = 0 */
     875        1694 :   bmk = gen_1; s  = powuu(m, n);
     876       20314 :   for (k = 1; k <= ((m-1)>>1); ++k)
     877             :   { /* bmk = binomial(m, k) */
     878             :     GEN c, kn, mkn;
     879       18620 :     bmk = diviuexact(mului(m-k+1, bmk), k);
     880       18620 :     kn  = powuu(k, n); mkn = powuu(m-k, n);
     881       18620 :     c = odd(m)? subii(mkn,kn): addii(mkn,kn);
     882       18620 :     c = mulii(bmk, c);
     883       18620 :     s = odd(k)? subii(s, c): addii(s, c);
     884       18620 :     if (gc_needed(av,2))
     885             :     {
     886           0 :       if(DEBUGMEM>1) pari_warn(warnmem,"stirling2");
     887           0 :       gerepileall(av, 2, &s, &bmk);
     888             :     }
     889             :   }
     890             :   /* k = m/2 */
     891        1694 :   if (!odd(m))
     892             :   {
     893             :     GEN c;
     894         805 :     bmk = diviuexact(mului(k+1, bmk), k);
     895         805 :     c = mulii(bmk, powuu(k,n));
     896         805 :     s = odd(k)? subii(s, c): addii(s, c);
     897             :   }
     898        1694 :   return gerepileuptoint(av, diviiexact(s, mpfact(m)));
     899             : }
     900             : 
     901             : /* Stirling number of the first kind. Up to the sign, the number of
     902             :    permutations of n symbols which have exactly m cycles. */
     903             : GEN
     904         154 : stirling1(ulong n, ulong m)
     905             : {
     906         154 :   pari_sp ltop=avma;
     907             :   ulong k;
     908             :   GEN s, t;
     909         154 :   if (n < m) return gen_0;
     910         154 :   else if (n==m) return gen_1;
     911             :   /* t = binomial(n-1+k, m-1) * binomial(2n-m, n-m-k) */
     912             :   /* k = n-m > 0 */
     913         154 :   t = binomialuu(2*n-m-1, m-1);
     914         154 :   s = mulii(t, stirling2(2*(n-m), n-m));
     915         154 :   if (odd(n-m)) togglesign(s);
     916        1547 :   for (k = n-m-1; k > 0; --k)
     917             :   {
     918             :     GEN c;
     919        1393 :     t = diviuuexact(muluui(n-m+k+1, n+k+1, t), n+k, n-m-k);
     920        1393 :     c = mulii(t, stirling2(n-m+k, k));
     921        1393 :     s = odd(k)? subii(s, c): addii(s, c);
     922        1393 :     if ((k & 0x1f) == 0) {
     923          21 :       t = gerepileuptoint(ltop, t);
     924          21 :       s = gerepileuptoint(avma, s);
     925             :     }
     926             :   }
     927         154 :   return gerepileuptoint(ltop, s);
     928             : }
     929             : 
     930             : GEN
     931         301 : stirling(long n, long m, long flag)
     932             : {
     933         301 :   if (n < 0) pari_err_DOMAIN("stirling", "n", "<", gen_0, stoi(n));
     934         301 :   if (m < 0) pari_err_DOMAIN("stirling", "m", "<", gen_0, stoi(m));
     935         301 :   switch (flag)
     936             :   {
     937         154 :     case 1: return stirling1((ulong)n,(ulong)m);
     938         147 :     case 2: return stirling2((ulong)n,(ulong)m);
     939           0 :     default: pari_err_FLAG("stirling");
     940             :   }
     941           0 :   return NULL; /*NOT REACHED*/
     942             : }
     943             : 
     944             : /***********************************************************************/
     945             : /**                          PERMUTATIONS                             **/
     946             : /***********************************************************************/
     947             : GEN
     948        9492 : Z_to_perm(long n, GEN x)
     949             : {
     950             :   pari_sp av;
     951             :   ulong i, r;
     952        9492 :   GEN v = cgetg(n+1, t_VECSMALL);
     953        9492 :   if (n==0) return v;
     954        9485 :   uel(v,n) = 1; av = avma;
     955        9485 :   if (signe(x) <= 0) x = modii(x, mpfact(n));
     956       38143 :   for (r=n-1; r>=1; r--)
     957             :   {
     958             :     ulong a;
     959       28658 :     x = diviu_rem(x, n+1-r, &a);
     960       87157 :     for (i=r+1; i<=(ulong)n; i++)
     961       58499 :       if (uel(v,i) > a) uel(v,i)++;
     962       28658 :     uel(v,r) = a+1;
     963             :   }
     964        9485 :   avma = av; return v;
     965             : }
     966             : GEN
     967        9492 : numtoperm(long n, GEN x)
     968             : {
     969             :   long i;
     970             :   GEN v;
     971             : 
     972        9492 :   if (n < 0) pari_err_DOMAIN("numtoperm", "n", "<", gen_0, stoi(n));
     973        9492 :   if (typ(x) != t_INT) pari_err_TYPE("numtoperm",x);
     974        9492 :   v = Z_to_perm(n, x); settyp(v, t_VEC);
     975        9492 :   for (i = 1; i <= n; i++) gel(v,i) = utoipos(uel(v,i));
     976        9492 :   return v;
     977             : }
     978             : 
     979             : /* destroys v */
     980             : static GEN
     981        6020 : perm_to_Z_inplace(GEN v)
     982             : {
     983        6020 :   long l = lg(v), i, r;
     984        6020 :   GEN x = gen_0;
     985       31430 :   for (i = 1; i < l; i++)
     986             :   {
     987       25417 :     long vi = v[i];
     988       25417 :     if (vi <= 0) return NULL;
     989       25410 :     x = i==1 ? utoi(vi-1): addiu(muliu(x,l-i), vi-1);
     990       67354 :     for (r = i+1; r < l; r++)
     991       41944 :       if (v[r] > vi) v[r]--;
     992             :   }
     993        6013 :   return x;
     994             : }
     995             : GEN
     996         840 : perm_to_Z(GEN v)
     997             : {
     998         840 :   pari_sp av = avma;
     999         840 :   GEN x = perm_to_Z_inplace(leafcopy(v));
    1000         840 :   if (!x) pari_err_TYPE("permtonum",v);
    1001         840 :   return gerepileuptoint(av, x);
    1002             : }
    1003             : GEN
    1004        6027 : permtonum(GEN p)
    1005             : {
    1006        6027 :   pari_sp av = avma;
    1007             :   GEN v, x;
    1008        6027 :   switch(typ(p))
    1009             :   {
    1010         840 :     case t_VECSMALL: return perm_to_Z(p);
    1011             :     case t_VEC: case t_COL:
    1012        5180 :       if (RgV_is_ZV(p)) { v = ZV_to_zv(p); break; }
    1013           7 :     default: pari_err_TYPE("permtonum",p); v = NULL;
    1014             :   }
    1015        5180 :   x = perm_to_Z_inplace(v);
    1016        5180 :   if (!x) pari_err_TYPE("permtonum",p);
    1017        5173 :   return gerepileuptoint(av, x);
    1018             : }
    1019             : 
    1020             : /*******************************************************************/
    1021             : /**                                                               **/
    1022             : /**                     RECIPROCAL POLYNOMIAL                     **/
    1023             : /**                                                               **/
    1024             : /*******************************************************************/
    1025             : /* return coefficients s.t x = x_0 X^n + ... + x_n */
    1026             : GEN
    1027       62118 : polrecip(GEN x)
    1028             : {
    1029       62118 :   if (typ(x) != t_POL) pari_err_TYPE("polrecip",x);
    1030       62118 :   return RgX_recip(x);
    1031             : }
    1032             : 
    1033             : /********************************************************************/
    1034             : /**                                                                **/
    1035             : /**                  POLYNOMIAL INTERPOLATION                      **/
    1036             : /**                                                                **/
    1037             : /********************************************************************/
    1038             : /* allow X = NULL for [1,...,n] */
    1039             : GEN
    1040          28 : RgV_polint(GEN X, GEN Y, long v)
    1041             : {
    1042          28 :   pari_sp av0 = avma, av;
    1043          28 :   GEN Q, P = NULL;
    1044          28 :   long i, l = lg(Y);
    1045          28 :   if (!X)
    1046             :   {
    1047          14 :     X = cgetg(l, t_VEC);
    1048          14 :     for (i=1; i<l; i++) gel(X,i) = utoipos(i);
    1049             :   }
    1050          28 :   Q = roots_to_pol(X, v); av = avma;
    1051         105 :   for (i=1; i<l; i++)
    1052             :   {
    1053             :     GEN inv, T, dP;
    1054          77 :     if (gequal0(gel(Y,i))) continue;
    1055          63 :     T = RgX_div_by_X_x(Q, gel(X,i), NULL);
    1056          63 :     inv = ginv(poleval(T,gel(X,i)));
    1057          63 :     dP = RgX_Rg_mul(T, gmul(gel(Y,i),inv));
    1058          63 :     P = P? RgX_add(P, dP): dP;
    1059          63 :     if (gc_needed(av,2))
    1060             :     {
    1061           0 :       if (DEBUGMEM>1) pari_warn(warnmem,"FpV_polint");
    1062           0 :       P = gerepileupto(av, P);
    1063             :     }
    1064             :   }
    1065          28 :   if (!P) { avma = av; return zeropol(v); }
    1066          21 :   return gerepileupto(av0, P);
    1067             : }
    1068             : /* X,Y are "spec" GEN vectors with n > 1 components ( at X[0], ... X[n-1] ) */
    1069             : GEN
    1070         417 : polint_i(GEN X, GEN Y, GEN x, long n, GEN *ptdy)
    1071             : {
    1072         417 :   long i, m, ns = 0;
    1073         417 :   pari_sp av = avma;
    1074         417 :   GEN y, c, d, dy = NULL; /* gcc -Wall */
    1075             : 
    1076         417 :   if (!X)
    1077             :   {
    1078           7 :     X = cgetg(n+1, t_VEC);
    1079           7 :     for (i=1; i<=n; i++) gel(X,i) = utoipos(i);
    1080           7 :     X++;
    1081             :   }
    1082         417 :   switch(typ(x)) {
    1083             :     case t_INT: case t_REAL: case t_FRAC: case t_COMPLEX: case t_QUAD:
    1084             :     {
    1085         417 :       GEN D = NULL;
    1086       16888 :       for (i=0; i<n; i++)
    1087             :       {
    1088       16478 :         GEN t = gsub(x,gel(X,i));
    1089       16478 :         switch(typ(t))
    1090             :         {
    1091             :           case t_INT: case t_REAL: case t_FRAC: case t_COMPLEX: case t_QUAD:
    1092       16471 :             t = gabs(t, DEFAULTPREC);
    1093       16471 :             if (!D || gcmp(t,D) < 0) { ns = i; D = t; }
    1094       16471 :             break;
    1095             :           default:
    1096           7 :             goto NODY;
    1097             :         }
    1098             :       }
    1099         410 :       break;
    1100             :       /* X[ns] is closest to x */
    1101             :     }
    1102             : NODY:
    1103             :     default:
    1104           7 :       if (ptdy) { *ptdy = gen_0; ptdy = NULL; }
    1105             :   }
    1106         417 :   c = new_chunk(n);
    1107         417 :   d = new_chunk(n); for (i=0; i<n; i++) gel(c,i) = gel(d,i) = gel(Y,i);
    1108         417 :   y = gel(d,ns--);
    1109             :   /* divided differences */
    1110       16485 :   for (m=1; m<n; m++)
    1111             :   {
    1112      455241 :     for (i=0; i<n-m; i++)
    1113             :     {
    1114      439173 :       GEN ho = gsub(gel(X,i),x), hp = gsub(gel(X,i+m),x), den = gsub(ho,hp);
    1115      439173 :       if (gequal0(den))
    1116             :       {
    1117           7 :         char *x1 = stack_sprintf("X[%ld]", i+1);
    1118           7 :         char *x2 = stack_sprintf("X[%ld]", i+m+1);
    1119           7 :         pari_err_DOMAIN("polinterpolate",x1,"=",strtoGENstr(x2), X);
    1120             :       }
    1121      439166 :       den = gdiv(gsub(gel(c,i+1),gel(d,i)), den);
    1122      439166 :       gel(c,i) = gmul(ho,den);
    1123      439166 :       gel(d,i) = gmul(hp,den);
    1124             :     }
    1125       16068 :     dy = (2*(ns+1) < n-m)? gel(c,ns+1): gel(d,ns--);
    1126       16068 :     y = gadd(y,dy);
    1127             :   }
    1128         410 :   if (!ptdy) return gerepileupto(av, y);
    1129          98 :   *ptdy = dy;
    1130          98 :   gerepileall(av, 2, &y, ptdy);
    1131          98 :   return y;
    1132             : }
    1133             : 
    1134             : GEN
    1135         368 : polint(GEN X, GEN Y, GEN x, GEN *ptdy)
    1136             : {
    1137         368 :   long lx = lg(X);
    1138             : 
    1139         368 :   if (! is_vec_t(typ(X))) pari_err_TYPE("polinterpolate",X);
    1140         368 :   if (Y)
    1141             :   {
    1142         347 :     if (! is_vec_t(typ(Y))) pari_err_TYPE("polinterpolate",Y);
    1143         347 :     if (lx != lg(Y)) pari_err_DIM("polinterpolate");
    1144             :   }
    1145             :   else
    1146             :   {
    1147          21 :     Y = X;
    1148          21 :     X = NULL;
    1149             :   }
    1150             : 
    1151         368 :   if (lx <= 2)
    1152             :   {
    1153          14 :     if (ptdy) *ptdy = gen_0;
    1154          14 :     if (lx == 1) return zeropol(0);
    1155           7 :     Y = gel(Y,1);
    1156           7 :     if (gvar(Y) == 0) pari_err_PRIORITY("polinterpolate", Y, "=", 0);
    1157           7 :     return scalarpol(Y, 0);
    1158             :   }
    1159         354 :   if (!x) return RgV_polint(X, Y, 0);
    1160         326 :   if (gequalX(x)) return RgV_polint(X, Y, varn(x));
    1161         326 :   return polint_i(X? X+1: NULL,Y+1,x,lx-1,ptdy);
    1162             : }
    1163             : 
    1164             : /********************************************************************/
    1165             : /**                                                                **/
    1166             : /**                       MODREVERSE                               **/
    1167             : /**                                                                **/
    1168             : /********************************************************************/
    1169             : static void
    1170           7 : err_reverse(GEN x, GEN T)
    1171             : {
    1172           7 :   pari_err_DOMAIN("modreverse","deg(minpoly(z))", "<", stoi(degpol(T)),
    1173             :                   mkpolmod(x,T));
    1174           0 : }
    1175             : 
    1176             : /* return y such that Mod(y, charpoly(Mod(a,T)) = Mod(a,T) */
    1177             : GEN
    1178          28 : RgXQ_reverse(GEN a, GEN T)
    1179             : {
    1180          28 :   pari_sp av = avma;
    1181          28 :   long n = degpol(T);
    1182             :   GEN y;
    1183             : 
    1184          28 :   if (n <= 1) {
    1185           7 :     if (n <= 0) return gcopy(a);
    1186           7 :     return gerepileupto(av, gneg(gdiv(gel(T,2), gel(T,3))));
    1187             :   }
    1188          21 :   if (typ(a) != t_POL || !signe(a)) err_reverse(a,T);
    1189          21 :   y = RgXV_to_RgM(RgXQ_powers(a,n-1,T), n);
    1190          21 :   y = RgM_solve(y, col_ei(n, 2));
    1191          21 :   if (!y) err_reverse(a,T);
    1192          14 :   return gerepilecopy(av, RgV_to_RgX(y, varn(T)));
    1193             : }
    1194             : GEN
    1195         560 : QXQ_reverse(GEN a, GEN T)
    1196             : {
    1197         560 :   pari_sp av = avma;
    1198         560 :   long n = degpol(T);
    1199             :   GEN y;
    1200             : 
    1201         560 :   if (n <= 1) {
    1202          14 :     if (n <= 0) return gcopy(a);
    1203          14 :     return gerepileupto(av, gneg(gdiv(gel(T,2), gel(T,3))));
    1204             :   }
    1205         546 :   if (typ(a) != t_POL || !signe(a)) err_reverse(a,T);
    1206         546 :   if (gequalX(a)) return gcopy(a);
    1207         546 :   y = RgXV_to_RgM(QXQ_powers(a,n-1,T), n);
    1208         546 :   y = RgM_solve(y, col_ei(n, 2));
    1209         546 :   if (!y) err_reverse(a,T);
    1210         546 :   return gerepilecopy(av, RgV_to_RgX(y, varn(T)));
    1211             : }
    1212             : 
    1213             : GEN
    1214          28 : modreverse(GEN x)
    1215             : {
    1216             :   long v, n;
    1217             :   GEN T, a, y;
    1218             : 
    1219          28 :   if (typ(x)!=t_POLMOD) pari_err_TYPE("modreverse",x);
    1220          28 :   T = gel(x,1); n = degpol(T); if (n <= 0) return gcopy(x);
    1221          21 :   a = gel(x,2);
    1222          21 :   v = varn(T);
    1223          21 :   y = cgetg(3,t_POLMOD);
    1224          21 :   gel(y,1) = (n==1)? gsub(pol_x(v), a): RgXQ_charpoly(a, T, v);
    1225          21 :   gel(y,2) = RgXQ_reverse(a, T); return y;
    1226             : }
    1227             : 
    1228             : /********************************************************************/
    1229             : /**                                                                **/
    1230             : /**                          MERGESORT                             **/
    1231             : /**                                                                **/
    1232             : /********************************************************************/
    1233             : static int
    1234    51147695 : cmp_small(GEN x, GEN y) {
    1235    51147695 :   long a = (long)x, b = (long)y;
    1236    51147695 :   return a>b? 1: (a<b? -1: 0);
    1237             : }
    1238             : 
    1239             : static int
    1240        1120 : veccmp(void *data, GEN x, GEN y)
    1241             : {
    1242        1120 :   GEN k = (GEN)data;
    1243        1120 :   long i, s, lk = lg(k), lx = minss(lg(x), lg(y));
    1244             : 
    1245        1120 :   if (!is_vec_t(typ(x))) pari_err_TYPE("lexicographic vecsort",x);
    1246        1120 :   if (!is_vec_t(typ(y))) pari_err_TYPE("lexicographic vecsort",y);
    1247        1197 :   for (i=1; i<lk; i++)
    1248             :   {
    1249        1148 :     long c = k[i];
    1250        1148 :     if (c >= lx)
    1251          14 :       pari_err_TYPE("lexicographic vecsort, index too large", stoi(c));
    1252        1134 :     s = lexcmp(gel(x,c), gel(y,c));
    1253        1134 :     if (s) return s;
    1254             :   }
    1255          49 :   return 0;
    1256             : }
    1257             : 
    1258             : /* return permutation sorting v[1..n], removing duplicates. Assume n > 0 */
    1259             : static GEN
    1260      973322 : gen_sortspec_uniq(GEN v, long n, void *E, int (*cmp)(void*,GEN,GEN))
    1261             : {
    1262             :   pari_sp av;
    1263             :   long NX, nx, ny, m, ix, iy, i;
    1264             :   GEN x, y, w, W;
    1265             :   int s;
    1266      973322 :   switch(n)
    1267             :   {
    1268        1043 :     case 1: return mkvecsmall(1);
    1269             :     case 2:
    1270      323953 :       s = cmp(E,gel(v,1),gel(v,2));
    1271      323953 :       if      (s < 0) return mkvecsmall2(1,2);
    1272      166761 :       else if (s > 0) return mkvecsmall2(2,1);
    1273        2100 :       return mkvecsmall(1);
    1274             :     case 3:
    1275      164101 :       s = cmp(E,gel(v,1),gel(v,2));
    1276      164101 :       if (s < 0) {
    1277       82509 :         s = cmp(E,gel(v,2),gel(v,3));
    1278       82509 :         if (s < 0) return mkvecsmall3(1,2,3);
    1279       53088 :         else if (s == 0) return mkvecsmall2(1,2);
    1280       52696 :         s = cmp(E,gel(v,1),gel(v,3));
    1281       52696 :         if      (s < 0) return mkvecsmall3(1,3,2);
    1282       26712 :         else if (s > 0) return mkvecsmall3(3,1,2);
    1283        1715 :         return mkvecsmall2(1,2);
    1284       81592 :       } else if (s > 0) {
    1285       79989 :         s = cmp(E,gel(v,1),gel(v,3));
    1286       79989 :         if (s < 0) return mkvecsmall3(2,1,3);
    1287       54866 :         else if (s == 0) return mkvecsmall2(2,1);
    1288       53788 :         s = cmp(E,gel(v,2),gel(v,3));
    1289       53788 :         if (s < 0) return mkvecsmall3(2,3,1);
    1290       28091 :         else if (s > 0) return mkvecsmall3(3,2,1);
    1291         630 :         return mkvecsmall2(2,1);
    1292             :       } else {
    1293        1603 :         s = cmp(E,gel(v,1),gel(v,3));
    1294        1603 :         if (s < 0) return mkvecsmall2(1,3);
    1295        1015 :         else if (s == 0) return mkvecsmall(1);
    1296         798 :         return mkvecsmall2(3,1);
    1297             :       }
    1298             :   }
    1299      484225 :   NX = nx = n>>1; ny = n-nx;
    1300      484225 :   av = avma;
    1301      484225 :   x = gen_sortspec_uniq(v,   nx,E,cmp); nx = lg(x)-1;
    1302      484225 :   y = gen_sortspec_uniq(v+NX,ny,E,cmp); ny = lg(y)-1;
    1303      484225 :   w = cgetg(n+1, t_VECSMALL);
    1304      484225 :   m = ix = iy = 1;
    1305     7710822 :   while (ix<=nx && iy<=ny)
    1306             :   {
    1307     6742372 :     s = cmp(E, gel(v,x[ix]), gel(v,y[iy]+NX));
    1308     6742372 :     if (s < 0)
    1309     2874354 :       w[m++] = x[ix++];
    1310     3868018 :     else if (s > 0)
    1311     2990659 :       w[m++] = y[iy++]+NX;
    1312             :     else {
    1313      877359 :       w[m++] = x[ix++];
    1314      877359 :       iy++;
    1315             :     }
    1316             :   }
    1317      484225 :   while (ix<=nx) w[m++] = x[ix++];
    1318      484225 :   while (iy<=ny) w[m++] = y[iy++]+NX;
    1319      484225 :   avma = av;
    1320      484225 :   W = cgetg(m, t_VECSMALL);
    1321      484225 :   for (i = 1; i < m; i++) W[i] = w[i];
    1322      484225 :   return W;
    1323             : }
    1324             : 
    1325             : /* return permutation sorting v[1..n]. Assume n > 0 */
    1326             : static GEN
    1327   113841834 : gen_sortspec(GEN v, long n, void *E, int (*cmp)(void*,GEN,GEN))
    1328             : {
    1329             :   long nx, ny, m, ix, iy;
    1330             :   GEN x, y, w;
    1331   113841834 :   switch(n)
    1332             :   {
    1333             :     case 1:
    1334     2425221 :       (void)cmp(E,gel(v,1),gel(v,1)); /* check for type error */
    1335     2425221 :       return mkvecsmall(1);
    1336             :     case 2:
    1337    93878364 :       return cmp(E,gel(v,1),gel(v,2)) <= 0? mkvecsmall2(1,2)
    1338    46939175 :                                           : mkvecsmall2(2,1);
    1339             :     case 3:
    1340    21978653 :       if (cmp(E,gel(v,1),gel(v,2)) <= 0) {
    1341    15739943 :         if (cmp(E,gel(v,2),gel(v,3)) <= 0) return mkvecsmall3(1,2,3);
    1342     9740128 :         return (cmp(E,gel(v,1),gel(v,3)) <= 0)? mkvecsmall3(1,3,2)
    1343     4870064 :                                               : mkvecsmall3(3,1,2);
    1344             :       } else {
    1345     6238710 :         if (cmp(E,gel(v,1),gel(v,3)) <= 0) return mkvecsmall3(2,1,3);
    1346     7720684 :         return (cmp(E,gel(v,2),gel(v,3)) <= 0)? mkvecsmall3(2,3,1)
    1347     3860342 :                                               : mkvecsmall3(3,2,1);
    1348             :       }
    1349             :   }
    1350    42498771 :   nx = n>>1; ny = n-nx;
    1351    42498771 :   w = cgetg(n+1,t_VECSMALL);
    1352    42498771 :   x = gen_sortspec(v,   nx,E,cmp);
    1353    42498757 :   y = gen_sortspec(v+nx,ny,E,cmp);
    1354    42498757 :   m = ix = iy = 1;
    1355   311504739 :   while (ix<=nx && iy<=ny)
    1356   226507225 :     if (cmp(E, gel(v,x[ix]), gel(v,y[iy]+nx))<=0)
    1357   125752625 :       w[m++] = x[ix++];
    1358             :     else
    1359   100754600 :       w[m++] = y[iy++]+nx;
    1360    42498757 :   while (ix<=nx) w[m++] = x[ix++];
    1361    42498757 :   while (iy<=ny) w[m++] = y[iy++]+nx;
    1362    42498757 :   avma = (pari_sp)w; return w;
    1363             : }
    1364             : 
    1365             : static void
    1366    26980508 : init_sort(GEN *x, long *tx, long *lx)
    1367             : {
    1368    26980508 :   *tx = typ(*x);
    1369    26980508 :   if (*tx == t_LIST)
    1370             :   {
    1371          35 :     if (list_typ(*x)!=t_LIST_RAW) pari_err_TYPE("sort",*x);
    1372          35 :     *x = list_data(*x);
    1373          35 :     *lx = *x? lg(*x): 1;
    1374             :   } else {
    1375    26980473 :     if (!is_matvec_t(*tx) && *tx != t_VECSMALL) pari_err_TYPE("gen_sort",*x);
    1376    26980473 :     *lx = lg(*x);
    1377             :   }
    1378    26980508 : }
    1379             : 
    1380             : /* (x o y)[1..lx-1], destroy y */
    1381             : INLINE GEN
    1382     1692520 : sort_extract(GEN x, GEN y, long tx, long lx)
    1383             : {
    1384             :   long i;
    1385     1692520 :   switch(tx)
    1386             :   {
    1387             :     case t_VECSMALL:
    1388          14 :       for (i=1; i<lx; i++) y[i] = x[y[i]];
    1389          14 :       break;
    1390             :     case t_LIST:
    1391           7 :       settyp(y,t_VEC);
    1392           7 :       for (i=1; i<lx; i++) gel(y,i) = gel(x,y[i]);
    1393           7 :       return gtolist(y);
    1394             :     default:
    1395     1692499 :       settyp(y,tx);
    1396     1692499 :       for (i=1; i<lx; i++) gel(y,i) = gcopy(gel(x,y[i]));
    1397             :   }
    1398     1692513 :   return y;
    1399             : }
    1400             : 
    1401             : /* Sort the vector x, using cmp to compare entries. */
    1402             : GEN
    1403        4914 : gen_sort_uniq(GEN x, void *E, int (*cmp)(void*,GEN,GEN))
    1404             : {
    1405             :   long tx, lx;
    1406             :   GEN y;
    1407             : 
    1408        4914 :   init_sort(&x, &tx, &lx);
    1409        4914 :   if (lx==1) return tx == t_LIST? listcreate(): cgetg(1, tx);
    1410        4746 :   y = gen_sortspec_uniq(x,lx-1,E,cmp);
    1411        4746 :   return sort_extract(x, y, tx, lg(y)); /* lg(y) <= lx */
    1412             : }
    1413             : /* Sort the vector x, using cmp to compare entries. */
    1414             : GEN
    1415     2411351 : gen_sort(GEN x, void *E, int (*cmp)(void*,GEN,GEN))
    1416             : {
    1417             :   long tx, lx;
    1418             :   GEN y;
    1419             : 
    1420     2411351 :   init_sort(&x, &tx, &lx);
    1421     2411351 :   if (lx==1) return tx == t_LIST? listcreate(): cgetg(1, tx);
    1422     1687788 :   y = gen_sortspec(x,lx-1,E,cmp);
    1423     1687774 :   return sort_extract(x, y, tx, lx);
    1424             : }
    1425             : /* indirect sort: return the permutation that would sort x */
    1426             : GEN
    1427         126 : gen_indexsort_uniq(GEN x, void *E, int (*cmp)(void*,GEN,GEN))
    1428             : {
    1429             :   long tx, lx;
    1430         126 :   init_sort(&x, &tx, &lx);
    1431         126 :   if (lx==1) return cgetg(1, t_VECSMALL);
    1432         126 :   return gen_sortspec_uniq(x,lx-1,E,cmp);
    1433             : }
    1434             : /* indirect sort: return the permutation that would sort x */
    1435             : GEN
    1436     1329290 : gen_indexsort(GEN x, void *E, int (*cmp)(void*,GEN,GEN))
    1437             : {
    1438             :   long tx, lx;
    1439     1329290 :   init_sort(&x, &tx, &lx);
    1440     1329290 :   if (lx==1) return cgetg(1, t_VECSMALL);
    1441     1245607 :   return gen_sortspec(x,lx-1,E,cmp);
    1442             : }
    1443             : 
    1444             : /* Sort the vector x in place, using cmp to compare entries */
    1445             : void
    1446    23234827 : gen_sort_inplace(GEN x, void *E, int (*cmp)(void*,GEN,GEN), GEN *perm)
    1447             : {
    1448             :   long tx, lx, i;
    1449    23234827 :   pari_sp av = avma;
    1450             :   GEN y;
    1451             : 
    1452    23234827 :   init_sort(&x, &tx, &lx);
    1453    23234827 :   if (lx<=2)
    1454             :   {
    1455      118517 :     if (perm) *perm = lx == 1? cgetg(1, t_VECSMALL): mkvecsmall(1);
    1456    23353344 :     return;
    1457             :   }
    1458    23116310 :   y = gen_sortspec(x,lx-1, E, cmp);
    1459    23116310 :   if (perm)
    1460             :   {
    1461        9296 :     GEN z = new_chunk(lx);
    1462        9296 :     for (i=1; i<lx; i++) gel(z,i) = gel(x,y[i]);
    1463        9296 :     for (i=1; i<lx; i++) gel(x,i) = gel(z,i);
    1464        9296 :     *perm = y;
    1465        9296 :     avma = (pari_sp)y;
    1466             :   } else {
    1467    23107014 :     for (i=1; i<lx; i++) gel(y,i) = gel(x,y[i]);
    1468    23107014 :     for (i=1; i<lx; i++) gel(x,i) = gel(y,i);
    1469    23107014 :     avma = av;
    1470             :   }
    1471             : }
    1472             : 
    1473             : static int
    1474        7889 : closurecmp(void *data, GEN x, GEN y)
    1475             : {
    1476        7889 :   pari_sp av = avma;
    1477        7889 :   GEN z = closure_callgen2((GEN)data, x,y);
    1478        7889 :   if (typ(z) != t_INT)
    1479           0 :     pari_err_TYPE("closurecmp, cmp. fun. must return an integer", z);
    1480        7889 :   avma = av; return signe(z);
    1481             : }
    1482             : 
    1483             : static void
    1484         126 : check_positive_entries(GEN k)
    1485             : {
    1486         126 :   long i, l = lg(k);
    1487         287 :   for (i=1; i<l; i++)
    1488         161 :     if (k[i] <= 0) pari_err_DOMAIN("sort_function", "index", "<", gen_0, stoi(k[i]));
    1489         126 : }
    1490             : 
    1491             : typedef int (*CMP_FUN)(void*,GEN,GEN);
    1492             : static CMP_FUN
    1493         742 : sort_function(void **E, GEN x, GEN k)
    1494             : {
    1495         742 :   int (*cmp)(GEN,GEN) = &lexcmp;
    1496         742 :   if (!k)
    1497             :   {
    1498         119 :     *E = (void*)((typ(x) == t_VECSMALL)? cmp_small: cmp);
    1499         119 :     return &cmp_nodata;
    1500             :   }
    1501         623 :   if (typ(x) == t_VECSMALL) pari_err_TYPE("sort_function", x);
    1502         623 :   switch(typ(k))
    1503             :   {
    1504          91 :     case t_INT: k = mkvecsmall(itos(k));  break;
    1505          35 :     case t_VEC: case t_COL: k = ZV_to_zv(k); break;
    1506           0 :     case t_VECSMALL: break;
    1507             :     case t_CLOSURE:
    1508         497 :      if (closure_arity(k) != 2 || closure_is_variadic(k))
    1509           0 :        pari_err_TYPE("sort_function, cmp. fun. needs exactly 2 arguments",k);
    1510         497 :      *E = (void*)k;
    1511         497 :      return &closurecmp;
    1512           0 :     default: pari_err_TYPE("sort_function",k);
    1513             :   }
    1514         126 :   check_positive_entries(k);
    1515         126 :   *E = (void*)k;
    1516         126 :   return &veccmp;
    1517             : }
    1518             : 
    1519             : #define cmp_IND 1
    1520             : #define cmp_LEX 2 /* FIXME: backward compatibility, ignored */
    1521             : #define cmp_REV 4
    1522             : #define cmp_UNIQ 8
    1523             : GEN
    1524         672 : vecsort0(GEN x, GEN k, long flag)
    1525             : {
    1526             :   void *E;
    1527         672 :   int (*CMP)(void*,GEN,GEN) = sort_function(&E, x, k);
    1528             : 
    1529         672 :   if (flag < 0 || flag > (cmp_REV|cmp_LEX|cmp_IND|cmp_UNIQ))
    1530           0 :     pari_err_FLAG("vecsort");
    1531         672 :   if (flag & cmp_UNIQ)
    1532          28 :     x = flag & cmp_IND? gen_indexsort_uniq(x, E, CMP): gen_sort_uniq(x, E, CMP);
    1533             :   else
    1534         644 :     x = flag & cmp_IND? gen_indexsort(x, E, CMP): gen_sort(x, E, CMP);
    1535         658 :   if (flag & cmp_REV) { /* reverse order */
    1536             :     long j, lx;
    1537             :     GEN y;
    1538          21 :     if (typ(x)==t_LIST)
    1539             :     {
    1540           7 :       y = list_data(x);
    1541           7 :       if (!y) return x;
    1542             :     }
    1543             :     else
    1544          14 :       y = x;
    1545          14 :     lx = lg(y);
    1546          14 :     for (j=1; j<=(lx-1)>>1; j++) swap(gel(y,j), gel(y,lx-j));
    1547             :   }
    1548         651 :   return x;
    1549             : }
    1550             : 
    1551             : GEN
    1552        5081 : indexsort(GEN x) { return gen_indexsort(x, (void*)&gcmp, cmp_nodata); }
    1553             : GEN
    1554           0 : indexlexsort(GEN x) { return gen_indexsort(x, (void*)&lexcmp, cmp_nodata); }
    1555             : GEN
    1556           0 : indexvecsort(GEN x, GEN k)
    1557             : {
    1558           0 :   if (typ(k) != t_VECSMALL) pari_err_TYPE("vecsort",k);
    1559           0 :   return gen_indexsort(x, (void*)k, &veccmp);
    1560             : }
    1561             : 
    1562             : GEN
    1563      589791 : sort(GEN x) { return gen_sort(x, (void*)gcmp, cmp_nodata); }
    1564             : GEN
    1565           0 : lexsort(GEN x) { return gen_sort(x, (void*)lexcmp, cmp_nodata); }
    1566             : GEN
    1567           0 : vecsort(GEN x, GEN k)
    1568             : {
    1569           0 :   if (typ(k) != t_VECSMALL) pari_err_TYPE("vecsort",k);
    1570           0 :   return gen_sort(x, (void*)k, &veccmp);
    1571             : }
    1572             : long
    1573          70 : vecsearch(GEN v, GEN x, GEN k)
    1574             : {
    1575          70 :   pari_sp av = avma;
    1576             :   void *E;
    1577          70 :   int (*CMP)(void*,GEN,GEN) = sort_function(&E, x, k);
    1578             :   long r;
    1579          70 :   if (!is_matvec_t(typ(v))) pari_err_TYPE("vecsearch", v);
    1580          70 :   r = gen_search(v, x, 0, E, CMP);
    1581          70 :   avma = av; return r;
    1582             : }
    1583             : 
    1584             : GEN
    1585         448 : ZV_indexsort(GEN L) { return gen_indexsort(L, (void*)&cmpii, &cmp_nodata); }
    1586             : GEN
    1587     1141189 : ZV_sort(GEN L) { return gen_sort(L, (void*)&cmpii, &cmp_nodata); }
    1588             : GEN
    1589        2184 : ZV_sort_uniq(GEN L) { return gen_sort_uniq(L, (void*)&cmpii, &cmp_nodata); }
    1590             : 
    1591             : /********************************************************************/
    1592             : /**                      SEARCH IN SORTED VECTOR                   **/
    1593             : /********************************************************************/
    1594             : /* index of x in table T, 0 otherwise */
    1595             : long
    1596    16224467 : tablesearch(GEN T, GEN x, int (*cmp)(GEN,GEN))
    1597             : {
    1598    16224467 :   long l = 1, u = lg(T)-1, i, s;
    1599             : 
    1600    86838303 :   while (u>=l)
    1601             :   {
    1602    57887953 :     i = (l+u)>>1; s = cmp(x, gel(T,i));
    1603    57887953 :     if (!s) return i;
    1604    54389369 :     if (s<0) u=i-1; else l=i+1;
    1605             :   }
    1606    12725883 :   return 0;
    1607             : }
    1608             : 
    1609             : /* looks if x belongs to the set T and returns the index if yes, 0 if no */
    1610             : long
    1611    14044674 : gen_search(GEN T, GEN x, long flag, void *data, int (*cmp)(void*,GEN,GEN))
    1612             : {
    1613    14044674 :   long lx = lg(T), i, l, u, s;
    1614             : 
    1615    14044674 :   if (lx==1) return flag? 1: 0;
    1616    14044674 :   l = 1; u = lx-1;
    1617             :   do
    1618             :   {
    1619    47536398 :     i = (l+u)>>1; s = cmp(data, x, gel(T,i));
    1620    47536398 :     if (!s) return flag? 0: i;
    1621    33570936 :     if (s<0) u=i-1; else l=i+1;
    1622    33570936 :   } while (u>=l);
    1623       79212 :   if (!flag) return 0;
    1624       41356 :   return (s<0)? i: i+1;
    1625             : }
    1626             : 
    1627             : long
    1628      894236 : ZV_search(GEN x, GEN y) { return tablesearch(x, y, cmpii); }
    1629             : 
    1630             : long
    1631    15327886 : zv_search(GEN x, long y) { return tablesearch(x, (GEN)y, cmp_small); }
    1632             : 
    1633             : /********************************************************************/
    1634             : /**                   COMPARISON FUNCTIONS                         **/
    1635             : /********************************************************************/
    1636             : int
    1637   379477994 : cmp_nodata(void *data, GEN x, GEN y)
    1638             : {
    1639   379477994 :   int (*cmp)(GEN,GEN)=(int (*)(GEN,GEN)) data;
    1640   379477994 :   return cmp(x,y);
    1641             : }
    1642             : 
    1643             : /* assume x and y come from the same idealprimedec call (uniformizer unique) */
    1644             : int
    1645      604397 : cmp_prime_over_p(GEN x, GEN y)
    1646             : {
    1647      604397 :   long k = pr_get_f(x) - pr_get_f(y); /* diff. between residue degree */
    1648      627066 :   return k? ((k > 0)? 1: -1)
    1649      627066 :           : ZV_cmp(pr_get_gen(x), pr_get_gen(y));
    1650             : }
    1651             : 
    1652             : int
    1653       38307 : cmp_prime_ideal(GEN x, GEN y)
    1654             : {
    1655       38307 :   int k = cmpii(pr_get_p(x), pr_get_p(y));
    1656       38307 :   return k? k: cmp_prime_over_p(x,y);
    1657             : }
    1658             : 
    1659             : /* assume x and y are t_POL in the same variable whose coeffs can be
    1660             :  * compared (used to sort polynomial factorizations) */
    1661             : int
    1662      940000 : gen_cmp_RgX(void *data, GEN x, GEN y)
    1663             : {
    1664      940000 :   int (*coeff_cmp)(GEN,GEN)=(int(*)(GEN,GEN))data;
    1665      940000 :   long i, lx = lg(x), ly = lg(y);
    1666             :   int fl;
    1667      940000 :   if (lx > ly) return  1;
    1668      874985 :   if (lx < ly) return -1;
    1669     2126019 :   for (i=lx-1; i>1; i--)
    1670     1934391 :     if ((fl = coeff_cmp(gel(x,i), gel(y,i)))) return fl;
    1671      191628 :   return 0;
    1672             : }
    1673             : 
    1674             : static int
    1675        2576 : cmp_RgX_Rg(GEN x, GEN y)
    1676             : {
    1677        2576 :   long lx = lgpol(x), ly;
    1678        2576 :   if (lx > 1) return  1;
    1679        2177 :   ly = gequal0(y) ? 0:1;
    1680        2177 :   if (lx > ly) return  1;
    1681        2177 :   if (lx < ly) return -1;
    1682        2177 :   if (lx==0) return 0;
    1683         427 :   return gcmp(gel(x,2), y);
    1684             : }
    1685             : int
    1686       90497 : cmp_RgX(GEN x, GEN y)
    1687             : {
    1688       90497 :   if (typ(x) == t_POLMOD) x = gel(x,2);
    1689       90497 :   if (typ(y) == t_POLMOD) y = gel(y,2);
    1690       90497 :   if (typ(x) == t_POL) {
    1691       43317 :     if (typ(y) != t_POL) return cmp_RgX_Rg(x, y);
    1692             :   } else {
    1693       47180 :     if (typ(y) != t_POL) return gcmp(x,y);
    1694        1225 :     return - cmp_RgX_Rg(y,x);
    1695             :   }
    1696       41966 :   return gen_cmp_RgX((void*)&gcmp,x,y);
    1697             : }
    1698             : 
    1699             : int
    1700      125368 : cmp_Flx(GEN x, GEN y)
    1701             : {
    1702      125368 :   long i, lx = lg(x), ly = lg(y);
    1703      125368 :   if (lx > ly) return  1;
    1704      114245 :   if (lx < ly) return -1;
    1705      243702 :   for (i=lx-1; i>1; i--)
    1706      242145 :     if (uel(x,i) != uel(y,i)) return uel(x,i)<uel(y,i)? -1: 1;
    1707        1557 :   return 0;
    1708             : }
    1709             : /********************************************************************/
    1710             : /**                   MERGE & SORT FACTORIZATIONS                  **/
    1711             : /********************************************************************/
    1712             : /* merge fx, fy two factorizations, whose 1st column is sorted in strictly
    1713             :  * increasing order wrt cmp */
    1714             : GEN
    1715      911795 : merge_factor(GEN fx, GEN fy, void *data, int (*cmp)(void *,GEN,GEN))
    1716             : {
    1717      911795 :   GEN x = gel(fx,1), e = gel(fx,2), M, E;
    1718      911795 :   GEN y = gel(fy,1), f = gel(fy,2);
    1719      911795 :   long ix, iy, m, lx = lg(x), ly = lg(y), l = lx+ly-1;
    1720             : 
    1721      911795 :   M = cgetg(l, t_COL);
    1722      911795 :   E = cgetg(l, t_COL);
    1723             : 
    1724      911795 :   m = ix = iy = 1;
    1725    11295925 :   while (ix<lx && iy<ly)
    1726             :   {
    1727     9472335 :     int s = cmp(data, gel(x,ix), gel(y,iy));
    1728     9472335 :     if (s < 0)
    1729     8612497 :     { gel(M,m) = gel(x,ix); gel(E,m) = gel(e,ix); ix++; }
    1730      859838 :     else if (s == 0)
    1731             :     {
    1732      162767 :       GEN z = gel(x,ix), g = addii(gel(e,ix), gel(f,iy));
    1733      162767 :       iy++; ix++; if (!signe(g)) continue;
    1734       85739 :       gel(M,m) = z; gel(E,m) = g;
    1735             :     }
    1736             :     else
    1737      697071 :     { gel(M,m) = gel(y,iy); gel(E,m) = gel(f,iy); iy++; }
    1738     9395307 :     m++;
    1739             :   }
    1740      911795 :   while (ix<lx) { gel(M,m) = gel(x,ix); gel(E,m) = gel(e,ix); ix++; m++; }
    1741      911795 :   while (iy<ly) { gel(M,m) = gel(y,iy); gel(E,m) = gel(f,iy); iy++; m++; }
    1742      911795 :   setlg(M, m);
    1743      911795 :   setlg(E, m); return mkmat2(M, E);
    1744             : }
    1745             : /* merge two sorted vectors, removing duplicates. Shallow */
    1746             : GEN
    1747       32403 : merge_sort_uniq(GEN x, GEN y, void *data, int (*cmp)(void *,GEN,GEN))
    1748             : {
    1749       32403 :   long ix, iy, m, lx = lg(x), ly = lg(y), l = lx+ly-1;
    1750             :   GEN M;
    1751             : 
    1752       32403 :   M = cgetg(l, t_COL);
    1753       32403 :   m = ix = iy = 1;
    1754       85659 :   while (ix<lx && iy<ly)
    1755             :   {
    1756       20853 :     int s = cmp(data, gel(x,ix), gel(y,iy));
    1757       20853 :     if (s < 0)
    1758       13321 :     { gel(M,m) = gel(x,ix); ix++; }
    1759        7532 :     else if (s == 0)
    1760           0 :     { gel(M,m) = gel(x,ix); iy++; ix++; }
    1761             :     else
    1762        7532 :     { gel(M,m) = gel(y,iy); iy++; }
    1763       20853 :     m++;
    1764             :   }
    1765       32403 :   while (ix<lx) { gel(M,m) = gel(x,ix); ix++; m++; }
    1766       32403 :   while (iy<ly) { gel(M,m) = gel(y,iy); iy++; m++; }
    1767       32403 :   setlg(M, m); return M;
    1768             : }
    1769             : 
    1770             : /* sort generic factorization, in place */
    1771             : GEN
    1772     2796911 : sort_factor(GEN y, void *data, int (*cmp)(void *,GEN,GEN))
    1773             : {
    1774             :   GEN a, b, A, B, w;
    1775             :   pari_sp av;
    1776             :   long n, i;
    1777             : 
    1778     2796911 :   a = gel(y,1); n = lg(a); if (n == 1) return y;
    1779     2794601 :   b = gel(y,2); av = avma;
    1780     2794601 :   A = new_chunk(n);
    1781     2794601 :   B = new_chunk(n);
    1782     2794601 :   w = gen_sortspec(a, n-1, data, cmp);
    1783     2794601 :   for (i=1; i<n; i++) { long k = w[i]; A[i] = a[k]; B[i] = b[k]; }
    1784     2794601 :   for (i=1; i<n; i++) { a[i] = A[i]; b[i] = B[i]; }
    1785     2794601 :   avma = av; return y;
    1786             : }
    1787             : /* sort polynomial factorization, in place */
    1788             : GEN
    1789      422541 : sort_factor_pol(GEN y,int (*cmp)(GEN,GEN))
    1790             : {
    1791      422541 :   (void)sort_factor(y,(void*)cmp, &gen_cmp_RgX);
    1792      422541 :   return y;
    1793             : }
    1794             : 
    1795             : /* assume f and g coprime integer factorizations */
    1796             : GEN
    1797         322 : merge_factor_i(GEN f, GEN g)
    1798             : {
    1799         322 :   if (lg(f) == 1) return g;
    1800         322 :   if (lg(g) == 1) return f;
    1801         322 :   return sort_factor(famat_mul_shallow(f,g), (void*)&cmpii, &cmp_nodata);
    1802             : }
    1803             : 
    1804             : /***********************************************************************/
    1805             : /*                                                                     */
    1806             : /*                          SET OPERATIONS                             */
    1807             : /*                                                                     */
    1808             : /***********************************************************************/
    1809             : GEN
    1810         896 : gtoset(GEN x)
    1811             : {
    1812             :   long lx;
    1813         896 :   if (!x) return cgetg(1, t_VEC);
    1814         896 :   switch(typ(x))
    1815             :   {
    1816             :     case t_VEC:
    1817         868 :     case t_COL: lx = lg(x); break;
    1818             :     case t_LIST:
    1819          14 :       if (list_typ(x)==t_LIST_MAP) return mapdomain(x);
    1820          14 :       x = list_data(x); lx = x? lg(x): 1; break;
    1821           7 :     case t_VECSMALL: lx = lg(x); x = zv_to_ZV(x); break;
    1822           7 :     default: return mkveccopy(x);
    1823             :   }
    1824         889 :   if (lx==1) return cgetg(1,t_VEC);
    1825         882 :   x = gen_sort_uniq(x, (void*)&cmp_universal, cmp_nodata);
    1826         882 :   settyp(x, t_VEC); /* it may be t_COL */
    1827         882 :   return x;
    1828             : }
    1829             : 
    1830             : long
    1831          14 : setisset(GEN x)
    1832             : {
    1833          14 :   long i, lx = lg(x);
    1834             : 
    1835          14 :   if (typ(x) != t_VEC) return 0;
    1836          14 :   if (lx == 1) return 1;
    1837          70 :   for (i=1; i<lx-1; i++)
    1838          63 :     if (cmp_universal(gel(x,i+1), gel(x,i)) <= 0) return 0;
    1839           7 :   return 1;
    1840             : }
    1841             : 
    1842             : long
    1843      105903 : setsearch(GEN T, GEN y, long flag)
    1844             : {
    1845             :   long lx;
    1846      105903 :   switch(typ(T))
    1847             :   {
    1848      105889 :     case t_VEC: lx = lg(T); break;
    1849             :     case t_LIST:
    1850           7 :     if (list_typ(T) != t_LIST_RAW) pari_err_TYPE("setsearch",T);
    1851           7 :     T = list_data(T); lx = T? lg(T): 1; break;
    1852           7 :     default: pari_err_TYPE("setsearch",T);
    1853           0 :       return 0; /*not reached*/
    1854             :   }
    1855      105896 :   if (lx==1) return flag? 1: 0;
    1856      105896 :   return gen_search(T,y,flag,(void*)cmp_universal,cmp_nodata);
    1857             : }
    1858             : 
    1859             : GEN
    1860          21 : setunion(GEN x, GEN y)
    1861             : {
    1862          21 :   pari_sp av = avma;
    1863          21 :   long i, j, k, lx = lg(x), ly = lg(y);
    1864          21 :   GEN z = cgetg(lx + ly - 1, t_VEC);
    1865          21 :   if (typ(x) != t_VEC) pari_err_TYPE("setunion",x);
    1866          21 :   if (typ(y) != t_VEC) pari_err_TYPE("setunion",y);
    1867          21 :   i = j = k = 1;
    1868         119 :   while (i<lx && j<ly)
    1869             :   {
    1870          77 :     int s = cmp_universal(gel(x,i), gel(y,j));
    1871          77 :     if (s < 0)
    1872          28 :       z[k++] = x[i++];
    1873          49 :     else if (s > 0)
    1874          14 :       z[k++] = y[j++];
    1875             :     else {
    1876          35 :       z[k++] = x[i++];
    1877          35 :       j++;
    1878             :     }
    1879             :   }
    1880          21 :   while (i<lx) z[k++] = x[i++];
    1881          21 :   while (j<ly) z[k++] = y[j++];
    1882          21 :   setlg(z, k);
    1883          21 :   return gerepilecopy(av, z);
    1884             : }
    1885             : /* in case of equal keys in x,y, take the key from x */
    1886             : GEN
    1887      295907 : ZV_union_shallow(GEN x, GEN y)
    1888             : {
    1889      295907 :   long i, j, k, lx = lg(x), ly = lg(y);
    1890      295907 :   GEN z = cgetg(lx + ly - 1, t_VEC);
    1891      295907 :   i = j = k = 1;
    1892     1127292 :   while (i<lx && j<ly)
    1893             :   {
    1894      535478 :     int s = cmpii(gel(x,i), gel(y,j));
    1895      535478 :     if (s < 0)
    1896        8097 :       z[k++] = x[i++];
    1897      527381 :     else if (s > 0)
    1898      196564 :       z[k++] = y[j++];
    1899             :     else {
    1900      330817 :       z[k++] = x[i++];
    1901      330817 :       j++;
    1902             :     }
    1903             :   }
    1904      295907 :   while (i<lx) z[k++] = x[i++];
    1905      295907 :   while (j<ly) z[k++] = y[j++];
    1906      295907 :   setlg(z, k); return z;
    1907             : }
    1908             : 
    1909             : 
    1910             : GEN
    1911           7 : setintersect(GEN x, GEN y)
    1912             : {
    1913           7 :   long ix = 1, iy = 1, iz = 1, lx = lg(x), ly = lg(y);
    1914           7 :   pari_sp av = avma;
    1915           7 :   GEN z = cgetg(lx,t_VEC);
    1916           7 :   if (typ(x) != t_VEC) pari_err_TYPE("setintersect",x);
    1917           7 :   if (typ(y) != t_VEC) pari_err_TYPE("setintersect",y);
    1918          77 :   while (ix < lx && iy < ly)
    1919             :   {
    1920          63 :     int c = cmp_universal(gel(x,ix), gel(y,iy));
    1921          63 :     if      (c < 0) ix++;
    1922          35 :     else if (c > 0) iy++;
    1923          21 :     else { gel(z, iz++) = gel(x,ix); ix++; iy++; }
    1924             :   }
    1925           7 :   setlg(z,iz); return gerepilecopy(av,z);
    1926             : }
    1927             : 
    1928             : GEN
    1929           7 : gen_setminus(GEN A, GEN B, int (*cmp)(GEN,GEN))
    1930             : {
    1931           7 :   pari_sp ltop = avma;
    1932           7 :   long i = 1, j = 1, k = 1, lx = lg(A), ly = lg(B);
    1933           7 :   GEN  diff = cgetg(lx,t_VEC);
    1934          91 :   while (i < lx && j < ly)
    1935          77 :     switch ( cmp(gel(A,i),gel(B,j)) )
    1936             :     {
    1937          28 :       case -1: gel(diff,k++) = gel(A,i++); break;
    1938          28 :       case 1: j++; break;
    1939          21 :       case 0: i++; break;
    1940             :     }
    1941           7 :   while (i < lx) gel(diff,k++) = gel(A,i++);
    1942           7 :   setlg(diff,k);
    1943           7 :   return gerepilecopy(ltop,diff);
    1944             : }
    1945             : 
    1946             : GEN
    1947           7 : setminus(GEN x, GEN y)
    1948             : {
    1949           7 :   if (typ(x) != t_VEC) pari_err_TYPE("setminus",x);
    1950           7 :   if (typ(y) != t_VEC) pari_err_TYPE("setminus",y);
    1951           7 :   return gen_setminus(x,y,cmp_universal);
    1952             : }
    1953             : 
    1954             : GEN
    1955          21 : setbinop(GEN f, GEN x, GEN y)
    1956             : {
    1957          21 :   pari_sp av = avma;
    1958          21 :   long i, j, lx, ly, k = 1;
    1959             :   GEN z;
    1960          21 :   if (typ(f) != t_CLOSURE || closure_arity(f) != 2 || closure_is_variadic(f))
    1961           7 :     pari_err_TYPE("setbinop [function needs exactly 2 arguments]",f);
    1962          14 :   lx = lg(x);
    1963          14 :   if (typ(x) != t_VEC) pari_err_TYPE("setbinop", x);
    1964          14 :   if (y == NULL) { /* assume x = y and f symmetric */
    1965           7 :     z = cgetg((((lx-1)*lx) >> 1) + 1, t_VEC);
    1966          28 :     for (i = 1; i < lx; i++)
    1967          63 :       for (j = i; j < lx; j++)
    1968          42 :         gel(z, k++) = closure_callgen2(f, gel(x,i),gel(x,j));
    1969             :   } else {
    1970           7 :     ly = lg(y);
    1971           7 :     if (typ(y) != t_VEC) pari_err_TYPE("setbinop", y);
    1972           7 :     z = cgetg((lx-1)*(ly-1) + 1, t_VEC);
    1973          28 :     for (i = 1; i < lx; i++)
    1974          84 :       for (j = 1; j < ly; j++)
    1975          63 :         gel(z, k++) = closure_callgen2(f, gel(x,i),gel(y,j));
    1976             :   }
    1977          14 :   return gerepileupto(av, gtoset(z));
    1978             : }

Generated by: LCOV version 1.11