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 to exceed 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 - ZX.c (source / functions) Hit Total Coverage
Test: PARI/GP v2.14.0 lcov report (development 27775-aca467eab2) Lines: 700 791 88.5 %
Date: 2022-07-03 07:33:15 Functions: 95 107 88.8 %
Legend: Lines: hit not hit

          Line data    Source code
       1             : /* Copyright (C) 2007  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; either version 2 of the License, or (at your option) any later
       8             : version. It is distributed in the hope that it will be useful, but WITHOUT
       9             : ANY WARRANTY WHATSOEVER.
      10             : 
      11             : Check the License for details. You should have received a copy of it, along
      12             : with the package; see the file 'COPYING'. If not, write to the Free Software
      13             : Foundation, Inc., 51 Franklin Street, Fifth Floor, Boston, MA 02110-1301 USA. */
      14             : 
      15             : #include "pari.h"
      16             : #include "paripriv.h"
      17             : 
      18             : /*******************************************************************/
      19             : /*                                                                 */
      20             : /*                                ZX                               */
      21             : /*                                                                 */
      22             : /*******************************************************************/
      23             : void
      24      297758 : RgX_check_QX(GEN x, const char *s)
      25      297758 : { if (!RgX_is_QX(x)) pari_err_TYPE(stack_strcat(s," [not in Q[X]]"), x); }
      26             : void
      27     1748615 : RgX_check_ZX(GEN x, const char *s)
      28     1748615 : { if (!RgX_is_ZX(x)) pari_err_TYPE(stack_strcat(s," [not in Z[X]]"), x); }
      29             : long
      30       92255 : ZX_max_lg(GEN x)
      31             : {
      32       92255 :   long i, prec = 0, lx = lg(x);
      33             : 
      34      573214 :   for (i=2; i<lx; i++) { long l = lgefint(gel(x,i)); if (l > prec) prec = l; }
      35       92255 :   return prec;
      36             : }
      37             : 
      38             : GEN
      39    38828904 : ZX_add(GEN x, GEN y)
      40             : {
      41             :   long lx,ly,i;
      42             :   GEN z;
      43    38828904 :   lx = lg(x); ly = lg(y); if (lx < ly) swapspec(x,y, lx,ly);
      44    38828904 :   z = cgetg(lx,t_POL); z[1] = x[1];
      45   243016177 :   for (i=2; i<ly; i++) gel(z,i) = addii(gel(x,i),gel(y,i));
      46    77451322 :   for (   ; i<lx; i++) gel(z,i) = icopy(gel(x,i));
      47    38623057 :   if (lx == ly) z = ZX_renormalize(z, lx);
      48    38640852 :   if (!lgpol(z)) { set_avma((pari_sp)(z + lx)); return pol_0(varn(x)); }
      49    37506284 :   return z;
      50             : }
      51             : 
      52             : GEN
      53    27829424 : ZX_sub(GEN x,GEN y)
      54             : {
      55    27829424 :   long i, lx = lg(x), ly = lg(y);
      56             :   GEN z;
      57    27829424 :   if (lx >= ly)
      58             :   {
      59    27407413 :     z = cgetg(lx,t_POL); z[1] = x[1];
      60   113732944 :     for (i=2; i<ly; i++) gel(z,i) = subii(gel(x,i),gel(y,i));
      61    27337239 :     if (lx == ly)
      62             :     {
      63    22057175 :       z = ZX_renormalize(z, lx);
      64    22069944 :       if (!lgpol(z)) { set_avma((pari_sp)(z + lx)); z = pol_0(varn(x)); }
      65             :     }
      66             :     else
      67    16915938 :       for (   ; i<lx; i++) gel(z,i) = icopy(gel(x,i));
      68             :   }
      69             :   else
      70             :   {
      71      422011 :     z = cgetg(ly,t_POL); z[1] = y[1];
      72     1136425 :     for (i=2; i<lx; i++) gel(z,i) = subii(gel(x,i),gel(y,i));
      73     1315245 :     for (   ; i<ly; i++) gel(z,i) = negi(gel(y,i));
      74             :   }
      75    27778512 :   return z;
      76             : }
      77             : 
      78             : GEN
      79      588793 : ZX_neg(GEN x)
      80             : {
      81      588793 :   long i, l = lg(x);
      82      588793 :   GEN y = cgetg(l,t_POL);
      83     2424303 :   y[1] = x[1]; for(i=2; i<l; i++) gel(y,i) = negi(gel(x,i));
      84      588844 :   return y;
      85             : }
      86             : GEN
      87     4918636 : ZX_copy(GEN x)
      88             : {
      89     4918636 :   long i, l = lg(x);
      90     4918636 :   GEN y = cgetg(l, t_POL);
      91     4918637 :   y[1] = x[1];
      92    20324119 :   for (i=2; i<l; i++)
      93             :   {
      94    15405499 :     GEN c = gel(x,i);
      95    15405499 :     gel(y,i) = lgefint(c) == 2? gen_0: icopy(c);
      96             :   }
      97     4918620 :   return y;
      98             : }
      99             : 
     100             : GEN
     101       55612 : scalar_ZX(GEN x, long v)
     102             : {
     103             :   GEN z;
     104       55612 :   if (!signe(x)) return pol_0(v);
     105        4406 :   z = cgetg(3, t_POL);
     106        4406 :   z[1] = evalsigne(1) | evalvarn(v);
     107        4406 :   gel(z,2) = icopy(x); return z;
     108             : }
     109             : 
     110             : GEN
     111       22296 : scalar_ZX_shallow(GEN x, long v)
     112             : {
     113             :   GEN z;
     114       22296 :   if (!signe(x)) return pol_0(v);
     115       19097 :   z = cgetg(3, t_POL);
     116       19097 :   z[1] = evalsigne(1) | evalvarn(v);
     117       19097 :   gel(z,2) = x; return z;
     118             : }
     119             : 
     120             : GEN
     121       87808 : ZX_Z_add(GEN y, GEN x)
     122             : {
     123             :   long lz, i;
     124       87808 :   GEN z = cgetg_copy(y, &lz);
     125       87809 :   if (lz == 2) { set_avma((pari_sp)(z + 2)); return scalar_ZX(x,varn(y)); }
     126       77071 :   z[1] = y[1];
     127       77071 :   gel(z,2) = addii(gel(y,2),x);
     128     1865791 :   for(i=3; i<lz; i++) gel(z,i) = icopy(gel(y,i));
     129       77746 :   if (lz==3) z = ZX_renormalize(z,lz);
     130       77071 :   return z;
     131             : }
     132             : GEN
     133       93165 : ZX_Z_add_shallow(GEN y, GEN x)
     134             : {
     135             :   long lz, i;
     136       93165 :   GEN z = cgetg_copy(y, &lz);
     137       93169 :   if (lz == 2) { set_avma((pari_sp)(z + 2)); return scalar_ZX_shallow(x,varn(y)); }
     138       92987 :   z[1] = y[1];
     139       92987 :   gel(z,2) = addii(gel(y,2),x);
     140      545856 :   for(i=3; i<lz; i++) gel(z,i) = gel(y,i);
     141       92964 :   if (lz==3) z = ZX_renormalize(z,lz);
     142       92970 :   return z;
     143             : }
     144             : 
     145             : GEN
     146       65606 : ZX_Z_sub(GEN y, GEN x)
     147             : {
     148             :   long lz, i;
     149       65606 :   GEN z = cgetg_copy(y, &lz);
     150       65606 :   if (lz == 2)
     151             :   { /* scalarpol(negi(x), v) */
     152         434 :     long v = varn(y);
     153         434 :     set_avma((pari_sp)(z + 2));
     154         434 :     if (!signe(x)) return pol_0(v);
     155         434 :     z = cgetg(3,t_POL);
     156         434 :     z[1] = evalvarn(v) | evalsigne(1);
     157         434 :     gel(z,2) = negi(x); return z;
     158             :   }
     159       65172 :   z[1] = y[1];
     160       65172 :   gel(z,2) = subii(gel(y,2),x);
     161      382439 :   for(i=3; i<lz; i++) gel(z,i) = icopy(gel(y,i));
     162       65172 :   if (lz==3) z = ZX_renormalize(z,lz);
     163       65172 :   return z;
     164             : }
     165             : 
     166             : GEN
     167     2654402 : Z_ZX_sub(GEN x, GEN y)
     168             : {
     169             :   long lz, i;
     170     2654402 :   GEN z = cgetg_copy(y, &lz);
     171     2654547 :   if (lz == 2) { set_avma((pari_sp)(z + 2)); return scalar_ZX(x,varn(y)); }
     172     2654547 :   z[1] = y[1];
     173     2654547 :   gel(z,2) = subii(x, gel(y,2));
     174    13367033 :   for(i=3; i<lz; i++) gel(z,i) = negi(gel(y,i));
     175     2655177 :   if (lz==3) z = ZX_renormalize(z,lz);
     176     2654627 :   return z;
     177             : }
     178             : 
     179             : GEN
     180     5021577 : ZX_Z_divexact(GEN y,GEN x)
     181             : {
     182     5021577 :   long i, l = lg(y);
     183     5021577 :   GEN z = cgetg(l,t_POL); z[1] = y[1];
     184    32362458 :   for(i=2; i<l; i++) gel(z,i) = diviiexact(gel(y,i),x);
     185     5017240 :   return z;
     186             : }
     187             : 
     188             : GEN
     189      977795 : ZX_divuexact(GEN y, ulong x)
     190             : {
     191      977795 :   long i, l = lg(y);
     192      977795 :   GEN z = cgetg(l,t_POL); z[1] = y[1];
     193     9684353 :   for(i=2; i<l; i++) gel(z,i) = diviuexact(gel(y,i),x);
     194      977795 :   return z;
     195             : }
     196             : 
     197             : GEN
     198         112 : zx_z_divexact(GEN y, long x)
     199             : {
     200         112 :   long i, l = lg(y);
     201         112 :   GEN z = cgetg(l,t_VECSMALL); z[1] = y[1];
     202         595 :   for (i=2; i<l; i++) z[i] = y[i]/x;
     203         112 :   return z;
     204             : }
     205             : 
     206             : GEN
     207    16224014 : ZX_Z_mul(GEN y,GEN x)
     208             : {
     209             :   GEN z;
     210             :   long i, l;
     211    16224014 :   if (!signe(x)) return pol_0(varn(y));
     212    13180210 :   l = lg(y); z = cgetg(l,t_POL); z[1] = y[1];
     213    90371155 :   for(i=2; i<l; i++) gel(z,i) = mulii(gel(y,i),x);
     214    13169368 :   return z;
     215             : }
     216             : 
     217             : GEN
     218      251775 : ZX_mulu(GEN y, ulong x)
     219             : {
     220             :   GEN z;
     221             :   long i, l;
     222      251775 :   if (!x) return pol_0(varn(y));
     223      251775 :   l = lg(y); z = cgetg(l,t_POL); z[1] = y[1];
     224     1020939 :   for(i=2; i<l; i++) gel(z,i) = mului(x,gel(y,i));
     225      251775 :   return z;
     226             : }
     227             : 
     228             : GEN
     229     6036943 : ZX_shifti(GEN y, long n)
     230             : {
     231             :   GEN z;
     232             :   long i, l;
     233     6036943 :   l = lg(y); z = cgetg(l,t_POL); z[1] = y[1];
     234    29707245 :   for(i=2; i<l; i++) gel(z,i) = shifti(gel(y,i),n);
     235     6005447 :   return ZX_renormalize(z,l);
     236             : }
     237             : 
     238             : GEN
     239      107334 : ZX_remi2n(GEN y, long n)
     240             : {
     241             :   GEN z;
     242             :   long i, l;
     243      107334 :   l = lg(y); z = cgetg(l,t_POL); z[1] = y[1];
     244     4753179 :   for(i=2; i<l; i++) gel(z,i) = remi2n(gel(y,i),n);
     245      100221 :   return ZX_renormalize(z,l);
     246             : }
     247             : 
     248             : GEN
     249       46488 : ZXT_remi2n(GEN z, long n)
     250             : {
     251       46488 :   if (typ(z) == t_POL)
     252       40061 :     return ZX_remi2n(z, n);
     253             :   else
     254             :   {
     255        6427 :     long i,l = lg(z);
     256        6427 :     GEN x = cgetg(l, t_VEC);
     257       19287 :     for (i=1; i<l; i++) gel(x,i) = ZXT_remi2n(gel(z,i), n);
     258        6426 :     return x;
     259             :   }
     260             : }
     261             : 
     262             : GEN
     263        1281 : zx_to_ZX(GEN z)
     264             : {
     265        1281 :   long i, l = lg(z);
     266        1281 :   GEN x = cgetg(l,t_POL);
     267      100660 :   for (i=2; i<l; i++) gel(x,i) = stoi(z[i]);
     268        1281 :   x[1] = evalsigne(l-2!=0)| z[1]; return x;
     269             : }
     270             : 
     271             : GEN
     272     4885424 : ZX_deriv(GEN x)
     273             : {
     274     4885424 :   long i,lx = lg(x)-1;
     275             :   GEN y;
     276             : 
     277     4885424 :   if (lx<3) return pol_0(varn(x));
     278     4803445 :   y = cgetg(lx,t_POL);
     279    28996590 :   for (i=2; i<lx ; i++) gel(y,i) = mului(i-1,gel(x,i+1));
     280     4801258 :   y[1] = x[1]; return y;
     281             : }
     282             : 
     283             : int
     284     1819066 : ZX_equal(GEN V, GEN W)
     285             : {
     286     1819066 :   long i, l = lg(V);
     287     1819066 :   if (lg(W) != l) return 0;
     288     5988139 :   for (i = 2; i < l; i++)
     289     4652926 :     if (!equalii(gel(V,i), gel(W,i))) return 0;
     290     1335213 :   return 1;
     291             : }
     292             : 
     293             : static long
     294   231403555 : ZX_valspec(GEN x, long nx)
     295             : {
     296             :   long vx;
     297   302695744 :   for (vx = 0; vx<nx ; vx++)
     298   302705378 :     if (signe(gel(x,vx))) break;
     299   231403555 :   return vx;
     300             : }
     301             : 
     302             : long
     303      327094 : ZX_val(GEN x)
     304             : {
     305             :   long vx;
     306      327094 :   if (!signe(x)) return LONG_MAX;
     307      212606 :   for (vx = 0;; vx++)
     308      262778 :     if (signe(gel(x,2+vx))) break;
     309      212606 :   return vx;
     310             : }
     311             : long
     312    22897808 : ZX_valrem(GEN x, GEN *Z)
     313             : {
     314             :   long vx;
     315    22897808 :   if (!signe(x)) { *Z = pol_0(varn(x)); return LONG_MAX; }
     316    22897808 :   for (vx = 0;; vx++)
     317    43507165 :     if (signe(gel(x,2+vx))) break;
     318    22897808 :   *Z = RgX_shift_shallow(x, -vx);
     319    22897840 :   return vx;
     320             : }
     321             : 
     322             : GEN
     323          21 : ZX_div_by_X_1(GEN a, GEN *r)
     324             : {
     325          21 :   long l = lg(a), i;
     326          21 :   GEN a0, z0, z = cgetg(l-1, t_POL);
     327          21 :   z[1] = a[1];
     328          21 :   a0 = a + l-1;
     329          21 :   z0 = z + l-2; *z0 = *a0--;
     330        4732 :   for (i=l-3; i>1; i--) /* z[i] = a[i+1] + z[i+1] */
     331             :   {
     332        4711 :     GEN t = addii(gel(a0--,0), gel(z0--,0));
     333        4711 :     gel(z0,0) = t;
     334             :   }
     335          21 :   if (r) *r = addii(gel(a0,0), gel(z0,0));
     336          21 :   return z;
     337             : }
     338             : 
     339             : /* return P(X + c) using destructive Horner, optimize for c = 1,-1 */
     340             : static GEN
     341     1242829 : ZX_translate_basecase(GEN P, GEN c)
     342             : {
     343     1242829 :   pari_sp av = avma;
     344             :   GEN Q, R;
     345             :   long i, k, n;
     346             : 
     347     1242829 :   if (!signe(P) || !signe(c)) return ZX_copy(P);
     348     1203353 :   Q = leafcopy(P);
     349     1203361 :   R = Q+2; n = degpol(P);
     350     1203361 :   if (equali1(c))
     351             :   {
     352     5887117 :     for (i=1; i<=n; i++)
     353             :     {
     354    35069890 :       for (k=n-i; k<n; k++) gel(R,k) = addii(gel(R,k), gel(R,k+1));
     355     4935619 :       if (gc_needed(av,2))
     356             :       {
     357           0 :         if(DEBUGMEM>1) pari_warn(warnmem,"ZX_translate(1), i = %ld/%ld", i,n);
     358           0 :         Q = gerepilecopy(av, Q); R = Q+2;
     359             :       }
     360             :     }
     361             :   }
     362      251864 :   else if (equalim1(c))
     363             :   {
     364       55098 :     for (i=1; i<=n; i++)
     365             :     {
     366      180940 :       for (k=n-i; k<n; k++) gel(R,k) = subii(gel(R,k), gel(R,k+1));
     367       45566 :       if (gc_needed(av,2))
     368             :       {
     369           0 :         if(DEBUGMEM>1) pari_warn(warnmem,"ZX_translate(-1), i = %ld/%ld", i,n);
     370           0 :         Q = gerepilecopy(av, Q); R = Q+2;
     371             :       }
     372             :     }
     373             :   }
     374             :   else
     375             :   {
     376     1494099 :     for (i=1; i<=n; i++)
     377             :     {
     378     6748548 :       for (k=n-i; k<n; k++) gel(R,k) = addmulii_inplace(gel(R,k), c, gel(R,k+1));
     379     1251767 :       if (gc_needed(av,2))
     380             :       {
     381           0 :         if(DEBUGMEM>1) pari_warn(warnmem,"ZX_translate, i = %ld/%ld", i,n);
     382           0 :         Q = gerepilecopy(av, Q); R = Q+2;
     383             :       }
     384             :     }
     385             :   }
     386     1203085 :   return gerepilecopy(av, Q);
     387             : }
     388             : 
     389             : static GEN
     390           0 : Z_Xpm1_powu(long n, long s, long v)
     391             : {
     392             :   long d, k;
     393             :   GEN C;
     394           0 :   if (!n) return pol_1(v);
     395           0 :   d = (n + 1) >> 1;
     396           0 :   C = cgetg(n+3, t_POL);
     397           0 :   C[1] = evalsigne(1)| evalvarn(v);
     398           0 :   gel(C,2) = gen_1;
     399           0 :   gel(C,3) = utoipos(n);
     400           0 :   for (k=2; k <= d; k++)
     401           0 :     gel(C,k+2) = diviuexact(mului(n-k+1, gel(C,k+1)), k);
     402           0 :   if (s < 0)
     403           0 :     for (k = odd(n)? 0: 1; k <= d; k += 2)
     404           0 :       togglesign_safe(&gel(C,k+2));
     405           0 :   if (s > 0 || !odd(n))
     406           0 :     for (k = d+1; k <= n; k++) gel(C,k+2) = gel(C,n-k+2);
     407             :   else
     408           0 :     for (k = d+1; k <= n; k++) gel(C,k+2) = negi(gel(C,n-k+2));
     409           0 :   return C;
     410             : }
     411             : /* return (x+u)^n */
     412             : static GEN
     413           0 : Z_XpN_powu(GEN u, long n, long v)
     414             : {
     415             :   pari_sp av;
     416             :   long k;
     417             :   GEN B, C, V;
     418           0 :   if (!n) return pol_1(v);
     419           0 :   if (is_pm1(u))
     420           0 :     return Z_Xpm1_powu(n, signe(u), v);
     421           0 :   av = avma;
     422           0 :   V = gpowers(u, n);
     423           0 :   B = vecbinomial(n);
     424           0 :   C = cgetg(n+3, t_POL);
     425           0 :   C[1] = evalsigne(1)| evalvarn(v);
     426           0 :   for (k=1; k <= n+1; k++)
     427           0 :     gel(C,k+1) = mulii(gel(V,n+2-k), gel(B,k));
     428           0 :   return gerepileupto(av, C);
     429             : }
     430             : 
     431             : GEN
     432     1242820 : ZX_translate(GEN P, GEN c)
     433             : {
     434     1242820 :   pari_sp av = avma;
     435     1242820 :   long n = degpol(P);
     436     1242826 :   if (n < 220)
     437     1242826 :     return ZX_translate_basecase(P, c);
     438             :   else
     439             :   {
     440           0 :     long d = n >> 1;
     441           0 :     GEN Q = ZX_translate(RgX_shift_shallow(P, -d), c);
     442           0 :     GEN R = ZX_translate(RgXn_red_shallow(P, d), c);
     443           0 :     GEN S = Z_XpN_powu(c, d, varn(P));
     444           0 :     return gerepileupto(av, ZX_add(ZX_mul(Q, S), R));
     445             :   }
     446             : }
     447             : 
     448             : /* P(ax + b) */
     449             : GEN
     450      975101 : ZX_affine(GEN P, GEN a, GEN b)
     451             : {
     452      975101 :   if (signe(b)) P = ZX_translate(P, b);
     453      975101 :   return ZX_unscale(P, a);
     454             : }
     455             : 
     456             : GEN
     457      478043 : ZX_Z_eval(GEN x, GEN y)
     458             : {
     459      478043 :   long i = lg(x)-1, j;
     460      478043 :   pari_sp av = avma;
     461             :   GEN t, r;
     462             : 
     463      478043 :   if (i<=2) return (i==2)? icopy(gel(x,2)): gen_0;
     464      442770 :   if (!signe(y)) return icopy(gel(x,2));
     465             : 
     466      442658 :   t = gel(x,i); i--;
     467             : #if 0 /* standard Horner's rule */
     468             :   for ( ; i>=2; i--)
     469             :     t = addii(mulii(t,y),gel(x,i));
     470             : #endif
     471             :   /* specific attention to sparse polynomials */
     472     2309968 :   for ( ; i>=2; i = j-1)
     473             :   {
     474     2017118 :     for (j = i; !signe(gel(x,j)); j--)
     475      149798 :       if (j==2)
     476             :       {
     477        2933 :         if (i != j) y = powiu(y, i-j+1);
     478        2933 :         return gerepileuptoint(av, mulii(t,y));
     479             :       }
     480     1867320 :     r = (i==j)? y: powiu(y, i-j+1);
     481     1867324 :     t = addii(mulii(t,r), gel(x,j));
     482     1867307 :     if (gc_needed(av,2))
     483             :     {
     484           0 :       if (DEBUGMEM>1) pari_warn(warnmem,"ZX_Z_eval: i = %ld",i);
     485           0 :       t = gerepileuptoint(av, t);
     486             :     }
     487             :   }
     488      439715 :   return gerepileuptoint(av, t);
     489             : }
     490             : 
     491             : /* Return 2^(n degpol(P))  P(x >> n) */
     492             : GEN
     493           0 : ZX_rescale2n(GEN P, long n)
     494             : {
     495           0 :   long i, l = lg(P), ni = n;
     496             :   GEN Q;
     497           0 :   if (l==2) return pol_0(varn(P));
     498           0 :   Q = cgetg(l,t_POL);
     499           0 :   gel(Q,l-1) = icopy(gel(P,l-1));
     500           0 :   for (i=l-2; i>=2; i--)
     501             :   {
     502           0 :     gel(Q,i) = shifti(gel(P,i), ni);
     503           0 :     ni += n;
     504             :   }
     505           0 :   Q[1] = P[1]; return Q;
     506             : }
     507             : 
     508             : /* Return h^deg(P) P(x / h), not memory clean. h integer, P ZX */
     509             : GEN
     510       36584 : ZX_rescale(GEN P, GEN h)
     511             : {
     512       36584 :   long l = lg(P);
     513       36584 :   GEN Q = cgetg(l,t_POL);
     514       36587 :   if (l != 2)
     515             :   {
     516       36586 :     long i = l-1;
     517       36586 :     GEN hi = h;
     518       36586 :     gel(Q,i) = gel(P,i);
     519       36586 :     if (l != 3) { i--; gel(Q,i) = mulii(gel(P,i), h); }
     520      113066 :     for (i--; i>=2; i--) { hi = mulii(hi,h); gel(Q,i) = mulii(gel(P,i), hi); }
     521             :   }
     522       36583 :   Q[1] = P[1]; return Q;
     523             : }
     524             : /* Return h^(deg(P)-1) P(x / h), P!=0, h=lt(P), memory unclean; monic result */
     525             : GEN
     526           0 : ZX_rescale_lt(GEN P)
     527             : {
     528           0 :   long l = lg(P);
     529           0 :   GEN Q = cgetg(l,t_POL);
     530           0 :   gel(Q,l-1) = gen_1;
     531           0 :   if (l != 3)
     532             :   {
     533           0 :     long i = l-1;
     534           0 :     GEN h = gel(P,i), hi = h;
     535           0 :     i--; gel(Q,i) = gel(P,i);
     536           0 :     if (l != 4) { i--; gel(Q,i) = mulii(gel(P,i), h); }
     537           0 :     for (i--; i>=2; i--) { hi = mulii(hi,h); gel(Q,i) = mulii(gel(P,i), hi); }
     538             :   }
     539           0 :   Q[1] = P[1]; return Q;
     540             : }
     541             : 
     542             : /*Eval x in 2^(k*BIL) in linear time*/
     543             : static GEN
     544    61286808 : ZX_eval2BILspec(GEN x, long k, long nx)
     545             : {
     546    61286808 :   pari_sp av = avma;
     547    61286808 :   long i,j, lz = k*nx, ki;
     548    61286808 :   GEN pz = cgetipos(2+lz);
     549    61275449 :   GEN nz = cgetipos(2+lz);
     550  2383166084 :   for(i=0; i < lz; i++)
     551             :   {
     552  2321870086 :     *int_W(pz,i) = 0UL;
     553  2321870086 :     *int_W(nz,i) = 0UL;
     554             :   }
     555   718910092 :   for(i=0, ki=0; i<nx; i++, ki+=k)
     556             :   {
     557   657614094 :     GEN c = gel(x,i);
     558   657614094 :     long lc = lgefint(c)-2;
     559   657614094 :     if (signe(c)==0) continue;
     560   479635282 :     if (signe(c) > 0)
     561  1280353818 :       for (j=0; j<lc; j++) *int_W(pz,j+ki) = *int_W(c,j);
     562             :     else
     563   289203610 :       for (j=0; j<lc; j++) *int_W(nz,j+ki) = *int_W(c,j);
     564             :   }
     565    61295998 :   pz = int_normalize(pz,0);
     566    61300605 :   nz = int_normalize(nz,0); return gerepileuptoint(av, subii(pz,nz));
     567             : }
     568             : 
     569             : static long
     570   113046660 : ZX_expispec(GEN x, long nx)
     571             : {
     572   113046660 :   long i, m = 0;
     573   916259932 :   for(i = 0; i < nx; i++)
     574             :   {
     575   803253005 :     long e = expi(gel(x,i));
     576   803213272 :     if (e > m) m = e;
     577             :   }
     578   113006927 :   return m;
     579             : }
     580             : 
     581             : static GEN
     582    38214731 : Z_mod2BIL_ZX(GEN x, long bs, long d, long vx)
     583             : {
     584    38214731 :   long i, offset, lm = lgefint(x)-2, l = d+vx+3, sx = signe(x);
     585    38214731 :   GEN s1 = int2n(bs*BITS_IN_LONG), pol = cgetg(l, t_POL);
     586    38198728 :   int carry = 0;
     587    38198728 :   pol[1] = evalsigne(1);
     588    58540139 :   for (i=0; i<vx; i++) gel(pol,i+2) = gen_0;
     589   723386914 :   for (offset=0; i <= d+vx; i++, offset += bs)
     590             :   {
     591   685660252 :     pari_sp av = avma;
     592   685660252 :     long lz = minss(bs, lm-offset);
     593   687372921 :     GEN z = lz > 0 ?adduispec_offset(carry, x, offset, lz): utoi(carry);
     594   684849428 :     if (lgefint(z) == 3+bs) { carry = 1; z = gen_0;}
     595             :     else
     596             :     {
     597   672597542 :       carry = (lgefint(z) == 2+bs && (HIGHBIT & *int_W(z,bs-1)));
     598   672597542 :       if (carry)
     599    49870981 :         z = gerepileuptoint(av, (sx==-1)? subii(s1,z): subii(z,s1));
     600   622726561 :       else if (sx==-1) togglesign(z);
     601             :     }
     602   685188186 :     gel(pol,i+2) = z;
     603             :   }
     604    37726662 :   return ZX_renormalize(pol,l);
     605             : }
     606             : 
     607             : static GEN
     608    13916442 : ZX_sqrspec_sqri(GEN x, long nx, long ex, long v)
     609             : {
     610    13916442 :   long e = 2*ex + expu(nx) + 3;
     611    13916305 :   long N = divsBIL(e)+1;
     612    13916260 :   GEN  z = sqri(ZX_eval2BILspec(x,N,nx));
     613    13914618 :   return Z_mod2BIL_ZX(z, N, nx*2-2, v);
     614             : }
     615             : 
     616             : static GEN
     617    23237965 : ZX_mulspec_mulii(GEN x, GEN y, long nx, long ny, long ex, long ey, long v)
     618             : {
     619    23237965 :   long e = ex + ey + expu(minss(nx,ny)) + 3;
     620    23237515 :   long N = divsBIL(e)+1;
     621    23237547 :   GEN  z = mulii(ZX_eval2BILspec(x,N,nx), ZX_eval2BILspec(y,N,ny));
     622    23236744 :   return Z_mod2BIL_ZX(z, N, nx+ny-2, v);
     623             : }
     624             : 
     625             : INLINE GEN
     626   239509770 : ZX_sqrspec_basecase_limb(GEN x, long a, long i)
     627             : {
     628   239509770 :   pari_sp av = avma;
     629   239509770 :   GEN s = gen_0;
     630   239509770 :   long j, l = (i+1)>>1;
     631   525404820 :   for (j=a; j<l; j++)
     632             :   {
     633   286746479 :     GEN xj = gel(x,j), xx = gel(x,i-j);
     634   286746479 :     if (signe(xj) && signe(xx))
     635   282769528 :       s = addii(s, mulii(xj, xx));
     636             :   }
     637   238658341 :   s = shifti(s,1);
     638   239019907 :   if ((i&1) == 0)
     639             :   {
     640   145719884 :     GEN t = gel(x, i>>1);
     641   145719884 :     if (signe(t))
     642   145466777 :       s = addii(s, sqri(t));
     643             :   }
     644   238531507 :   return gerepileuptoint(av,s);
     645             : }
     646             : 
     647             : static GEN
     648    51741084 : ZX_sqrspec_basecase(GEN x, long nx, long v)
     649             : {
     650             :   long i, lz, nz;
     651             :   GEN z;
     652             : 
     653    51741084 :   lz = (nx << 1) + 1; nz = lz-2;
     654    51741084 :   lz += v;
     655    51741084 :   z = cgetg(lz,t_POL); z[1] = evalsigne(1); z += 2;
     656    51924694 :   for (i=0; i<v; i++) gel(z++, 0) = gen_0;
     657   197412109 :   for (i=0; i<nx; i++)
     658   145771775 :     gel(z,i) = ZX_sqrspec_basecase_limb(x, 0, i);
     659   145745955 :   for (  ; i<nz; i++) gel(z,i) = ZX_sqrspec_basecase_limb(x, i-nx+1, i);
     660    51633698 :   z -= v+2; return z;
     661             : }
     662             : 
     663             : static GEN
     664    27748026 : Z_sqrshiftspec_ZX(GEN x, long vx)
     665             : {
     666    27748026 :   long i, nz = 2*vx+1;
     667    27748026 :   GEN z = cgetg(2+nz, t_POL);
     668    27741696 :   z[1] = evalsigne(1);
     669    29985642 :   for(i=2;i<nz+1;i++) gel(z,i) = gen_0;
     670    27741696 :   gel(z,nz+1) = sqri(x);
     671    27680742 :   return z;
     672             : }
     673             : 
     674             : static GEN
     675    34358009 : Z_ZX_mulshiftspec(GEN x, GEN y, long ny, long vz)
     676             : {
     677    34358009 :   long i, nz = vz+ny;
     678    34358009 :   GEN z = cgetg(2+nz, t_POL);
     679    34357793 :   z[1] = evalsigne(1);
     680    82166425 :   for (i=0; i<vz; i++)   gel(z,i+2)    = gen_0;
     681   103791309 :   for (i=0; i<ny; i++) gel(z,i+vz+2) = mulii(x, gel(y,i));
     682    34345531 :   return z;
     683             : }
     684             : 
     685             : GEN
     686    93781409 : ZX_sqrspec(GEN x, long nx)
     687             : {
     688             : #ifdef PARI_KERNEL_GMP
     689    53968007 :   const long low[]={ 17, 32, 96, 112, 160, 128, 128, 160, 160, 160, 160, 160, 176, 192, 192, 192, 192, 192, 224, 224, 224, 240, 240, 240, 272, 288, 288, 240, 288, 304, 304, 304, 304, 304, 304, 352, 352, 368, 352, 352, 352, 368, 368, 432, 432, 496, 432, 496, 496};
     690    53968007 :   const long high[]={ 102860, 70254, 52783, 27086, 24623, 18500, 15289, 13899, 12635, 11487, 10442, 9493, 8630, 7845, 7132, 7132, 6484, 6484, 5894, 5894, 4428, 4428, 3660, 4428, 3660, 3660, 2749, 2499, 2272, 2066, 1282, 1282, 1166, 1166, 1166, 1166, 1166, 1166, 1166, 963, 963, 724, 658, 658, 658, 528, 528, 528, 528};
     691             : #else
     692    39813402 :   const long low[]={ 17, 17, 32, 32, 96, 96, 96, 96, 96, 96, 96, 96, 96, 96, 112, 112, 128, 112, 112, 112, 112, 112, 128, 128, 160, 160, 112, 128, 128, 160, 160, 160, 160, 160, 160, 160, 160, 160, 160, 160, 160, 160, 160, 160, 160, 160, 176, 160, 160, 176, 160, 160, 176, 176, 208, 176, 176, 176, 192, 192, 176, 176, 224, 176, 224, 224, 176, 224, 224, 224, 176, 176, 176, 176, 176, 176, 176, 176, 224, 176, 176, 224, 224, 224, 224, 224, 224, 224, 240, 288, 240, 288, 288, 240, 288, 288, 240, 240, 304, 304};
     693    39813402 :   const long high[]={ 165657, 85008, 52783, 43622, 32774, 27086, 22385, 15289, 13899, 12635, 11487, 10442, 9493, 7845, 6484, 6484, 5894, 5894, 4871, 4871, 4428, 4026, 3660, 3660, 3660, 3327, 3327, 3024, 2749, 2749, 2272, 2749, 2499, 2499, 2272, 1878, 1878, 1878, 1707, 1552, 1552, 1552, 1552, 1552, 1411, 1411, 1411, 1282, 1282, 1282, 1282, 1282, 1166, 1166, 1166, 1166, 1166, 1166, 1166, 1166, 1060, 1060, 963, 963, 963, 963, 963, 963, 963, 963, 963, 963, 963, 876, 876, 876, 876, 796, 658, 724, 658, 724, 658, 658, 658, 658, 658, 658, 658, 658, 658, 658, 658, 658, 336, 658, 658, 592, 336, 336};
     694             : #endif
     695    93781409 :   const long nblow = numberof(low);
     696             :   pari_sp av;
     697             :   long ex, vx;
     698             :   GEN z;
     699    93781409 :   if (!nx) return pol_0(0);
     700    93301084 :   vx = ZX_valspec(x,nx); nx-=vx; x+=vx;
     701    93420476 :   if (nx==1) return Z_sqrshiftspec_ZX(gel(x, 0), vx);
     702    65671285 :   av = avma;
     703    65671285 :   ex = ZX_expispec(x,nx);
     704    65645214 :   if (nx-2 < nblow && low[nx-2]<=ex && ex<=high[nx-2])
     705    51743997 :     z = ZX_sqrspec_basecase(x, nx, 2*vx);
     706             :   else
     707    13901217 :     z = ZX_sqrspec_sqri(x, nx, ex, 2*vx);
     708    65604867 :   return gerepileupto(av, z);
     709             : }
     710             : 
     711             : GEN
     712    93784819 : ZX_sqr(GEN x)
     713             : {
     714    93784819 :   GEN z = ZX_sqrspec(x+2, lgpol(x));
     715    93612528 :   z[1] = x[1];
     716    93612528 :   return z;
     717             : }
     718             : 
     719             : GEN
     720    78617174 : ZX_mulspec(GEN x, GEN y, long nx, long ny)
     721             : {
     722             :   pari_sp av;
     723             :   long ex, ey, vx, vy, v;
     724    78617174 :   if (!nx || !ny) return pol_0(0);
     725    68999551 :   vx = ZX_valspec(x,nx); nx-=vx; x += vx;
     726    68999535 :   vy = ZX_valspec(y,ny); ny-=vy; y += vy;
     727    69003453 :   v = vx + vy;
     728    69003453 :   if (nx==1) return Z_ZX_mulshiftspec(gel(x,0), y, ny, v);
     729    44378724 :   if (ny==1) return Z_ZX_mulshiftspec(gel(y,0), x, nx, v);
     730    34645241 :   if (nx == 2 && ny == 2)
     731             :   {
     732    11407862 :     GEN a0 = gel(x,0), a1 = gel(x,1), A0, A1, A2;
     733    11407862 :     GEN b0 = gel(y,0), b1 = gel(y,1), z = cgetg(5 + v, t_POL);
     734             :     long i;
     735    11407864 :     z[1] = evalvarn(0) | evalsigne(1);
     736    11407864 :     A0 = mulii(a0, b0);
     737    11407503 :     A2 = mulii(a1, b1); av = avma;
     738    11407406 :     A1 = gerepileuptoint(av, subii(addii(A0,A2),
     739             :                                    mulii(subii(a1,a0), subii(b1,b0))));
     740    11407815 :     i = 4 + v;
     741    11407815 :     gel(z,i--) = A2;
     742    11407815 :     gel(z,i--) = A1;
     743    13484434 :     gel(z,i--) = A0; while (i > 1) gel(z, i--) = gen_0;
     744    11407815 :     return z;
     745             :   }
     746             : #if 0
     747             :   /* generically slower even when degrees differ a lot; sometimes about twice
     748             :    * faster when bitsize is moderate */
     749             :   if (DEBUGVAR)
     750             :     return RgX_mulspec(x - vx, y - vy, nx + vx, ny + vy);
     751             : #endif
     752    23237379 :   av = avma;
     753    23237379 :   ex = ZX_expispec(x, nx);
     754    23236820 :   ey = ZX_expispec(y, ny);
     755    23237476 :   return gerepileupto(av,  ZX_mulspec_mulii(x,y,nx,ny,ex,ey,v));
     756             : }
     757             : GEN
     758    77341106 : ZX_mul(GEN x, GEN y)
     759             : {
     760             :   GEN z;
     761    77341106 :   if (x == y) return ZX_sqr(x);
     762    76605143 :   z = ZX_mulspec(x+2,y+2,lgpol(x),lgpol(y));
     763    76595079 :   z[1] = x[1];
     764    76595079 :   if (!signe(y)) z[1] &= VARNBITS;
     765    76595079 :   return z;
     766             : }
     767             : 
     768             : /* x,y two ZX in the same variable; assume y is monic */
     769             : GEN
     770     9409036 : ZX_rem(GEN x, GEN y)
     771             : {
     772             :   long vx, dx, dy, dz, i, j, sx, lr;
     773             :   pari_sp av0, av;
     774             :   GEN z,p1,rem;
     775             : 
     776     9409036 :   vx = varn(x);
     777     9409036 :   dy = degpol(y);
     778     9408994 :   dx = degpol(x);
     779     9409025 :   if (dx < dy) return ZX_copy(x);
     780     4820971 :   if (!dy) return pol_0(vx); /* y is constant */
     781     4818346 :   av0 = avma; dz = dx-dy;
     782     4818346 :   z=cgetg(dz+3,t_POL); z[1] = x[1];
     783     4818389 :   x += 2; y += 2; z += 2;
     784             : 
     785     4818389 :   p1 = gel(x,dx);
     786     4818389 :   gel(z,dz) = icopy(p1);
     787    26033355 :   for (i=dx-1; i>=dy; i--)
     788             :   {
     789    21214990 :     av=avma; p1=gel(x,i);
     790   207786766 :     for (j=i-dy+1; j<=i && j<=dz; j++)
     791   186581492 :       p1 = subii(p1, mulii(gel(z,j),gel(y,i-j)));
     792    21205274 :     gel(z,i-dy) = avma == av? icopy(p1): gerepileuptoint(av, p1);
     793             :   }
     794     4818365 :   rem = (GEN)avma; av = (pari_sp)new_chunk(dx+3);
     795     4818337 :   for (sx=0; ; i--)
     796             :   {
     797     2447071 :     p1 = gel(x,i);
     798    51851524 :     for (j=0; j<=i && j<=dz; j++)
     799    44586282 :       p1 = subii(p1, mulii(gel(z,j),gel(y,i-j)));
     800     7265242 :     if (signe(p1)) { sx = 1; break; }
     801     2700142 :     if (!i) break;
     802     2447151 :     set_avma(av);
     803             :   }
     804     4818091 :   lr=i+3; rem -= lr;
     805     4818091 :   rem[0] = evaltyp(t_POL) | evallg(lr);
     806     4818039 :   rem[1] = z[-1];
     807     4818039 :   p1 = gerepileuptoint((pari_sp)rem, p1);
     808     4818841 :   rem += 2; gel(rem,i) = p1;
     809    30185923 :   for (i--; i>=0; i--)
     810             :   {
     811    25367646 :     av=avma; p1 = gel(x,i);
     812   242254867 :     for (j=0; j<=i && j<=dz; j++)
     813   216911250 :       p1 = subii(p1, mulii(gel(z,j),gel(y,i-j)));
     814    25343617 :     gel(rem,i) = avma == av? icopy(p1): gerepileuptoint(av, p1);
     815             :   }
     816     4818277 :   rem -= 2;
     817     4818277 :   if (!sx) (void)ZX_renormalize(rem, lr);
     818     4818278 :   return gerepileupto(av0,rem);
     819             : }
     820             : 
     821             : /* return x(1) */
     822             : GEN
     823        4718 : ZX_eval1(GEN x)
     824             : {
     825        4718 :   pari_sp av = avma;
     826        4718 :   long i = lg(x)-1;
     827             :   GEN s;
     828        4718 :   if (i < 2) return gen_0;
     829        4193 :   s = gel(x,i); i--;
     830        4193 :   if (i == 1) return icopy(s);
     831       34608 :   for ( ; i>=2; i--)
     832             :   {
     833       31423 :     GEN c = gel(x,i);
     834       31423 :     if (signe(c)) s = addii(s, c);
     835             :   }
     836        3185 :   return gerepileuptoint(av,s);
     837             : }
     838             : 
     839             : /* reduce T mod X^n - 1. Shallow function */
     840             : GEN
     841     1800507 : ZX_mod_Xnm1(GEN T, ulong n)
     842             : {
     843     1800507 :   long i, j, L = lg(T), l = n+2;
     844             :   GEN S;
     845     1800507 :   if (L <= l) return T;
     846     1609470 :   S = cgetg(l, t_POL);
     847     1608599 :   S[1] = T[1];
     848    12370611 :   for (i = 2; i < l; i++) gel(S,i) = gel(T,i);
     849     7513682 :   for (j = 2; i < L; i++) {
     850     5911214 :     gel(S,j) = addii(gel(S,j), gel(T,i));
     851     5905083 :     if (++j == l) j = 2;
     852             :   }
     853     1602468 :   return normalizepol_lg(S, l);
     854             : }
     855             : 
     856             : static GEN
     857           0 : _ZX_mul(void* E, GEN x, GEN y)
     858           0 : { (void) E; return ZX_mul(x, y); }
     859             : static GEN
     860           0 : _ZX_sqr(void *E, GEN x)
     861           0 : { (void) E; return ZX_sqr(x); }
     862             : 
     863             : static GEN
     864           0 : _ZX_divrem(void * E, GEN x, GEN y, GEN *r)
     865           0 : { (void) E; return RgX_divrem(x, y, r); }
     866             : 
     867             : static GEN
     868           0 : _ZX_add(void * E, GEN x, GEN y)
     869           0 : { (void) E; return ZX_add(x, y); }
     870             : 
     871             : static struct bb_ring ZX_ring = { _ZX_add,_ZX_mul,_ZX_sqr };
     872             : 
     873             : GEN
     874           0 : ZX_digits(GEN x, GEN T)
     875             : {
     876           0 :   long d = degpol(T), n = (lgpol(x)+d-1)/d;
     877           0 :   return gen_digits(x, T, n, NULL, &ZX_ring, _ZX_divrem);
     878             : }
     879             : 
     880             : GEN
     881           0 : ZXV_ZX_fromdigits(GEN x, GEN T)
     882           0 : { return gen_fromdigits(x,T, NULL, &ZX_ring); }
     883             : 
     884             : /*******************************************************************/
     885             : /*                                                                 */
     886             : /*                                ZXV                              */
     887             : /*                                                                 */
     888             : /*******************************************************************/
     889             : 
     890             : int
     891          35 : ZXV_equal(GEN V, GEN W)
     892             : {
     893          35 :   long l = lg(V);
     894          35 :   if (l!=lg(W)) return 0;
     895          35 :   while (--l > 0)
     896          35 :     if (!ZX_equal(gel(V,l), gel(W,l))) return 0;
     897           0 :   return 1;
     898             : }
     899             : 
     900             : GEN
     901        7476 : ZXV_Z_mul(GEN x, GEN y)
     902       29904 : { pari_APPLY_same(ZX_Z_mul(gel(x,i), y)) }
     903             : 
     904             : GEN
     905           0 : ZXV_remi2n(GEN x, long N)
     906           0 : { pari_APPLY_same(ZX_remi2n(gel(x,i), N)) }
     907             : 
     908             : GEN
     909       82859 : ZXV_dotproduct(GEN x, GEN y)
     910             : {
     911       82859 :   pari_sp av = avma;
     912       82859 :   long i, lx = lg(x);
     913             :   GEN c;
     914       82859 :   if (lx == 1) return pol_0(varn(x));
     915       82859 :   c = ZX_mul(gel(x,1), gel(y,1));
     916      391475 :   for (i = 2; i < lx; i++)
     917             :   {
     918      308616 :     GEN t = ZX_mul(gel(x,i), gel(y,i));
     919      308616 :     if (signe(t)) c = ZX_add(c, t);
     920             :   }
     921       82859 :   return gerepileupto(av, c);
     922             : }
     923             : 
     924             : GEN
     925      122417 : ZXC_to_FlxC(GEN x, ulong p, long sv)
     926     5048669 : { pari_APPLY_type(t_COL,typ(gel(x,i))==t_INT ?
     927             :                   Z_to_Flx(gel(x,i), p, sv): ZX_to_Flx(gel(x,i), p)) }
     928             : 
     929             : GEN
     930       11498 : ZXM_to_FlxM(GEN x, ulong p, long sv)
     931      101979 : { pari_APPLY_same(ZXC_to_FlxC(gel(x,i), p, sv)) }
     932             : 
     933             : /*******************************************************************/
     934             : /*                                                                 */
     935             : /*                                ZXQM                             */
     936             : /*                                                                 */
     937             : /*******************************************************************/
     938             : 
     939             : GEN
     940     2470100 : ZXn_mul(GEN x, GEN y, long n)
     941     2470100 : { return RgXn_red_shallow(ZX_mul(x, y), n); }
     942             : 
     943             : GEN
     944        1407 : ZXn_sqr(GEN x, long n)
     945        1407 : { return RgXn_red_shallow(ZX_sqr(x), n); }
     946             : 
     947             : /*******************************************************************/
     948             : /*                                                                 */
     949             : /*                                ZXQM                             */
     950             : /*                                                                 */
     951             : /*******************************************************************/
     952             : 
     953             : static long
     954     3073024 : ZX_expi(GEN x)
     955             : {
     956     3073024 :   if (signe(x)==0) return 0;
     957     1389503 :   if (typ(x)==t_INT) return expi(x);
     958      911725 :   return ZX_expispec(x+2, lgpol(x));
     959             : }
     960             : 
     961             : static long
     962      713697 : ZXC_expi(GEN x)
     963             : {
     964      713697 :   long i, l = lg(x), m=0;
     965     3786721 :   for(i = 1; i < l; i++)
     966             :   {
     967     3073024 :     long e = ZX_expi(gel(x,i));
     968     3073024 :     if (e > m) m = e;
     969             :   }
     970      713697 :   return m;
     971             : }
     972             : 
     973             : static long
     974      137354 : ZXM_expi(GEN x)
     975             : {
     976      137354 :   long i, l = lg(x), m=0;
     977      851051 :   for(i = 1; i < l; i++)
     978             :   {
     979      713697 :     long e = ZXC_expi(gel(x,i));
     980      713697 :     if (e > m) m = e;
     981             :   }
     982      137354 :   return m;
     983             : }
     984             : 
     985             : static GEN
     986     3073024 : ZX_eval2BIL(GEN x, long k)
     987             : {
     988     3073024 :   if (signe(x)==0) return gen_0;
     989     1389503 :   if (typ(x)==t_INT) return x;
     990      911725 :   return ZX_eval2BILspec(x+2, k, lgpol(x));
     991             : }
     992             : 
     993             : /*Eval x in 2^(k*BIL) in linear time*/
     994             : static GEN
     995      713697 : ZXC_eval2BIL(GEN x, long k)
     996             : {
     997      713697 :   long i, lx = lg(x);
     998      713697 :   GEN A = cgetg(lx, t_COL);
     999     3786721 :   for (i=1; i<lx; i++) gel(A,i) = ZX_eval2BIL(gel(x,i), k);
    1000      713697 :   return A;
    1001             : }
    1002             : 
    1003             : static GEN
    1004      137354 : ZXM_eval2BIL(GEN x, long k)
    1005             : {
    1006      137354 :   long i, lx = lg(x);
    1007      137354 :   GEN A = cgetg(lx, t_MAT);
    1008      851051 :   for (i=1; i<lx; i++) gel(A,i) = ZXC_eval2BIL(gel(x,i), k);
    1009      137354 :   return A;
    1010             : }
    1011             : 
    1012             : static GEN
    1013       61124 : Z_mod2BIL_ZXQ(GEN x, long bs, GEN T)
    1014             : {
    1015       61124 :   pari_sp av = avma;
    1016       61124 :   long v = varn(T), d = 2*(degpol(T)-1);
    1017       61124 :   GEN z = Z_mod2BIL_ZX(x, bs, d, 0);
    1018       61124 :   setvarn(z, v);
    1019       61124 :   return gerepileupto(av, ZX_rem(z, T));
    1020             : }
    1021             : 
    1022             : static GEN
    1023       12684 : ZC_mod2BIL_ZXQC(GEN x, long bs, GEN T)
    1024             : {
    1025       12684 :   long i, lx = lg(x);
    1026       12684 :   GEN A = cgetg(lx, t_COL);
    1027       73808 :   for (i=1; i<lx; i++) gel(A,i) = Z_mod2BIL_ZXQ(gel(x,i), bs, T);
    1028       12684 :   return A;
    1029             : }
    1030             : 
    1031             : static GEN
    1032        5943 : ZM_mod2BIL_ZXQM(GEN x, long bs, GEN T)
    1033             : {
    1034        5943 :   long i, lx = lg(x);
    1035        5943 :   GEN A = cgetg(lx, t_MAT);
    1036       18627 :   for (i=1; i<lx; i++) gel(A,i) = ZC_mod2BIL_ZXQC(gel(x,i), bs, T);
    1037        5943 :   return A;
    1038             : }
    1039             : 
    1040             : GEN
    1041        5803 : ZXQM_mul(GEN x, GEN y, GEN T)
    1042             : {
    1043        5803 :   long d = degpol(T);
    1044             :   GEN z;
    1045        5803 :   pari_sp av = avma;
    1046        5803 :   if (d == 0)
    1047           0 :     z = ZM_mul(simplify_shallow(x),simplify_shallow(y));
    1048             :   else
    1049             :   {
    1050        5803 :     long e, N, ex = ZXM_expi(x), ey = ZXM_expi(y), n = lg(x)-1;
    1051        5803 :     e = ex + ey + expu(d) + expu(n) + 4;
    1052        5803 :     N = divsBIL(e)+1;
    1053        5803 :     z = ZM_mul(ZXM_eval2BIL(x,N), ZXM_eval2BIL(y,N));
    1054        5803 :     z = ZM_mod2BIL_ZXQM(z, N, T);
    1055             :   }
    1056        5803 :   return gerepileupto(av, z);
    1057             : }
    1058             : 
    1059             : GEN
    1060         140 : ZXQM_sqr(GEN x, GEN T)
    1061             : {
    1062         140 :   long d = degpol(T);
    1063             :   GEN z;
    1064         140 :   pari_sp av = avma;
    1065         140 :   if (d == 0)
    1066           0 :     z = ZM_sqr(simplify_shallow(x));
    1067             :   else
    1068             :   {
    1069         140 :     long ex = ZXM_expi(x), d = degpol(T), n = lg(x)-1;
    1070         140 :     long e = 2*ex + expu(d) + expu(n) + 4;
    1071         140 :     long N = divsBIL(e)+1;
    1072         140 :     z = ZM_sqr(ZXM_eval2BIL(x,N));
    1073         140 :     z = ZM_mod2BIL_ZXQM(z, N, T);
    1074             :   }
    1075         140 :   return gerepileupto(av, z);
    1076             : }
    1077             : 
    1078             : GEN
    1079        4389 : QXQM_mul(GEN x, GEN y, GEN T)
    1080             : {
    1081        4389 :   GEN dx, nx = Q_primitive_part(x, &dx);
    1082        4389 :   GEN dy, ny = Q_primitive_part(y, &dy);
    1083        4389 :   GEN z = ZXQM_mul(nx, ny, T);
    1084        4389 :   if (dx || dy)
    1085             :   {
    1086        4389 :     GEN d = dx ? dy ? gmul(dx, dy): dx : dy;
    1087        4389 :     if (!gequal1(d)) z = RgM_Rg_mul(z, d);
    1088             :   }
    1089        4389 :   return z;
    1090             : }
    1091             : 
    1092             : GEN
    1093           7 : QXQM_sqr(GEN x, GEN T)
    1094             : {
    1095           7 :   GEN dx, nx = Q_primitive_part(x, &dx);
    1096           7 :   GEN z = ZXQM_sqr(nx, T);
    1097           7 :   if (dx) z = RgM_Rg_mul(z, gsqr(dx));
    1098           7 :   return z;
    1099             : }
    1100             : 
    1101             : static GEN
    1102     1004988 : Z_mod2BIL_Fq(GEN x, long bs, GEN T, GEN p)
    1103             : {
    1104     1004988 :   pari_sp av = avma;
    1105     1004988 :   long v = get_FpX_var(T), d = 2*(get_FpX_degree(T)-1);
    1106     1004988 :   GEN z = Z_mod2BIL_ZX(x, bs, d, 0);
    1107     1004988 :   setvarn(z, v);
    1108     1004988 :   return gerepileupto(av, FpX_rem(z, T, p));
    1109             : }
    1110             : 
    1111             : static GEN
    1112      240757 : ZC_mod2BIL_FqC(GEN x, long bs, GEN T, GEN p)
    1113             : {
    1114      240757 :   long i, lx = lg(x);
    1115      240757 :   GEN A = cgetg(lx, t_COL);
    1116     1245745 :   for (i=1; i<lx; i++) gel(A,i) = Z_mod2BIL_Fq(gel(x,i), bs, T, p);
    1117      240757 :   return A;
    1118             : }
    1119             : 
    1120             : static GEN
    1121       62804 : ZM_mod2BIL_FqM(GEN x, long bs, GEN T, GEN p)
    1122             : {
    1123       62804 :   long i, lx = lg(x);
    1124       62804 :   GEN A = cgetg(lx, t_MAT);
    1125      303561 :   for (i=1; i<lx; i++) gel(A,i) = ZC_mod2BIL_FqC(gel(x,i), bs, T, p);
    1126       62804 :   return A;
    1127             : }
    1128             : 
    1129             : GEN
    1130       62804 : FqM_mul_Kronecker(GEN x, GEN y, GEN T, GEN p)
    1131             : {
    1132       62804 :   pari_sp av = avma;
    1133       62804 :   long ex = ZXM_expi(x), ey = ZXM_expi(y), d= degpol(T), n = lg(x)-1;
    1134       62804 :   long e = ex + ey + expu(d) + expu(n) + 4;
    1135       62804 :   long N = divsBIL(e)+1;
    1136       62804 :   GEN  z = ZM_mul(ZXM_eval2BIL(x,N), ZXM_eval2BIL(y,N));
    1137       62804 :   return gerepileupto(av, ZM_mod2BIL_FqM(z, N, T, p));
    1138             : }
    1139             : 
    1140             : /*******************************************************************/
    1141             : /*                                                                 */
    1142             : /*                                ZXX                              */
    1143             : /*                                                                 */
    1144             : /*******************************************************************/
    1145             : 
    1146             : void
    1147         896 : RgX_check_ZXX(GEN x, const char *s)
    1148             : {
    1149         896 :   long k = lg(x)-1;
    1150        9324 :   for ( ; k>1; k--) {
    1151        8428 :     GEN t = gel(x,k);
    1152        8428 :     switch(typ(t)) {
    1153        7490 :       case t_INT: break;
    1154         938 :       case t_POL: if (RgX_is_ZX(t)) break;
    1155             :       /* fall through */
    1156           0 :       default: pari_err_TYPE(stack_strcat(s, " not in Z[X,Y]"),x);
    1157             :     }
    1158             :   }
    1159         896 : }
    1160             : 
    1161             : /*Renormalize (in place) polynomial with t_INT or ZX coefficients.*/
    1162             : GEN
    1163   227837970 : ZXX_renormalize(GEN x, long lx)
    1164             : {
    1165             :   long i;
    1166   292340846 :   for (i = lx-1; i>1; i--)
    1167   279324410 :     if (signe(gel(x,i))) break;
    1168   227837970 :   stackdummy((pari_sp)(x + lg(x)), (pari_sp)(x + (i+1)));
    1169   228044430 :   setlg(x, i+1); setsigne(x, i!=1); return x;
    1170             : }
    1171             : 
    1172             : GEN
    1173       20629 : ZXX_evalx0(GEN y)
    1174             : {
    1175       20629 :   long i, l = lg(y);
    1176       20629 :   GEN z = cgetg(l,t_POL); z[1] = y[1];
    1177      400260 :   for(i=2; i<l; i++)
    1178             :   {
    1179      379631 :     GEN yi = gel(y,i);
    1180      379631 :     gel(z,i) = typ(yi)==t_INT? yi: constant_coeff(yi);
    1181             :   }
    1182       20629 :   return ZX_renormalize(z,l);
    1183             : }
    1184             : 
    1185             : long
    1186        2723 : ZXX_max_lg(GEN x)
    1187             : {
    1188        2723 :   long i, prec = 0, lx = lg(x);
    1189       19425 :   for (i=2; i<lx; i++)
    1190             :   {
    1191       16702 :     GEN p = gel(x,i);
    1192       16702 :     long l = (typ(p) == t_INT)? lgefint(p): ZX_max_lg(p);
    1193       16702 :     if (l > prec) prec = l;
    1194             :   }
    1195        2723 :   return prec;
    1196             : }
    1197             : 
    1198             : GEN
    1199        8785 : ZXX_Z_mul(GEN y, GEN x)
    1200             : {
    1201        8785 :   long i, l = lg(y);
    1202        8785 :   GEN z = cgetg(l,t_POL); z[1] = y[1];
    1203      230587 :   for(i=2; i<l; i++)
    1204      221802 :     if(typ(gel(y,i))==t_INT)
    1205           0 :       gel(z,i) = mulii(gel(y,i),x);
    1206             :     else
    1207      221802 :       gel(z,i) = ZX_Z_mul(gel(y,i),x);
    1208        8785 :   return z;
    1209             : }
    1210             : 
    1211             : GEN
    1212           0 : ZXX_Z_add_shallow(GEN x, GEN y)
    1213             : {
    1214           0 :   long i, l = lg(x);
    1215             :   GEN z, a;
    1216           0 :   if (signe(x)==0) return scalarpol(y,varn(x));
    1217           0 :   z = cgetg(l,t_POL); z[1] = x[1];
    1218           0 :   a = gel(x,2);
    1219           0 :   gel(z, 2) = typ(a)==t_INT? addii(a,y): ZX_Z_add(a,y);
    1220           0 :   for(i=3; i<l; i++)
    1221           0 :     gel(z,i) = gel(x,i);
    1222           0 :   return z;
    1223             : }
    1224             : 
    1225             : GEN
    1226       55209 : ZXX_Z_divexact(GEN y, GEN x)
    1227             : {
    1228       55209 :   long i, l = lg(y);
    1229       55209 :   GEN z = cgetg(l,t_POL); z[1] = y[1];
    1230      463008 :   for(i=2; i<l; i++)
    1231      407799 :     if(typ(gel(y,i))==t_INT)
    1232        2996 :       gel(z,i) = diviiexact(gel(y,i),x);
    1233             :     else
    1234      404803 :       gel(z,i) = ZX_Z_divexact(gel(y,i),x);
    1235       55209 :   return z;
    1236             : }
    1237             : 
    1238             : /* Kronecker substitution, ZXX -> ZX:
    1239             :  * P(X,Y) = sum_{0<=i<lP} P_i(X) * Y^i, where deg P_i < n.
    1240             :  * Returns P(X,X^(2n-1)) */
    1241             : GEN
    1242     9579910 : RgXX_to_Kronecker_spec(GEN P, long lP, long n)
    1243             : {
    1244     9579910 :   long i, k, N = (n<<1) + 1;
    1245             :   GEN y;
    1246     9579910 :   if (!lP) return pol_0(0);
    1247     9562181 :   y = cgetg((N-2)*lP + 2, t_POL) + 2;
    1248    34510118 :   for (k=i=0; i<lP; i++)
    1249             :   {
    1250             :     long j;
    1251    34517877 :     GEN c = gel(P,i);
    1252    34517877 :     if (typ(c)!=t_POL)
    1253             :     {
    1254     2955983 :       gel(y,k++) = c;
    1255     2955983 :       j = 3;
    1256             :     }
    1257             :     else
    1258             :     {
    1259    31561894 :       long l = lg(c);
    1260    31561894 :       if (l-3 >= n)
    1261           0 :         pari_err_BUG("RgXX_to_Kronecker, P is not reduced mod Q");
    1262   192173887 :       for (j=2; j < l; j++) gel(y,k++) = gel(c,j);
    1263             :     }
    1264    34518211 :     if (i == lP-1) break;
    1265   149192009 :     for (   ; j < N; j++) gel(y,k++) = gen_0;
    1266             :   }
    1267     9558315 :   y-=2; setlg(y, k+2); y[1] = evalsigne(1); return y;
    1268             : }
    1269             : 
    1270             : /* shallow, n = deg(T) */
    1271             : GEN
    1272     1762450 : Kronecker_to_ZXX(GEN z, long n, long v)
    1273             : {
    1274     1762450 :   long i,j,lx,l, N = (n<<1)+1;
    1275             :   GEN x, t;
    1276     1762450 :   l = lg(z); lx = (l-2) / (N-2);
    1277     1762450 :   x = cgetg(lx+3,t_POL);
    1278     1762350 :   x[1] = z[1];
    1279     8899797 :   for (i=2; i<lx+2; i++)
    1280             :   {
    1281     7137321 :     t = cgetg(N,t_POL); t[1] = evalvarn(v);
    1282   106609115 :     for (j=2; j<N; j++) gel(t,j) = gel(z,j);
    1283     7136524 :     z += (N-2);
    1284     7136524 :     gel(x,i) = ZX_renormalize(t,N);
    1285             :   }
    1286     1762476 :   N = (l-2) % (N-2) + 2;
    1287     1762476 :   t = cgetg(N,t_POL); t[1] = evalvarn(v);
    1288     5118548 :   for (j=2; j<N; j++) gel(t,j) = gel(z,j);
    1289     1762734 :   gel(x,i) = ZX_renormalize(t,N);
    1290     1762636 :   return ZXX_renormalize(x, i+1);
    1291             : }
    1292             : 
    1293             : GEN
    1294     6565725 : RgXX_to_Kronecker(GEN P, long n)
    1295             : {
    1296     6565725 :   GEN z = RgXX_to_Kronecker_spec(P+2, lgpol(P), n);
    1297     6565698 :   setvarn(z,varn(P)); return z;
    1298             : }
    1299             : GEN
    1300     3017416 : ZXX_mul_Kronecker(GEN x, GEN y, long n)
    1301     3017416 : { return ZX_mul(RgXX_to_Kronecker(x,n), RgXX_to_Kronecker(y,n)); }
    1302             : 
    1303             : GEN
    1304        5194 : ZXX_sqr_Kronecker(GEN x, long n)
    1305        5194 : { return ZX_sqr(RgXX_to_Kronecker(x,n)); }
    1306             : 
    1307             : /* shallow, n = deg(T) */
    1308             : GEN
    1309       23443 : Kronecker_to_ZXQX(GEN z, GEN T)
    1310             : {
    1311       23443 :   long i,j,lx,l, N = (degpol(T)<<1)+1;
    1312             :   GEN x, t;
    1313       23443 :   l = lg(z); lx = (l-2) / (N-2);
    1314       23443 :   x = cgetg(lx+3,t_POL);
    1315       23443 :   x[1] = z[1];
    1316      234005 :   for (i=2; i<lx+2; i++)
    1317             :   {
    1318      210562 :     t = cgetg(N,t_POL); t[1] = T[1];
    1319     1888224 :     for (j=2; j<N; j++) gel(t,j) = gel(z,j);
    1320      210559 :     z += (N-2);
    1321      210559 :     gel(x,i) = ZX_rem(ZX_renormalize(t,N), T);
    1322             :   }
    1323       23443 :   N = (l-2) % (N-2) + 2;
    1324       23443 :   t = cgetg(N,t_POL); t[1] = T[1];
    1325       54467 :   for (j=2; j<N; j++) gel(t,j) = gel(z,j);
    1326       23443 :   gel(x,i) = ZX_rem(ZX_renormalize(t,N), T);
    1327       23443 :   return ZXX_renormalize(x, i+1);
    1328             : }
    1329             : 
    1330             : GEN
    1331        1792 : ZXQX_sqr(GEN x, GEN T)
    1332             : {
    1333        1792 :   pari_sp av = avma;
    1334        1792 :   long n = degpol(T);
    1335        1792 :   GEN z = ZXX_sqr_Kronecker(x, n);
    1336        1792 :   z = Kronecker_to_ZXQX(z, T);
    1337        1792 :   return gerepilecopy(av, z);
    1338             : }
    1339             : 
    1340             : GEN
    1341       21651 : ZXQX_mul(GEN x, GEN y, GEN T)
    1342             : {
    1343       21651 :   pari_sp av = avma;
    1344       21651 :   long n = degpol(T);
    1345       21651 :   GEN z = ZXX_mul_Kronecker(x, y, n);
    1346       21651 :   z = Kronecker_to_ZXQX(z, T);
    1347       21651 :   return gerepilecopy(av, z);
    1348             : }
    1349             : 
    1350             : GEN
    1351        8862 : ZXQX_ZXQ_mul(GEN P, GEN U, GEN T)
    1352             : {
    1353             :   long i, lP;
    1354             :   GEN res;
    1355        8862 :   res = cgetg_copy(P, &lP); res[1] = P[1];
    1356       46615 :   for(i=2; i<lP; i++)
    1357       37755 :     gel(res,i) = typ(gel(P,i))==t_POL? ZXQ_mul(U, gel(P,i), T)
    1358       37755 :                                      : gmul(U, gel(P,i));
    1359        8860 :   return ZXX_renormalize(res,lP);
    1360             : }
    1361             : 
    1362             : GEN
    1363     2936422 : QX_mul(GEN x, GEN y)
    1364             : {
    1365     2936422 :   GEN dx, nx = Q_primitive_part(x, &dx);
    1366     2936422 :   GEN dy, ny = Q_primitive_part(y, &dy);
    1367     2936421 :   GEN z = ZX_mul(nx, ny);
    1368     2936422 :   if (dx || dy)
    1369             :   {
    1370     2847855 :     GEN d = dx ? dy ? gmul(dx, dy): dx : dy;
    1371     2847855 :     return ZX_Q_mul(z, d);
    1372             :   } else
    1373       88567 :     return z;
    1374             : }
    1375             : 
    1376             : GEN
    1377       30155 : QX_sqr(GEN x)
    1378             : {
    1379       30155 :   GEN dx, nx = Q_primitive_part(x, &dx);
    1380       30155 :   GEN z = ZX_sqr(nx);
    1381       30155 :   if (dx)
    1382       19738 :     return ZX_Q_mul(z, gsqr(dx));
    1383             :   else
    1384       10417 :     return z;
    1385             : }
    1386             : 
    1387             : GEN
    1388     1656866 : QX_ZX_rem(GEN x, GEN y)
    1389             : {
    1390     1656866 :   pari_sp av = avma;
    1391     1656866 :   GEN d, nx = Q_primitive_part(x, &d);
    1392     1656866 :   GEN r = ZX_rem(nx, y);
    1393     1656866 :   if (d) r = ZX_Q_mul(r, d);
    1394     1656866 :   return gerepileupto(av, r);
    1395             : }
    1396             : 
    1397             : GEN
    1398       12789 : QXQX_mul(GEN x, GEN y, GEN T)
    1399             : {
    1400       12789 :   GEN dx, nx = Q_primitive_part(x, &dx);
    1401       12789 :   GEN dy, ny = Q_primitive_part(y, &dy);
    1402       12789 :   GEN z = ZXQX_mul(nx, ny, T);
    1403       12789 :   if (dx || dy)
    1404             :   {
    1405        8190 :     GEN d = dx ? dy ? gmul(dx, dy): dx : dy;
    1406        8190 :     return ZXX_Q_mul(z, d);
    1407             :   } else
    1408        4599 :     return z;
    1409             : }
    1410             : 
    1411             : GEN
    1412        1792 : QXQX_sqr(GEN x, GEN T)
    1413             : {
    1414        1792 :   GEN dx, nx = Q_primitive_part(x, &dx);
    1415        1792 :   GEN z = ZXQX_sqr(nx, T);
    1416        1792 :   if (dx)
    1417         588 :     return ZXX_Q_mul(z, gsqr(dx));
    1418             :   else
    1419        1204 :     return z;
    1420             : }
    1421             : 
    1422             : GEN
    1423         448 : QXQX_powers(GEN P, long n, GEN T)
    1424             : {
    1425         448 :   GEN v = cgetg(n+2, t_VEC);
    1426             :   long i;
    1427         448 :   gel(v, 1) = pol_1(varn(T));
    1428         448 :   if (n==0) return v;
    1429         448 :   gel(v, 2) = gcopy(P);
    1430        1512 :   for (i = 2; i <= n; i++) gel(v,i+1) = QXQX_mul(P, gel(v,i), T);
    1431         448 :   return v;
    1432             : }
    1433             : 
    1434             : GEN
    1435        3458 : QXQX_QXQ_mul(GEN P, GEN U, GEN T)
    1436             : {
    1437             :   long i, lP;
    1438             :   GEN res;
    1439        3458 :   res = cgetg_copy(P, &lP); res[1] = P[1];
    1440       22687 :   for(i=2; i<lP; i++)
    1441       19229 :     gel(res,i) = typ(gel(P,i))==t_POL? QXQ_mul(U, gel(P,i), T)
    1442       19229 :                                      : gmul(U, gel(P,i));
    1443        3458 :   return ZXX_renormalize(res,lP);
    1444             : }

Generated by: LCOV version 1.13