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 - FpXX.c (source / functions) Hit Total Coverage
Test: PARI/GP v2.12.0 lcov report (development 23036-b751c0af5) Lines: 830 983 84.4 %
Date: 2018-09-26 05:46:06 Functions: 101 114 88.6 %
Legend: Lines: hit not hit

          Line data    Source code
       1             : /* Copyright (C) 2012  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             : /* Not so fast arithmetic with polynomials over FpX */
      18             : 
      19             : /*******************************************************************/
      20             : /*                                                                 */
      21             : /*                             FpXX                                */
      22             : /*                                                                 */
      23             : /*******************************************************************/
      24             : /*Polynomials whose coefficients are either polynomials or integers*/
      25             : 
      26             : static ulong
      27      215839 : to_FlxqX(GEN P, GEN Q, GEN T, GEN p, GEN *pt_P, GEN *pt_Q, GEN *pt_T)
      28             : {
      29      215839 :   ulong pp = uel(p,2);
      30      215839 :   long v = get_FpX_var(T);
      31      215839 :   *pt_P = ZXX_to_FlxX(P, pp, v);
      32      215839 :   if (pt_Q) *pt_Q = ZXX_to_FlxX(Q, pp, v);
      33      215839 :   *pt_T = ZXT_to_FlxT(T, pp);
      34      215839 :   return pp;
      35             : }
      36             : 
      37             : static GEN
      38          81 : ZXX_copy(GEN a) { return gcopy(a); }
      39             : 
      40             : GEN
      41       32606 : FpXX_red(GEN z, GEN p)
      42             : {
      43             :   GEN res;
      44       32606 :   long i, l = lg(z);
      45       32606 :   res = cgetg(l,t_POL); res[1] = z[1];
      46      212247 :   for (i=2; i<l; i++)
      47             :   {
      48      179641 :     GEN zi = gel(z,i), c;
      49      179641 :     if (typ(zi)==t_INT)
      50        1799 :       c = modii(zi,p);
      51             :     else
      52             :     {
      53      177842 :       pari_sp av = avma;
      54      177842 :       c = FpX_red(zi,p);
      55      177842 :       switch(lg(c)) {
      56          14 :         case 2: set_avma(av); c = gen_0; break;
      57       16114 :         case 3: c = gerepilecopy(av, gel(c,2)); break;
      58             :       }
      59             :     }
      60      179641 :     gel(res,i) = c;
      61             :   }
      62       32606 :   return FpXX_renormalize(res,lg(res));
      63             : }
      64             : GEN
      65      397520 : FpXX_add(GEN x, GEN y, GEN p)
      66             : {
      67             :   long i,lz;
      68             :   GEN z;
      69      397520 :   long lx=lg(x);
      70      397520 :   long ly=lg(y);
      71      397520 :   if (ly>lx) swapspec(x,y, lx,ly);
      72      397520 :   lz = lx; z = cgetg(lz, t_POL); z[1]=x[1];
      73      397520 :   for (i=2; i<ly; i++) gel(z,i) = Fq_add(gel(x,i), gel(y,i), NULL, p);
      74      397520 :   for (   ; i<lx; i++) gel(z,i) = gcopy(gel(x,i));
      75      397520 :   return FpXX_renormalize(z, lz);
      76             : }
      77             : GEN
      78       13294 : FpXX_sub(GEN x, GEN y, GEN p)
      79             : {
      80             :   long i,lz;
      81             :   GEN z;
      82       13294 :   long lx=lg(x);
      83       13294 :   long ly=lg(y);
      84       13294 :   if (ly <= lx)
      85             :   {
      86       11579 :     lz = lx; z = cgetg(lz, t_POL); z[1]=x[1];
      87       11579 :     for (i=2; i<ly; i++) gel(z,i) = Fq_sub(gel(x,i), gel(y,i), NULL, p);
      88       11579 :     for (   ; i<lx; i++) gel(z,i) = gcopy(gel(x,i));
      89             :   }
      90             :   else
      91             :   {
      92        1715 :     lz = ly; z = cgetg(lz, t_POL); z[1]=x[1];
      93        1715 :     for (i=2; i<lx; i++) gel(z,i) = Fq_sub(gel(x,i), gel(y,i), NULL, p);
      94        1715 :     for (   ; i<ly; i++) gel(z,i) = Fq_neg(gel(y,i), NULL, p);
      95             :   }
      96       13294 :   return FpXX_renormalize(z, lz);
      97             : }
      98             : 
      99             : static GEN
     100      118022 : FpXX_subspec(GEN x, GEN y, GEN p, long nx, long ny)
     101             : {
     102             :   long i,lz;
     103             :   GEN z;
     104      118022 :   if (ny <= nx)
     105             :   {
     106      118022 :     lz = nx+2; z = cgetg(lz, t_POL)+2;
     107      118022 :     for (i=0; i<ny; i++) gel(z,i) = Fq_sub(gel(x,i), gel(y,i), NULL, p);
     108      118022 :     for (   ; i<nx; i++) gel(z,i) = gcopy(gel(x,i));
     109             :   }
     110             :   else
     111             :   {
     112           0 :     lz = ny+2; z = cgetg(lz, t_POL)+2;
     113           0 :     for (i=0; i<nx; i++) gel(z,i) = Fq_sub(gel(x,i), gel(y,i), NULL, p);
     114           0 :     for (   ; i<ny; i++) gel(z,i) = Fq_neg(gel(y,i), NULL, p);
     115             :   }
     116      118022 :   return FpXX_renormalize(z-2, lz);
     117             : }
     118             : 
     119             : GEN
     120         471 : FpXX_neg(GEN x, GEN p)
     121             : {
     122         471 :   long i, lx = lg(x);
     123         471 :   GEN y = cgetg(lx,t_POL);
     124         471 :   y[1] = x[1];
     125         471 :   for(i=2; i<lx; i++) gel(y,i) = Fq_neg(gel(x,i), NULL, p);
     126         471 :   return FpXX_renormalize(y, lx);
     127             : }
     128             : 
     129             : GEN
     130       64248 : FpXX_Fp_mul(GEN P, GEN u, GEN p)
     131             : {
     132             :   long i, lP;
     133       64248 :   GEN res = cgetg_copy(P, &lP); res[1] = P[1];
     134      473492 :   for(i=2; i<lP; i++)
     135             :   {
     136      409244 :     GEN x = gel(P,i);
     137      409244 :     gel(res,i) = typ(x)==t_INT? Fp_mul(x,u,p): FpX_Fp_mul(x,u,p);
     138             :   }
     139       64248 :   return FpXX_renormalize(res,lP);
     140             : }
     141             : 
     142             : GEN
     143       50368 : FpXX_mulu(GEN P, ulong u, GEN p)
     144             : {
     145             :   long i, lP;
     146       50368 :   GEN res = cgetg_copy(P, &lP); res[1] = P[1];
     147      269214 :   for(i=2; i<lP; i++)
     148             :   {
     149      218846 :     GEN x = gel(P,i);
     150      218846 :     gel(res,i) = typ(x)==t_INT? Fp_mulu(x,u,p): FpX_mulu(x,u,p);
     151             :   }
     152       50368 :   return FpXX_renormalize(res,lP);
     153             : }
     154             : 
     155             : GEN
     156           0 : FpXX_halve(GEN P, GEN p)
     157             : {
     158             :   long i, lP;
     159           0 :   GEN res = cgetg_copy(P, &lP); res[1] = P[1];
     160           0 :   for(i=2; i<lP; i++)
     161             :   {
     162           0 :     GEN x = gel(P,i);
     163           0 :     gel(res,i) = typ(x)==t_INT? Fp_halve(x,p): FpX_halve(x,p);
     164             :   }
     165           0 :   return FpXX_renormalize(res,lP);
     166             : }
     167             : 
     168             : GEN
     169        9502 : FpXX_deriv(GEN P, GEN p)
     170             : {
     171        9502 :   long i, l = lg(P)-1;
     172             :   GEN res;
     173             : 
     174        9502 :   if (l < 3) return pol_0(varn(P));
     175        9285 :   res = cgetg(l, t_POL);
     176        9285 :   res[1] = P[1];
     177       50223 :   for (i=2; i<l ; i++)
     178             :   {
     179       40938 :     GEN x = gel(P,i+1);
     180       40938 :     gel(res,i) = typ(x)==t_INT? Fp_mulu(x,i-1,p): FpX_mulu(x,i-1,p);
     181             :   }
     182        9285 :   return FpXX_renormalize(res, l);
     183             : }
     184             : 
     185             : GEN
     186           0 : FpXX_integ(GEN P, GEN p)
     187             : {
     188           0 :   long i, l = lg(P);
     189             :   GEN res;
     190             : 
     191           0 :   if (l == 2) return pol_0(varn(P));
     192           0 :   res = cgetg(l+1, t_POL);
     193           0 :   res[1] = P[1];
     194           0 :   gel(res,2) = gen_0;
     195           0 :   for (i=3; i<=l ; i++)
     196             :   {
     197           0 :     GEN x = gel(P,i-1);
     198           0 :     GEN i1 = Fp_inv(utoi(i-2), p);
     199           0 :     gel(res,i) = typ(x)==t_INT? Fp_mul(x,i1,p): FpX_Fp_mul(x,i1,p);
     200             :   }
     201           0 :   return FpXX_renormalize(res, l+1);
     202             : }
     203             : 
     204             : /*******************************************************************/
     205             : /*                                                                 */
     206             : /*                             (Fp[X]/(Q))[Y]                      */
     207             : /*                                                                 */
     208             : /*******************************************************************/
     209             : 
     210             : static GEN
     211      392991 : get_FpXQX_red(GEN T, GEN *B)
     212             : {
     213      392991 :   if (typ(T)!=t_VEC) { *B=NULL; return T; }
     214       60120 :   *B = gel(T,1); return gel(T,2);
     215             : }
     216             : 
     217             : GEN
     218          16 : random_FpXQX(long d1, long v, GEN T, GEN p)
     219             : {
     220          16 :   long dT = get_FpX_degree(T), vT = get_FpX_var(T);
     221          16 :   long i, d = d1+2;
     222          16 :   GEN y = cgetg(d,t_POL); y[1] = evalsigne(1) | evalvarn(v);
     223          16 :   for (i=2; i<d; i++) gel(y,i) = random_FpX(dT, vT, p);
     224          16 :   return FpXQX_renormalize(y,d);
     225             : }
     226             : 
     227             : /*Not stack clean*/
     228             : GEN
     229      672349 : Kronecker_to_FpXQX(GEN Z, GEN T, GEN p)
     230             : {
     231      672349 :   long i,j,lx,l, N = (get_FpX_degree(T)<<1) + 1;
     232      672349 :   GEN x, t = cgetg(N,t_POL), z = FpX_red(Z, p);
     233      672349 :   t[1] = evalvarn(get_FpX_var(T));
     234      672349 :   l = lg(z); lx = (l-2) / (N-2);
     235      672349 :   x = cgetg(lx+3,t_POL);
     236      672349 :   x[1] = z[1];
     237    14471747 :   for (i=2; i<lx+2; i++)
     238             :   {
     239    13799398 :     for (j=2; j<N; j++) gel(t,j) = gel(z,j);
     240    13799398 :     z += (N-2);
     241    13799398 :     gel(x,i) = FpX_rem(FpX_renormalize(t,N), T,p);
     242             :   }
     243      672349 :   N = (l-2) % (N-2) + 2;
     244      672349 :   for (j=2; j<N; j++) gel(t,j) = gel(z,j);
     245      672349 :   gel(x,i) = FpX_rem(FpX_renormalize(t,N), T,p);
     246      672349 :   return FpXQX_renormalize(x, i+1);
     247             : }
     248             : 
     249             : /* shallow, n = deg(T) */
     250             : GEN
     251      100282 : Kronecker_to_ZXX(GEN z, long n, long v)
     252             : {
     253      100282 :   long i,j,lx,l, N = (n<<1)+1;
     254             :   GEN x, t;
     255      100282 :   l = lg(z); lx = (l-2) / (N-2);
     256      100282 :   x = cgetg(lx+3,t_POL);
     257      100282 :   x[1] = z[1];
     258     1699138 :   for (i=2; i<lx+2; i++)
     259             :   {
     260     1598856 :     t = cgetg(N,t_POL); t[1] = evalvarn(v);
     261     1598856 :     for (j=2; j<N; j++) gel(t,j) = gel(z,j);
     262     1598856 :     z += (N-2);
     263     1598856 :     gel(x,i) = ZX_renormalize(t,N);
     264             :   }
     265      100282 :   N = (l-2) % (N-2) + 2;
     266      100282 :   t = cgetg(N,t_POL); t[1] = evalvarn(v);
     267      100282 :   for (j=2; j<N; j++) gel(t,j) = gel(z,j);
     268      100282 :   gel(x,i) = ZX_renormalize(t,N);
     269      100282 :   return ZXX_renormalize(x, i+1);
     270             : }
     271             : /* shallow */
     272             : GEN
     273      159012 : ZXX_mul_Kronecker(GEN x, GEN y, long n)
     274      159012 : { return ZX_mul(ZXX_to_Kronecker(x,n), ZXX_to_Kronecker(y,n)); }
     275             : 
     276             : GEN
     277        2086 : ZXX_sqr_Kronecker(GEN x, long n)
     278        2086 : { return ZX_sqr(ZXX_to_Kronecker(x,n)); }
     279             : 
     280             : GEN
     281      178513 : FpXQX_red(GEN z, GEN T, GEN p)
     282             : {
     283      178513 :   long i, l = lg(z);
     284      178513 :   GEN res = cgetg(l,t_POL); res[1] = z[1];
     285     1812013 :   for(i=2;i<l;i++)
     286     1633500 :     if (typ(gel(z,i)) == t_INT)
     287      107725 :       gel(res,i) = modii(gel(z,i),p);
     288             :     else
     289     1525775 :       gel(res,i) = FpXQ_red(gel(z,i),T,p);
     290      178513 :   return FpXQX_renormalize(res,l);
     291             : }
     292             : 
     293             : /* z in Z[X], return z * Mod(1,p), normalized*/
     294             : GEN
     295          77 : FpXQX_to_mod(GEN z, GEN T, GEN p)
     296             : {
     297          77 :   long i,l = lg(z);
     298          77 :   GEN x = cgetg(l, t_POL);
     299          77 :   x[1] = z[1];
     300          77 :   if (l == 2) return x;
     301          77 :   T = FpX_to_mod(T, p);
     302        1148 :   for (i=2; i<l; i++)
     303             :   {
     304        1071 :     GEN zi = gel(z,i);
     305        1071 :     gel(x,i) = typ(zi) == t_POL? mkpolmod(FpX_to_mod(zi, p), T): icopy(zi);
     306             :   }
     307          77 :   return normalizepol_lg(x,l);
     308             : }
     309             : 
     310             : static int
     311     1228875 : ZXX_is_ZX_spec(GEN a,long na)
     312             : {
     313             :   long i;
     314     1334658 :   for(i=0;i<na;i++)
     315     1320118 :     if(typ(gel(a,i))!=t_INT) return 0;
     316       14540 :   return 1;
     317             : }
     318             : 
     319             : static int
     320      141051 : ZXX_is_ZX(GEN a) { return ZXX_is_ZX_spec(a+2,lgpol(a)); }
     321             : 
     322             : static GEN
     323       45619 : FpXX_FpX_mulspec(GEN P, GEN U, GEN p, long v, long lU)
     324             : {
     325       45619 :   long i, lP =lg(P);
     326             :   GEN res;
     327       45619 :   res = cgetg(lP, t_POL); res[1] = P[1];
     328     1788324 :   for(i=2; i<lP; i++)
     329             :   {
     330     1742705 :     GEN Pi = gel(P,i);
     331     3478829 :     gel(res,i) = typ(Pi)==t_INT? FpX_Fp_mulspec(U, Pi, p, lU):
     332     1736124 :                                  FpX_mulspec(U, Pi+2, p, lU, lgpol(Pi));
     333     1742705 :     setvarn(gel(res,i),v);
     334             :   }
     335       45619 :   return FpXQX_renormalize(res,lP);
     336             : }
     337             : 
     338             : GEN
     339       36636 : FpXX_FpX_mul(GEN P, GEN U, GEN p)
     340       36636 : { return FpXX_FpX_mulspec(P,U+2,p,varn(U),lgpol(U)); }
     341             : 
     342             : static GEN
     343        8983 : FpXY_FpY_mulspec(GEN x, GEN y, GEN T, GEN p, long lx, long ly)
     344             : {
     345        8983 :   pari_sp av = avma;
     346        8983 :   long v = fetch_var();
     347        8983 :   GEN z = RgXY_swapspec(x,get_FpX_degree(T)-1,v,lx);
     348        8983 :   z = FpXX_FpX_mulspec(z,y,p,v,ly);
     349        8983 :   z = RgXY_swapspec(z+2,lx+ly+3,get_FpX_var(T),lgpol(z));
     350        8983 :   (void)delete_var(); return gerepilecopy(av,z);
     351             : }
     352             : 
     353             : static GEN
     354      543912 : FpXQX_mulspec(GEN x, GEN y, GEN T, GEN p, long lx, long ly)
     355             : {
     356      543912 :   pari_sp av = avma;
     357             :   GEN z, kx, ky;
     358             :   long n;
     359      543912 :   if (ZXX_is_ZX_spec(y,ly))
     360             :   {
     361        6457 :     if (ZXX_is_ZX_spec(x,lx))
     362        1926 :       return FpX_mulspec(x,y,p,lx,ly);
     363             :     else
     364        4531 :       return FpXY_FpY_mulspec(x,y,T,p,lx,ly);
     365      537455 :   } else if (ZXX_is_ZX_spec(x,lx))
     366        4452 :       return FpXY_FpY_mulspec(y,x,T,p,ly,lx);
     367      533003 :   n = get_FpX_degree(T);
     368      533003 :   kx = ZXX_to_Kronecker_spec(x, lx, n);
     369      533003 :   ky = ZXX_to_Kronecker_spec(y, ly, n);
     370      533003 :   z = Kronecker_to_FpXQX(ZX_mul(ky,kx), T, p);
     371      533003 :   return gerepileupto(av, z);
     372             : }
     373             : 
     374             : GEN
     375      305558 : FpXQX_mul(GEN x, GEN y, GEN T, GEN p)
     376             : {
     377      305558 :   GEN z = FpXQX_mulspec(x+2,y+2,T,p,lgpol(x),lgpol(y));
     378      305558 :   setvarn(z,varn(x)); return z;
     379             : }
     380             : 
     381             : GEN
     382      141019 : FpXQX_sqr(GEN x, GEN T, GEN p)
     383             : {
     384      141019 :   pari_sp av = avma;
     385             :   GEN z, kx;
     386      141019 :   if (ZXX_is_ZX(x)) return ZX_sqr(x);
     387      139314 :   kx= ZXX_to_Kronecker(x, get_FpX_degree(T));
     388      139314 :   z = Kronecker_to_FpXQX(ZX_sqr(kx), T, p);
     389      139314 :   return gerepileupto(av, z);
     390             : }
     391             : 
     392             : GEN
     393      240234 : FpXQX_FpXQ_mul(GEN P, GEN U, GEN T, GEN p)
     394             : {
     395             :   long i, lP;
     396             :   GEN res;
     397      240234 :   res = cgetg_copy(P, &lP); res[1] = P[1];
     398     1113248 :   for(i=2; i<lP; i++)
     399     1652358 :     gel(res,i) = typ(gel(P,i))==t_INT? FpX_Fp_mul(U, gel(P,i), p):
     400      779344 :                                        FpXQ_mul(U, gel(P,i), T,p);
     401      240234 :   return FpXQX_renormalize(res,lP);
     402             : }
     403             : 
     404             : /* x and y in Z[Y][X]. Assume T irreducible mod p */
     405             : static GEN
     406      125775 : FpXQX_divrem_basecase(GEN x, GEN y, GEN T, GEN p, GEN *pr)
     407             : {
     408             :   long vx, dx, dy, dy1, dz, i, j, sx, lr;
     409             :   pari_sp av0, av, tetpil;
     410             :   GEN z,p1,rem,lead;
     411             : 
     412      125775 :   if (!signe(y)) pari_err_INV("FpX_divrem",y);
     413      125775 :   vx=varn(x); dy=degpol(y); dx=degpol(x);
     414      125775 :   if (dx < dy)
     415             :   {
     416         143 :     if (pr)
     417             :     {
     418         135 :       av0 = avma; x = FpXQX_red(x, T, p);
     419         135 :       if (pr == ONLY_DIVIDES) { set_avma(av0); return signe(x)? NULL: pol_0(vx); }
     420         135 :       if (pr == ONLY_REM) return x;
     421         135 :       *pr = x;
     422             :     }
     423         143 :     return pol_0(vx);
     424             :   }
     425      125632 :   lead = leading_coeff(y);
     426      125632 :   if (!dy) /* y is constant */
     427             :   {
     428         828 :     if (pr && pr != ONLY_DIVIDES)
     429             :     {
     430         570 :       if (pr == ONLY_REM) return pol_0(vx);
     431           4 :       *pr = pol_0(vx);
     432             :     }
     433         262 :     if (gequal1(lead)) return FpXQX_red(x,T,p);
     434         246 :     av0 = avma; x = FqX_Fq_mul(x, Fq_inv(lead, T,p), T,p);
     435         246 :     return gerepileupto(av0,x);
     436             :   }
     437      124804 :   av0 = avma; dz = dx-dy;
     438      124804 :   lead = gequal1(lead)? NULL: gclone(Fq_inv(lead,T,p));
     439      124804 :   set_avma(av0);
     440      124804 :   z = cgetg(dz+3,t_POL); z[1] = x[1];
     441      124804 :   x += 2; y += 2; z += 2;
     442      124804 :   for (dy1=dy-1; dy1>=0 && !signe(gel(y, dy1)); dy1--);
     443             : 
     444      124804 :   p1 = gel(x,dx); av = avma;
     445      124804 :   gel(z,dz) = lead? gerepileupto(av, Fq_mul(p1,lead, T, p)): gcopy(p1);
     446      400047 :   for (i=dx-1; i>=dy; i--)
     447             :   {
     448      275243 :     av=avma; p1=gel(x,i);
     449      923742 :     for (j=i-dy1; j<=i && j<=dz; j++)
     450      648499 :       p1 = Fq_sub(p1, Fq_mul(gel(z,j),gel(y,i-j),NULL,p),NULL,p);
     451      275243 :     if (lead) p1 = Fq_mul(p1, lead, NULL,p);
     452      275243 :     tetpil=avma; gel(z,i-dy) = gerepile(av,tetpil,Fq_red(p1,T,p));
     453             :   }
     454      124804 :   if (!pr) { if (lead) gunclone(lead); return z-2; }
     455             : 
     456      121703 :   rem = (GEN)avma; av = (pari_sp)new_chunk(dx+3);
     457      135764 :   for (sx=0; ; i--)
     458             :   {
     459      149825 :     p1 = gel(x,i);
     460      531770 :     for (j=maxss(0,i-dy1); j<=i && j<=dz; j++)
     461      396006 :       p1 = Fq_sub(p1, Fq_mul(gel(z,j),gel(y,i-j),NULL,p),NULL,p);
     462      135764 :     tetpil=avma; p1 = Fq_red(p1, T, p); if (signe(p1)) { sx = 1; break; }
     463       15849 :     if (!i) break;
     464       14061 :     set_avma(av);
     465             :   }
     466      121703 :   if (pr == ONLY_DIVIDES)
     467             :   {
     468           0 :     if (lead) gunclone(lead);
     469           0 :     if (sx) return gc_NULL(av0);
     470           0 :     avma = (pari_sp)rem; return z-2;
     471             :   }
     472      121703 :   lr=i+3; rem -= lr;
     473      121703 :   rem[0] = evaltyp(t_POL) | evallg(lr);
     474      121703 :   rem[1] = z[-1];
     475      121703 :   p1 = gerepile((pari_sp)rem,tetpil,p1);
     476      121703 :   rem += 2; gel(rem,i) = p1;
     477      825845 :   for (i--; i>=0; i--)
     478             :   {
     479      704142 :     av=avma; p1 = gel(x,i);
     480     2285159 :     for (j=maxss(0,i-dy1); j<=i && j<=dz; j++)
     481     1581017 :       p1 = Fq_sub(p1, Fq_mul(gel(z,j),gel(y,i-j), NULL,p), NULL,p);
     482      704142 :     tetpil=avma; gel(rem,i) = gerepile(av,tetpil, Fq_red(p1, T, p));
     483             :   }
     484      121703 :   rem -= 2;
     485      121703 :   if (lead) gunclone(lead);
     486      121703 :   if (!sx) (void)FpXQX_renormalize(rem, lr);
     487      121703 :   if (pr == ONLY_REM) return gerepileupto(av0,rem);
     488        7176 :   *pr = rem; return z-2;
     489             : }
     490             : 
     491             : static GEN
     492          76 : FpXQX_halfgcd_basecase(GEN a, GEN b, GEN T, GEN p)
     493             : {
     494          76 :   pari_sp av=avma;
     495             :   GEN u,u1,v,v1;
     496          76 :   long vx = varn(a);
     497          76 :   long n = lgpol(a)>>1;
     498          76 :   u1 = v = pol_0(vx);
     499          76 :   u = v1 = pol_1(vx);
     500         903 :   while (lgpol(b)>n)
     501             :   {
     502         751 :     GEN r, q = FpXQX_divrem(a,b, T, p, &r);
     503         751 :     a = b; b = r; swap(u,u1); swap(v,v1);
     504         751 :     u1 = FpXX_sub(u1, FpXQX_mul(u, q, T, p), p);
     505         751 :     v1 = FpXX_sub(v1, FpXQX_mul(v, q ,T, p), p);
     506         751 :     if (gc_needed(av,2))
     507             :     {
     508           0 :       if (DEBUGMEM>1) pari_warn(warnmem,"FpXQX_halfgcd (d = %ld)",degpol(b));
     509           0 :       gerepileall(av,6, &a,&b,&u1,&v1,&u,&v);
     510             :     }
     511             :   }
     512          76 :   return gerepilecopy(av, mkmat2(mkcol2(u,u1), mkcol2(v,v1)));
     513             : }
     514             : static GEN
     515         136 : FpXQX_addmulmul(GEN u, GEN v, GEN x, GEN y, GEN T, GEN p)
     516             : {
     517         136 :   return FpXX_add(FpXQX_mul(u, x, T, p),FpXQX_mul(v, y, T, p), p);
     518             : }
     519             : 
     520             : static GEN
     521          68 : FpXQXM_FpXQX_mul2(GEN M, GEN x, GEN y, GEN T, GEN p)
     522             : {
     523          68 :   GEN res = cgetg(3, t_COL);
     524          68 :   gel(res, 1) = FpXQX_addmulmul(gcoeff(M,1,1), gcoeff(M,1,2), x, y, T, p);
     525          68 :   gel(res, 2) = FpXQX_addmulmul(gcoeff(M,2,1), gcoeff(M,2,2), x, y, T, p);
     526          68 :   return res;
     527             : }
     528             : 
     529             : static GEN
     530          56 : FpXQXM_mul2(GEN A, GEN B, GEN T, GEN p)
     531             : {
     532          56 :   GEN A11=gcoeff(A,1,1),A12=gcoeff(A,1,2), B11=gcoeff(B,1,1),B12=gcoeff(B,1,2);
     533          56 :   GEN A21=gcoeff(A,2,1),A22=gcoeff(A,2,2), B21=gcoeff(B,2,1),B22=gcoeff(B,2,2);
     534          56 :   GEN M1 = FpXQX_mul(FpXX_add(A11,A22, p), FpXX_add(B11,B22, p), T, p);
     535          56 :   GEN M2 = FpXQX_mul(FpXX_add(A21,A22, p), B11, T, p);
     536          56 :   GEN M3 = FpXQX_mul(A11, FpXX_sub(B12,B22, p), T, p);
     537          56 :   GEN M4 = FpXQX_mul(A22, FpXX_sub(B21,B11, p), T, p);
     538          56 :   GEN M5 = FpXQX_mul(FpXX_add(A11,A12, p), B22, T, p);
     539          56 :   GEN M6 = FpXQX_mul(FpXX_sub(A21,A11, p), FpXX_add(B11,B12, p), T, p);
     540          56 :   GEN M7 = FpXQX_mul(FpXX_sub(A12,A22, p), FpXX_add(B21,B22, p), T, p);
     541          56 :   GEN T1 = FpXX_add(M1,M4, p), T2 = FpXX_sub(M7,M5, p);
     542          56 :   GEN T3 = FpXX_sub(M1,M2, p), T4 = FpXX_add(M3,M6, p);
     543          56 :   retmkmat2(mkcol2(FpXX_add(T1,T2, p), FpXX_add(M2,M4, p)),
     544             :             mkcol2(FpXX_add(M3,M5, p), FpXX_add(T3,T4, p)));
     545             : }
     546             : /* Return [0,1;1,-q]*M */
     547             : static GEN
     548          56 : FpXQX_FpXQXM_qmul(GEN q, GEN M, GEN T, GEN p)
     549             : {
     550          56 :   GEN u, v, res = cgetg(3, t_MAT);
     551          56 :   u = FpXX_sub(gcoeff(M,1,1), FpXQX_mul(gcoeff(M,2,1), q, T, p), p);
     552          56 :   gel(res,1) = mkcol2(gcoeff(M,2,1), u);
     553          56 :   v = FpXX_sub(gcoeff(M,1,2), FpXQX_mul(gcoeff(M,2,2), q, T, p), p);
     554          56 :   gel(res,2) = mkcol2(gcoeff(M,2,2), v);
     555          56 :   return res;
     556             : }
     557             : 
     558             : static GEN
     559           0 : matid2_FpXQXM(long v)
     560             : {
     561           0 :   retmkmat2(mkcol2(pol_1(v),pol_0(v)),
     562             :             mkcol2(pol_0(v),pol_1(v)));
     563             : }
     564             : 
     565             : static GEN
     566          56 : FpXQX_halfgcd_split(GEN x, GEN y, GEN T, GEN p)
     567             : {
     568          56 :   pari_sp av=avma;
     569             :   GEN R, S, V;
     570             :   GEN y1, r, q;
     571          56 :   long l = lgpol(x), n = l>>1, k;
     572          56 :   if (lgpol(y)<=n) return matid2_FpXQXM(varn(x));
     573          56 :   R = FpXQX_halfgcd(RgX_shift_shallow(x,-n),RgX_shift_shallow(y,-n), T, p);
     574          56 :   V = FpXQXM_FpXQX_mul2(R,x,y, T, p); y1 = gel(V,2);
     575          56 :   if (lgpol(y1)<=n) return gerepilecopy(av, R);
     576          56 :   q = FpXQX_divrem(gel(V,1), y1, T, p, &r);
     577          56 :   k = 2*n-degpol(y1);
     578          56 :   S = FpXQX_halfgcd(RgX_shift_shallow(y1,-k), RgX_shift_shallow(r,-k), T, p);
     579          56 :   return gerepileupto(av, FpXQXM_mul2(S,FpXQX_FpXQXM_qmul(q,R, T, p), T, p));
     580             : }
     581             : 
     582             : /* Return M in GL_2(Fp[X]) such that:
     583             : if [a',b']~=M*[a,b]~ then degpol(a')>= (lgpol(a)>>1) >degpol(b')
     584             : */
     585             : 
     586             : static GEN
     587         132 : FpXQX_halfgcd_i(GEN x, GEN y, GEN T, GEN p)
     588             : {
     589         132 :   if (lg(x)<=FpXQX_HALFGCD_LIMIT) return FpXQX_halfgcd_basecase(x, y, T, p);
     590          56 :   return FpXQX_halfgcd_split(x, y, T, p);
     591             : }
     592             : 
     593             : GEN
     594         132 : FpXQX_halfgcd(GEN x, GEN y, GEN T, GEN p)
     595             : {
     596         132 :   pari_sp av = avma;
     597             :   GEN M,q,r;
     598         132 :   if (lgefint(p)==3)
     599             :   {
     600           0 :     ulong pp = to_FlxqX(x, y, T, p, &x, &y, &T);
     601           0 :     M = FlxXM_to_ZXXM(FlxqX_halfgcd(x, y, T, pp));
     602             :   }
     603             :   else
     604             :   {
     605         132 :     if (!signe(x))
     606             :     {
     607           0 :       long v = varn(x);
     608           0 :       retmkmat2(mkcol2(pol_0(v),pol_1(v)),
     609             :                 mkcol2(pol_1(v),pol_0(v)));
     610             :     }
     611         132 :     if (degpol(y)<degpol(x)) return FpXQX_halfgcd_i(x, y, T, p);
     612           8 :     q = FpXQX_divrem(y, x, T, p, &r);
     613           8 :     M = FpXQX_halfgcd_i(x, r, T, p);
     614           8 :     gcoeff(M,1,1) = FpXX_sub(gcoeff(M,1,1), FpXQX_mul(q, gcoeff(M,1,2), T, p), p);
     615           8 :     gcoeff(M,2,1) = FpXX_sub(gcoeff(M,2,1), FpXQX_mul(q, gcoeff(M,2,2), T, p), p);
     616             :   }
     617           8 :   return gerepilecopy(av, M);
     618             : }
     619             : 
     620             : static GEN
     621        3290 : FpXQX_gcd_basecase(GEN a, GEN b, GEN T, GEN p)
     622             : {
     623        3290 :   pari_sp av = avma, av0=avma;
     624       22783 :   while (signe(b))
     625             :   {
     626             :     GEN c;
     627       16203 :     if (gc_needed(av0,2))
     628             :     {
     629           0 :       if (DEBUGMEM>1) pari_warn(warnmem,"FpXQX_gcd (d = %ld)",degpol(b));
     630           0 :       gerepileall(av0,2, &a,&b);
     631             :     }
     632       16203 :     av = avma; c = FpXQX_rem(a, b, T, p); a=b; b=c;
     633             :   }
     634        3290 :   set_avma(av); return a;
     635             : }
     636             : 
     637             : GEN
     638       12772 : FpXQX_gcd(GEN x, GEN y, GEN T, GEN p)
     639             : {
     640       12772 :   pari_sp av = avma;
     641       12772 :   if (lgefint(p) == 3)
     642             :   {
     643             :     GEN Pl, Ql, Tl, U;
     644        9405 :     ulong pp = to_FlxqX(x, y, T, p, &Pl, &Ql, &Tl);
     645        9405 :     U  = FlxqX_gcd(Pl, Ql, Tl, pp);
     646        9405 :     return gerepileupto(av, FlxX_to_ZXX(U));
     647             :   }
     648        3367 :   x = FpXQX_red(x, T, p);
     649        3367 :   y = FpXQX_red(y, T, p);
     650        3367 :   if (!signe(x)) return gerepileupto(av, y);
     651        6592 :   while (lg(y)>FpXQX_GCD_LIMIT)
     652             :   {
     653             :     GEN c;
     654          12 :     if (lgpol(y)<=(lgpol(x)>>1))
     655             :     {
     656           0 :       GEN r = FpXQX_rem(x, y, T, p);
     657           0 :       x = y; y = r;
     658             :     }
     659          12 :     c = FpXQXM_FpXQX_mul2(FpXQX_halfgcd(x,y, T, p), x, y, T, p);
     660          12 :     x = gel(c,1); y = gel(c,2);
     661          12 :     gerepileall(av,2,&x,&y);
     662             :   }
     663        3290 :   return gerepileupto(av, FpXQX_gcd_basecase(x, y, T, p));
     664             : }
     665             : 
     666             : static GEN
     667           0 : FpXQX_extgcd_basecase(GEN a, GEN b, GEN T, GEN p, GEN *ptu, GEN *ptv)
     668             : {
     669           0 :   pari_sp av=avma;
     670             :   GEN u,v,d,d1,v1;
     671           0 :   long vx = varn(a);
     672           0 :   d = a; d1 = b;
     673           0 :   v = pol_0(vx); v1 = pol_1(vx);
     674           0 :   while (signe(d1))
     675             :   {
     676           0 :     GEN r, q = FpXQX_divrem(d, d1, T, p, &r);
     677           0 :     v = FpXX_sub(v,FpXQX_mul(q,v1,T, p),p);
     678           0 :     u=v; v=v1; v1=u;
     679           0 :     u=r; d=d1; d1=u;
     680           0 :     if (gc_needed(av,2))
     681             :     {
     682           0 :       if (DEBUGMEM>1) pari_warn(warnmem,"FpXQX_extgcd (d = %ld)",degpol(d));
     683           0 :       gerepileall(av,5, &d,&d1,&u,&v,&v1);
     684             :     }
     685             :   }
     686           0 :   if (ptu) *ptu = FpXQX_div(FpXX_sub(d,FpXQX_mul(b,v, T, p), p), a, T, p);
     687           0 :   *ptv = v; return d;
     688             : }
     689             : 
     690             : static GEN
     691           0 : FpXQX_extgcd_halfgcd(GEN x, GEN y, GEN T, GEN p, GEN *ptu, GEN *ptv)
     692             : {
     693           0 :   pari_sp av=avma;
     694           0 :   GEN u,v,R = matid2_FpXQXM(varn(x));
     695           0 :   while (lg(y)>FpXQX_EXTGCD_LIMIT)
     696             :   {
     697             :     GEN M, c;
     698           0 :     if (lgpol(y)<=(lgpol(x)>>1))
     699             :     {
     700           0 :       GEN r, q = FpXQX_divrem(x, y, T, p, &r);
     701           0 :       x = y; y = r;
     702           0 :       R = FpXQX_FpXQXM_qmul(q, R, T, p);
     703             :     }
     704           0 :     M = FpXQX_halfgcd(x,y, T, p);
     705           0 :     c = FpXQXM_FpXQX_mul2(M, x,y, T, p);
     706           0 :     R = FpXQXM_mul2(M, R, T, p);
     707           0 :     x = gel(c,1); y = gel(c,2);
     708           0 :     gerepileall(av,3,&x,&y,&R);
     709             :   }
     710           0 :   y = FpXQX_extgcd_basecase(x,y, T, p, &u,&v);
     711           0 :   if (ptu) *ptu = FpXQX_addmulmul(u,v,gcoeff(R,1,1),gcoeff(R,2,1), T, p);
     712           0 :   *ptv = FpXQX_addmulmul(u,v,gcoeff(R,1,2),gcoeff(R,2,2), T, p);
     713           0 :   return y;
     714             : }
     715             : 
     716             : /* x and y in Z[Y][X], return lift(gcd(x mod T,p, y mod T,p)). Set u and v st
     717             :  * ux + vy = gcd (mod T,p) */
     718             : GEN
     719        6118 : FpXQX_extgcd(GEN x, GEN y, GEN T, GEN p, GEN *ptu, GEN *ptv)
     720             : {
     721             :   GEN d;
     722        6118 :   pari_sp ltop=avma;
     723        6118 :   if (lgefint(p) == 3)
     724             :   {
     725             :     GEN Pl, Ql, Tl, Dl;
     726        6118 :     ulong pp = to_FlxqX(x, y, T, p, &Pl, &Ql, &Tl);
     727        6118 :     Dl = FlxqX_extgcd(Pl, Ql, Tl, pp, ptu, ptv);
     728        6118 :     if (ptu) *ptu = FlxX_to_ZXX(*ptu);
     729        6118 :     *ptv = FlxX_to_ZXX(*ptv);
     730        6118 :     d = FlxX_to_ZXX(Dl);
     731             :   }
     732             :   else
     733             :   {
     734           0 :     x = FpXQX_red(x, T, p);
     735           0 :     y = FpXQX_red(y, T, p);
     736           0 :     if (lg(y)>FpXQX_EXTGCD_LIMIT)
     737           0 :       d = FpXQX_extgcd_halfgcd(x, y, T, p, ptu, ptv);
     738             :     else
     739           0 :       d = FpXQX_extgcd_basecase(x, y, T, p, ptu, ptv);
     740             :   }
     741        6118 :   gerepileall(ltop,ptu?3:2,&d,ptv,ptu);
     742        6118 :   return d;
     743             : }
     744             : 
     745             : GEN
     746          68 : FpXQX_dotproduct(GEN x, GEN y, GEN T, GEN p)
     747             : {
     748          68 :   long i, l = minss(lg(x), lg(y));
     749             :   pari_sp av;
     750             :   GEN c;
     751          68 :   if (l == 2) return gen_0;
     752          68 :   av = avma; c = gmul(gel(x,2),gel(y,2));
     753          68 :   for (i=3; i<l; i++) c = gadd(c, gmul(gel(x,i),gel(y,i)));
     754          68 :   return gerepileupto(av, Fq_red(c,T,p));
     755             : }
     756             : 
     757             : /***********************************************************************/
     758             : /**                                                                   **/
     759             : /**                       Barrett reduction                           **/
     760             : /**                                                                   **/
     761             : /***********************************************************************/
     762             : 
     763             : /* Return new lgpol */
     764             : static long
     765      239509 : ZXX_lgrenormalizespec(GEN x, long lx)
     766             : {
     767             :   long i;
     768      240019 :   for (i = lx-1; i>=0; i--)
     769      240019 :     if (signe(gel(x,i))) break;
     770      239509 :   return i+1;
     771             : }
     772             : 
     773             : static GEN
     774        2865 : FpXQX_invBarrett_basecase(GEN S, GEN T, GEN p)
     775             : {
     776        2865 :   long i, l=lg(S)-1, lr = l-1, k;
     777        2865 :   GEN r=cgetg(lr, t_POL); r[1]=S[1];
     778        2865 :   gel(r,2) = gen_1;
     779       22687 :   for (i=3; i<lr; i++)
     780             :   {
     781       19822 :     pari_sp av = avma;
     782       19822 :     GEN u = gel(S,l-i+2);
     783      185711 :     for (k=3; k<i; k++)
     784      165889 :       u = Fq_add(u, Fq_mul(gel(S,l-i+k), gel(r,k), NULL, p), NULL, p);
     785       19822 :     gel(r,i) = gerepileupto(av, Fq_red(Fq_neg(u, NULL, p), T, p));
     786             :   }
     787        2865 :   return FpXQX_renormalize(r,lr);
     788             : }
     789             : 
     790             : INLINE GEN
     791      236044 : FpXX_recipspec(GEN x, long l, long n)
     792             : {
     793      236044 :   return RgX_recipspec_shallow(x, l, n);
     794             : }
     795             : 
     796             : static GEN
     797         182 : FpXQX_invBarrett_Newton(GEN S, GEN T, GEN p)
     798             : {
     799         182 :   pari_sp av = avma;
     800         182 :   long nold, lx, lz, lq, l = degpol(S), i, lQ;
     801         182 :   GEN q, y, z, x = cgetg(l+2, t_POL) + 2;
     802         182 :   ulong mask = quadratic_prec_mask(l-2); /* assume l > 2 */
     803         182 :   for (i=0;i<l;i++) gel(x,i) = gen_0;
     804         182 :   q = RgX_recipspec_shallow(S+2,l+1,l+1); lQ = lgpol(q); q+=2;
     805             :   /* We work on _spec_ FpX's, all the l[xzq] below are lgpol's */
     806             : 
     807             :   /* initialize */
     808         182 :   gel(x,0) = Fq_inv(gel(q,0), T, p);
     809         182 :   if (lQ>1) gel(q,1) = Fq_red(gel(q,1), T, p);
     810         182 :   if (lQ>1 && signe(gel(q,1)))
     811         182 :   {
     812         182 :     GEN u = gel(q, 1);
     813         182 :     if (!gequal1(gel(x,0))) u = Fq_mul(u, Fq_sqr(gel(x,0), T, p), T, p);
     814         182 :     gel(x,1) = Fq_neg(u, T, p); lx = 2;
     815             :   }
     816             :   else
     817           0 :     lx = 1;
     818         182 :   nold = 1;
     819        1519 :   for (; mask > 1; )
     820             :   { /* set x -= x(x*q - 1) + O(t^(nnew + 1)), knowing x*q = 1 + O(t^(nold+1)) */
     821        1155 :     long i, lnew, nnew = nold << 1;
     822             : 
     823        1155 :     if (mask & 1) nnew--;
     824        1155 :     mask >>= 1;
     825             : 
     826        1155 :     lnew = nnew + 1;
     827        1155 :     lq = ZXX_lgrenormalizespec(q, minss(lQ,lnew));
     828        1155 :     z = FpXQX_mulspec(x, q, T, p, lx, lq); /* FIXME: high product */
     829        1155 :     lz = lgpol(z); if (lz > lnew) lz = lnew;
     830        1155 :     z += 2;
     831             :     /* subtract 1 [=>first nold words are 0]: renormalize so that z(0) != 0 */
     832        1155 :     for (i = nold; i < lz; i++) if (signe(gel(z,i))) break;
     833        1155 :     nold = nnew;
     834        1155 :     if (i >= lz) continue; /* z-1 = 0(t^(nnew + 1)) */
     835             : 
     836             :     /* z + i represents (x*q - 1) / t^i */
     837        1155 :     lz = ZXX_lgrenormalizespec (z+i, lz-i);
     838        1155 :     z = FpXQX_mulspec(x, z+i, T, p, lx, lz); /* FIXME: low product */
     839        1155 :     lz = lgpol(z); z += 2;
     840        1155 :     if (lz > lnew-i) lz = ZXX_lgrenormalizespec(z, lnew-i);
     841             : 
     842        1155 :     lx = lz+ i;
     843        1155 :     y  = x + i; /* x -= z * t^i, in place */
     844        1155 :     for (i = 0; i < lz; i++) gel(y,i) = Fq_neg(gel(z,i), T, p);
     845             :   }
     846         182 :   x -= 2; setlg(x, lx + 2); x[1] = S[1];
     847         182 :   return gerepilecopy(av, x);
     848             : }
     849             : 
     850             : /* 1/polrecip(S)+O(x^(deg(S)-1)) */
     851             : GEN
     852        3068 : FpXQX_invBarrett(GEN S, GEN T, GEN p)
     853             : {
     854        3068 :   pari_sp ltop = avma;
     855        3068 :   long l = lg(S);
     856             :   GEN r;
     857        3068 :   if (l<5) return pol_0(varn(S));
     858        3047 :   if (l<=FpXQX_INVBARRETT_LIMIT)
     859             :   {
     860        2865 :     GEN c = gel(S,l-1), ci=gen_1;
     861        2865 :     if (!gequal1(c))
     862             :     {
     863        1862 :       ci = Fq_inv(c, T, p);
     864        1862 :       S = FqX_Fq_mul(S, ci, T, p);
     865        1862 :       r = FpXQX_invBarrett_basecase(S, T, p);
     866        1862 :       r = FqX_Fq_mul(r, ci, T, p);
     867             :     } else
     868        1003 :       r = FpXQX_invBarrett_basecase(S, T, p);
     869             :   }
     870             :   else
     871         182 :     r = FpXQX_invBarrett_Newton(S, T, p);
     872        3047 :   return gerepileupto(ltop, r);
     873             : }
     874             : 
     875             : GEN
     876        9162 : FpXQX_get_red(GEN S, GEN T, GEN p)
     877             : {
     878        9162 :   if (typ(S)==t_POL && lg(S)>FpXQX_BARRETT_LIMIT)
     879         752 :     retmkvec2(FpXQX_invBarrett(S,T,p),S);
     880        8410 :   return S;
     881             : }
     882             : 
     883             : /* Compute x mod S where 2 <= degpol(S) <= l+1 <= 2*(degpol(S)-1)
     884             :  * and mg is the Barrett inverse of S. */
     885             : static GEN
     886      118022 : FpXQX_divrem_Barrettspec(GEN x, long l, GEN mg, GEN S, GEN T, GEN p, GEN *pr)
     887             : {
     888             :   GEN q, r;
     889      118022 :   long lt = degpol(S); /*We discard the leading term*/
     890             :   long ld, lm, lT, lmg;
     891      118022 :   ld = l-lt;
     892      118022 :   lm = minss(ld, lgpol(mg));
     893      118022 :   lT  = ZXX_lgrenormalizespec(S+2,lt);
     894      118022 :   lmg = ZXX_lgrenormalizespec(mg+2,lm);
     895      118022 :   q = FpXX_recipspec(x+lt,ld,ld);                 /* q = rec(x)     lq<=ld*/
     896      118022 :   q = FpXQX_mulspec(q+2,mg+2,T,p,lgpol(q),lmg);    /* q = rec(x) * mg lq<=ld+lm*/
     897      118022 :   q = FpXX_recipspec(q+2,minss(ld,lgpol(q)),ld);  /* q = rec (rec(x) * mg) lq<=ld*/
     898      118022 :   if (!pr) return q;
     899      118022 :   r = FpXQX_mulspec(q+2,S+2,T,p,lgpol(q),lT);      /* r = q*pol        lr<=ld+lt*/
     900      118022 :   r = FpXX_subspec(x,r+2,p,lt,minss(lt,lgpol(r))); /* r = x - r   lr<=lt */
     901      118022 :   if (pr == ONLY_REM) return r;
     902      118022 :   *pr = r; return q;
     903             : }
     904             : 
     905             : static GEN
     906       50088 : FpXQX_divrem_Barrett_noGC(GEN x, GEN mg, GEN S, GEN T, GEN p, GEN *pr)
     907             : {
     908       50088 :   long l = lgpol(x), lt = degpol(S), lm = 2*lt-1;
     909       50088 :   GEN q = NULL, r;
     910             :   long i;
     911       50088 :   if (l <= lt)
     912             :   {
     913           0 :     if (pr == ONLY_REM) return ZXX_copy(x);
     914           0 :     if (pr == ONLY_DIVIDES) return signe(x)? NULL: pol_0(varn(x));
     915           0 :     if (pr) *pr =  ZXX_copy(x);
     916           0 :     return pol_0(varn(x));
     917             :   }
     918       50088 :   if (lt <= 1)
     919          21 :     return FpXQX_divrem_basecase(x,S,T,p,pr);
     920       50067 :   if (pr != ONLY_REM && l>lm)
     921             :   {
     922        1341 :     q = cgetg(l-lt+2, t_POL);
     923        1341 :     for (i=0;i<l-lt;i++) gel(q+2,i) = gen_0;
     924             :   }
     925       50067 :   r = l>lm ? shallowcopy(x): x;
     926      168170 :   while (l>lm)
     927             :   {
     928       68036 :     GEN zr, zq = FpXQX_divrem_Barrettspec(r+2+l-lm,lm,mg,S,T,p,&zr);
     929       68036 :     long lz = lgpol(zr);
     930       68036 :     if (pr != ONLY_REM)
     931             :     {
     932       34893 :       long lq = lgpol(zq);
     933       34893 :       for(i=0; i<lq; i++) gel(q+2+l-lm,i) = gel(zq,2+i);
     934             :     }
     935       68036 :     for(i=0; i<lz; i++) gel(r+2+l-lm,i) = gel(zr,2+i);
     936       68036 :     l = l-lm+lz;
     937             :   }
     938       50067 :   if (pr != ONLY_REM)
     939             :   {
     940        1362 :     if (l > lt)
     941             :     {
     942        1281 :       GEN zq = FpXQX_divrem_Barrettspec(r+2,l,mg,S,T,p,&r);
     943        1281 :       if (!q) q = zq;
     944             :       else
     945             :       {
     946        1260 :         long lq = lgpol(zq);
     947        1260 :         for(i=0; i<lq; i++) gel(q+2,i) = gel(zq,2+i);
     948             :       }
     949             :     }
     950             :     else
     951          81 :     { setlg(r, l+2); r = ZXX_copy(r); }
     952             :   }
     953             :   else
     954             :   {
     955       48705 :     if (l > lt)
     956       48705 :       (void) FpXQX_divrem_Barrettspec(r+2,l,mg,S,T,p,&r);
     957             :     else
     958           0 :     { setlg(r, l+2); r = ZXX_copy(r); }
     959       48705 :     r[1] = x[1]; return FpXQX_renormalize(r, lg(r));
     960             :   }
     961        1362 :   if (pr) { r[1] = x[1]; r = FpXQX_renormalize(r, lg(r)); }
     962        1362 :   q[1] = x[1]; q = FpXQX_renormalize(q, lg(q));
     963        1362 :   if (pr == ONLY_DIVIDES) return signe(r)? NULL: q;
     964        1362 :   if (pr) *pr = r;
     965        1362 :   return q;
     966             : }
     967             : 
     968             : GEN
     969      202653 : FpXQX_divrem(GEN x, GEN S, GEN T, GEN p, GEN *pr)
     970             : {
     971      202653 :   GEN B, y = get_FpXQX_red(S, &B);
     972      202653 :   long dy = degpol(y), dx = degpol(x), d = dx-dy;
     973      202653 :   if (pr==ONLY_REM) return FpXQX_rem(x, y, T, p);
     974      202653 :   if (lgefint(p) == 3)
     975             :   {
     976             :     GEN a, b, t, z;
     977      190609 :     pari_sp tetpil, av = avma;
     978      190609 :     ulong pp = to_FlxqX(x, y, T, p, &a, &b, &t);
     979      190609 :     z = FlxqX_divrem(a, b, t, pp, pr);
     980      190609 :     if (pr == ONLY_DIVIDES && !z) return gc_NULL(av);
     981      190609 :     tetpil=avma;
     982      190609 :     z = FlxX_to_ZXX(z);
     983      190609 :     if (pr && pr != ONLY_DIVIDES && pr != ONLY_REM)
     984      189132 :       *pr = FlxX_to_ZXX(*pr);
     985        1477 :     else return gerepile(av, tetpil, z);
     986      189132 :     gerepileallsp(av,tetpil,2, pr, &z);
     987      189132 :     return z;
     988             :   }
     989       12044 :   if (!B && d+3 < FpXQX_DIVREM_BARRETT_LIMIT)
     990       10664 :     return FpXQX_divrem_basecase(x,y,T,p,pr);
     991             :   else
     992             :   {
     993        1380 :     pari_sp av=avma;
     994        1380 :     GEN mg = B? B: FpXQX_invBarrett(y, T, p);
     995        1380 :     GEN q = FpXQX_divrem_Barrett_noGC(x,mg,y,T,p,pr);
     996        1380 :     if (!q) return gc_NULL(av);
     997        1380 :     if (!pr || pr==ONLY_DIVIDES) return gerepilecopy(av, q);
     998         393 :     gerepileall(av,2,&q,pr);
     999         393 :     return q;
    1000             :   }
    1001             : }
    1002             : 
    1003             : GEN
    1004      190322 : FpXQX_rem(GEN x, GEN S, GEN T, GEN p)
    1005             : {
    1006      190322 :   GEN B, y = get_FpXQX_red(S, &B);
    1007      190322 :   long dy = degpol(y), dx = degpol(x), d = dx-dy;
    1008      190322 :   if (d < 0) return FpXQX_red(x, T, p);
    1009      173505 :   if (lgefint(p) == 3)
    1010             :   {
    1011        9707 :     pari_sp av = avma;
    1012             :     GEN a, b, t, z;
    1013        9707 :     ulong pp = to_FlxqX(x, y, T, p, &a, &b, &t);
    1014        9707 :     z = FlxqX_rem(a, b, t, pp);
    1015        9707 :     z = FlxX_to_ZXX(z);
    1016        9707 :     return gerepileupto(av, z);
    1017             :   }
    1018      163798 :   if (!B && d+3 < FpXQX_REM_BARRETT_LIMIT)
    1019      115090 :     return FpXQX_divrem_basecase(x,y, T, p, ONLY_REM);
    1020             :   else
    1021             :   {
    1022       48708 :     pari_sp av=avma;
    1023       48708 :     GEN mg = B? B: FpXQX_invBarrett(y, T, p);
    1024       48708 :     GEN r = FpXQX_divrem_Barrett_noGC(x, mg, y, T, p, ONLY_REM);
    1025       48708 :     return gerepileupto(av, r);
    1026             :   }
    1027             : }
    1028             : 
    1029             : /* x + y*z mod p */
    1030             : INLINE GEN
    1031       11872 : Fq_addmul(GEN x, GEN y, GEN z, GEN T, GEN p)
    1032             : {
    1033             :   pari_sp av;
    1034       11872 :   if (!signe(y) || !signe(z)) return Fq_red(x, T, p);
    1035       11760 :   if (!signe(x)) return Fq_mul(z,y, T, p);
    1036       11746 :   av = avma;
    1037       11746 :   return gerepileupto(av, Fq_add(x, Fq_mul(y, z, T, p), T, p));
    1038             : }
    1039             : 
    1040             : GEN
    1041        5936 : FpXQX_div_by_X_x(GEN a, GEN x, GEN T, GEN p, GEN *r)
    1042             : {
    1043        5936 :   long l = lg(a)-1, i;
    1044        5936 :   GEN z = cgetg(l, t_POL);
    1045        5936 :   z[1] = evalsigne(1) | evalvarn(0);
    1046        5936 :   gel(z, l-1) = gel(a,l);
    1047       17808 :   for (i=l-2; i>1; i--) /* z[i] = a[i+1] + x*z[i+1] */
    1048       11872 :     gel(z, i) = Fq_addmul(gel(a,i+1), x, gel(z,i+1), T, p);
    1049        5936 :   if (r) *r = Fq_addmul(gel(a,2), x, gel(z,2), T, p);
    1050        5936 :   return z;
    1051             : }
    1052             : 
    1053             : struct _FpXQXQ {
    1054             :   GEN T, S;
    1055             :   GEN p;
    1056             : };
    1057             : 
    1058        2766 : static GEN _FpXQX_mul(void *data, GEN a,GEN b)
    1059             : {
    1060        2766 :   struct _FpXQXQ *d=(struct _FpXQXQ*)data;
    1061        2766 :   return FpXQX_mul(a,b,d->T,d->p);
    1062             : }
    1063             : 
    1064        1925 : static GEN _FpXQX_sqr(void *data, GEN a)
    1065             : {
    1066        1925 :   struct _FpXQXQ *d=(struct _FpXQXQ*)data;
    1067        1925 :   return FpXQX_sqr(a, d->T, d->p);
    1068             : }
    1069             : 
    1070             : GEN
    1071          35 : FpXQX_powu(GEN x, ulong n, GEN T, GEN p)
    1072             : {
    1073             :   struct _FpXQXQ D;
    1074          35 :   if (n==0) return pol_1(varn(x));
    1075          35 :   D.T = T; D.p = p;
    1076          35 :   return gen_powu(x, n, (void *)&D, _FpXQX_sqr, _FpXQX_mul);
    1077             : }
    1078             : 
    1079             : GEN
    1080         158 : FpXQXV_prod(GEN V, GEN T, GEN p)
    1081             : {
    1082         158 :   if (lgefint(p) == 3)
    1083             :   {
    1084          12 :     pari_sp av = avma;
    1085          12 :     ulong pp = p[2];
    1086          12 :     GEN Tl = ZXT_to_FlxT(T, pp);
    1087          12 :     GEN Vl = ZXXV_to_FlxXV(V, pp, get_FpX_var(T));
    1088          12 :     Tl = FlxqXV_prod(Vl, Tl, pp);
    1089          12 :     return gerepileupto(av, FlxX_to_ZXX(Tl));
    1090             :   }
    1091             :   else
    1092             :   {
    1093             :     struct _FpXQXQ d;
    1094         146 :     d.T=T; d.p=p;
    1095         146 :     return gen_product(V, (void*)&d, &_FpXQX_mul);
    1096             :   }
    1097             : }
    1098             : 
    1099             : static GEN
    1100        9905 : _FpXQX_divrem(void * E, GEN x, GEN y, GEN *r)
    1101             : {
    1102        9905 :   struct _FpXQXQ *d = (struct _FpXQXQ *) E;
    1103        9905 :   return FpXQX_divrem(x, y, d->T, d->p, r);
    1104             : }
    1105             : 
    1106             : static GEN
    1107       34846 : _FpXQX_add(void * E, GEN x, GEN y)
    1108             : {
    1109       34846 :   struct _FpXQXQ *d = (struct _FpXQXQ *) E;
    1110       34846 :   return FpXX_add(x, y, d->p);
    1111             : }
    1112             : 
    1113             : static GEN
    1114        4088 : _FpXQX_sub(void * E, GEN x, GEN y) {
    1115        4088 :   struct _FpXQXQ *d = (struct _FpXQXQ*) E;
    1116        4088 :   return FpXX_sub(x,y, d->p);
    1117             : }
    1118             : 
    1119             : static struct bb_ring FpXQX_ring = { _FpXQX_add, _FpXQX_mul, _FpXQX_sqr };
    1120             : 
    1121             : GEN
    1122         588 : FpXQX_digits(GEN x, GEN B, GEN T, GEN p)
    1123             : {
    1124         588 :   pari_sp av = avma;
    1125         588 :   long d = degpol(B), n = (lgpol(x)+d-1)/d;
    1126             :   GEN z;
    1127             :   struct _FpXQXQ D;
    1128         588 :   D.T = T; D.p = p;
    1129         588 :   z = gen_digits(x, B, n, (void *)&D, &FpXQX_ring, _FpXQX_divrem);
    1130         588 :   return gerepileupto(av, z);
    1131             : }
    1132             : 
    1133             : GEN
    1134         189 : FpXQXV_FpXQX_fromdigits(GEN x, GEN B, GEN T, GEN p)
    1135             : {
    1136         189 :   pari_sp av = avma;
    1137             :   struct _FpXQXQ D;
    1138             :   GEN z;
    1139         189 :   D.T = T; D.p = p;
    1140         189 :   z = gen_fromdigits(x,B,(void *)&D, &FpXQX_ring);
    1141         189 :   return gerepileupto(av, z);
    1142             : }
    1143             : 
    1144             : /* Q an FpXY (t_POL with FpX coeffs), evaluate at X = x */
    1145             : GEN
    1146       48083 : FpXY_evalx(GEN Q, GEN x, GEN p)
    1147             : {
    1148       48083 :   long i, lb = lg(Q);
    1149             :   GEN z;
    1150       48083 :   z = cgetg(lb, t_POL); z[1] = Q[1];
    1151      436296 :   for (i=2; i<lb; i++)
    1152             :   {
    1153      388213 :     GEN q = gel(Q,i);
    1154      388213 :     gel(z,i) = typ(q) == t_INT? modii(q,p): FpX_eval(q, x, p);
    1155             :   }
    1156       48083 :   return FpX_renormalize(z, lb);
    1157             : }
    1158             : /* Q an FpXY, evaluate at Y = y */
    1159             : GEN
    1160       18834 : FpXY_evaly(GEN Q, GEN y, GEN p, long vx)
    1161             : {
    1162       18834 :   pari_sp av = avma;
    1163       18834 :   long i, lb = lg(Q);
    1164             :   GEN z;
    1165       18834 :   if (!signe(Q)) return pol_0(vx);
    1166       18806 :   if (lb == 3 || !signe(y)) {
    1167          84 :     z = gel(Q, 2);
    1168          84 :     return typ(z)==t_INT? scalar_ZX(z, vx): ZX_copy(z);
    1169             :   }
    1170       18722 :   z = gel(Q, lb-1);
    1171       18722 :   if (typ(z) == t_INT) z = scalar_ZX_shallow(z, vx);
    1172       18722 :   for (i=lb-2; i>=2; i--) z = Fq_add(gel(Q,i), FpX_Fp_mul(z, y, p), NULL, p);
    1173       18722 :   return gerepileupto(av, z);
    1174             : }
    1175             : /* Q an FpXY, evaluate at (X,Y) = (x,y) */
    1176             : GEN
    1177       13909 : FpXY_eval(GEN Q, GEN y, GEN x, GEN p)
    1178             : {
    1179       13909 :   pari_sp av = avma;
    1180       13909 :   return gerepileuptoint(av, FpX_eval(FpXY_evalx(Q, x, p), y, p));
    1181             : }
    1182             : 
    1183             : GEN
    1184        2825 : FpXY_FpXQV_evalx(GEN P, GEN x, GEN T, GEN p)
    1185             : {
    1186        2825 :   long i, lP = lg(P);
    1187        2825 :   GEN res = cgetg(lP,t_POL);
    1188        2825 :   res[1] = P[1];
    1189       32908 :   for(i=2; i<lP; i++)
    1190       60166 :     gel(res,i) = typ(gel(P,i))==t_INT? icopy(gel(P,i)):
    1191       30083 :                                        FpX_FpXQV_eval(gel(P,i), x, T, p);
    1192        2825 :   return FlxX_renormalize(res, lP);
    1193             : }
    1194             : 
    1195             : GEN
    1196         147 : FpXY_FpXQ_evalx(GEN P, GEN x, GEN T, GEN p)
    1197             : {
    1198         147 :   pari_sp av = avma;
    1199         147 :   long n = brent_kung_optpow(get_FpX_degree(T)-1,lgpol(P),1);
    1200         147 :   GEN xp = FpXQ_powers(x, n, T, p);
    1201         147 :   return gerepileupto(av, FpXY_FpXQV_evalx(P, xp, T, p));
    1202             : }
    1203             : 
    1204             : /*******************************************************************/
    1205             : /*                                                                 */
    1206             : /*                       (Fp[X]/T(X))[Y] / S(Y)                    */
    1207             : /*                                                                 */
    1208             : /*******************************************************************/
    1209             : 
    1210             : /*Preliminary implementation to speed up FpX_ffisom*/
    1211             : typedef struct {
    1212             :   GEN S, T, p;
    1213             : } FpXYQQ_muldata;
    1214             : 
    1215             : /* reduce x in Fp[X, Y] in the algebra Fp[X,Y]/ (S(X),T(Y)) */
    1216             : static GEN
    1217         476 : FpXYQQ_redswap(GEN x, GEN S, GEN T, GEN p)
    1218             : {
    1219         476 :   pari_sp ltop=avma;
    1220         476 :   long n = get_FpX_degree(S);
    1221         476 :   long m = get_FpX_degree(T);
    1222         476 :   long v = get_FpX_var(T);
    1223         476 :   GEN V = RgXY_swap(x,m,v);
    1224         476 :   V = FpXQX_red(V,S,p);
    1225         476 :   V = RgXY_swap(V,n,v);
    1226         476 :   return gerepilecopy(ltop,V);
    1227             : }
    1228             : static GEN
    1229         280 : FpXYQQ_sqr(void *data, GEN x)
    1230             : {
    1231         280 :   FpXYQQ_muldata *D = (FpXYQQ_muldata*)data;
    1232         280 :   return FpXYQQ_redswap(FpXQX_sqr(x, D->T, D->p),D->S,D->T,D->p);
    1233             : 
    1234             : }
    1235             : static GEN
    1236         196 : FpXYQQ_mul(void *data, GEN x, GEN y)
    1237             : {
    1238         196 :   FpXYQQ_muldata *D = (FpXYQQ_muldata*)data;
    1239         196 :   return FpXYQQ_redswap(FpXQX_mul(x,y, D->T, D->p),D->S,D->T,D->p);
    1240             : }
    1241             : 
    1242             : /* x in Z[X,Y], S in Z[X] over Fq = Z[Y]/(p,T); compute lift(x^n mod (S,T,p)) */
    1243             : GEN
    1244         182 : FpXYQQ_pow(GEN x, GEN n, GEN S, GEN T, GEN p)
    1245             : {
    1246         182 :   pari_sp av = avma;
    1247             :   FpXYQQ_muldata D;
    1248             :   GEN y;
    1249         182 :   if (lgefint(p) == 3)
    1250             :   {
    1251           0 :     ulong pp = to_FlxqX(x, NULL, T, p, &x, NULL, &T);
    1252           0 :     S = ZX_to_Flx(S, pp);
    1253           0 :     y = FlxX_to_ZXX( FlxYqq_pow(x, n, S, T, pp) );
    1254             :   }
    1255             :   else
    1256             :   {
    1257         182 :     D.S = S;
    1258         182 :     D.T = T;
    1259         182 :     D.p = p;
    1260         182 :     y = gen_pow(x, n, (void*)&D, &FpXYQQ_sqr, &FpXYQQ_mul);
    1261             :   }
    1262         182 :   return gerepileupto(av, y);
    1263             : }
    1264             : 
    1265             : GEN
    1266       32647 : FpXQXQ_mul(GEN x, GEN y, GEN S, GEN T, GEN p) {
    1267       32647 :   return FpXQX_rem(FpXQX_mul(x, y, T, p), S, T, p);
    1268             : }
    1269             : 
    1270             : GEN
    1271      137897 : FpXQXQ_sqr(GEN x, GEN S, GEN T, GEN p) {
    1272      137897 :   return FpXQX_rem(FpXQX_sqr(x, T, p), S, T, p);
    1273             : }
    1274             : 
    1275             : /* Inverse of x in Z/pZ[X]/(pol) or NULL if inverse doesn't exist
    1276             :  * return lift(1 / (x mod (p,pol))) */
    1277             : GEN
    1278          21 : FpXQXQ_invsafe(GEN x, GEN S, GEN T, GEN p)
    1279             : {
    1280          21 :   GEN V, z = FpXQX_extgcd(get_FpXQX_mod(S), x, T, p, NULL, &V);
    1281          21 :   if (degpol(z)) return NULL;
    1282          21 :   z = gel(z,2);
    1283          21 :   z = typ(z)==t_INT ? Fp_invsafe(z,p) : FpXQ_invsafe(z,T,p);
    1284          21 :   if (!z) return NULL;
    1285          21 :   return typ(z)==t_INT ? FpXX_Fp_mul(V, z, p): FpXQX_FpXQ_mul(V, z, T, p);
    1286             : }
    1287             : 
    1288             : GEN
    1289          21 : FpXQXQ_inv(GEN x, GEN S, GEN T,GEN p)
    1290             : {
    1291          21 :   pari_sp av = avma;
    1292          21 :   GEN U = FpXQXQ_invsafe(x, S, T, p);
    1293          21 :   if (!U) pari_err_INV("FpXQXQ_inv",x);
    1294          21 :   return gerepileupto(av, U);
    1295             : }
    1296             : 
    1297             : GEN
    1298           0 : FpXQXQ_div(GEN x,GEN y,GEN S, GEN T,GEN p)
    1299             : {
    1300           0 :   pari_sp av = avma;
    1301           0 :   return gerepileupto(av, FpXQXQ_mul(x, FpXQXQ_inv(y,S,T,p),S,T,p));
    1302             : }
    1303             : 
    1304             : static GEN
    1305       36659 : _FpXQXQ_cmul(void *data, GEN P, long a, GEN x) {
    1306       36659 :   struct _FpXQXQ *d = (struct _FpXQXQ*) data;
    1307       36659 :   GEN y = gel(P,a+2);
    1308       73295 :   return typ(y)==t_INT ? FpXX_Fp_mul(x,y, d->p):
    1309       36636 :                          FpXX_FpX_mul(x,y,d->p);
    1310             : }
    1311             : static GEN
    1312       12420 : _FpXQXQ_red(void *data, GEN x) {
    1313       12420 :   struct _FpXQXQ *d = (struct _FpXQXQ*) data;
    1314       12420 :   return FpXQX_red(x, d->T, d->p);
    1315             : }
    1316             : static GEN
    1317       30537 : _FpXQXQ_mul(void *data, GEN x, GEN y) {
    1318       30537 :   struct _FpXQXQ *d = (struct _FpXQXQ*) data;
    1319       30537 :   return FpXQXQ_mul(x,y, d->S,d->T, d->p);
    1320             : }
    1321             : static GEN
    1322      137897 : _FpXQXQ_sqr(void *data, GEN x) {
    1323      137897 :   struct _FpXQXQ *d = (struct _FpXQXQ*) data;
    1324      137897 :   return FpXQXQ_sqr(x, d->S,d->T, d->p);
    1325             : }
    1326             : 
    1327             : static GEN
    1328        9729 : _FpXQXQ_one(void *data) {
    1329        9729 :   struct _FpXQXQ *d = (struct _FpXQXQ*) data;
    1330        9729 :   return pol_1(get_FpXQX_var(d->S));
    1331             : }
    1332             : 
    1333             : static GEN
    1334          68 : _FpXQXQ_zero(void *data) {
    1335          68 :   struct _FpXQXQ *d = (struct _FpXQXQ*) data;
    1336          68 :   return pol_0(get_FpXQX_var(d->S));
    1337             : }
    1338             : 
    1339             : static struct bb_algebra FpXQXQ_algebra = { _FpXQXQ_red, _FpXQX_add,
    1340             :        _FpXQX_sub, _FpXQXQ_mul, _FpXQXQ_sqr, _FpXQXQ_one, _FpXQXQ_zero };
    1341             : 
    1342             : const struct bb_algebra *
    1343         628 : get_FpXQXQ_algebra(void **E, GEN S, GEN T, GEN p)
    1344             : {
    1345         628 :   GEN z = new_chunk(sizeof(struct _FpXQXQ));
    1346         628 :   struct _FpXQXQ *e = (struct _FpXQXQ *) z;
    1347         628 :   e->T = FpX_get_red(T, p);
    1348         628 :   e->S = FpXQX_get_red(S, e->T, p);
    1349         628 :   e->p  = p; *E = (void*)e;
    1350         628 :   return &FpXQXQ_algebra;
    1351             : }
    1352             : 
    1353             : static struct bb_algebra FpXQX_algebra = { _FpXQXQ_red, _FpXQX_add,
    1354             :        _FpXQX_sub, _FpXQX_mul, _FpXQX_sqr, _FpXQXQ_one, _FpXQXQ_zero };
    1355             : 
    1356             : const struct bb_algebra *
    1357          42 : get_FpXQX_algebra(void **E, GEN T, GEN p, long v)
    1358             : {
    1359          42 :   GEN z = new_chunk(sizeof(struct _FpXQXQ));
    1360          42 :   struct _FpXQXQ *e = (struct _FpXQXQ *) z;
    1361          42 :   e->T = FpX_get_red(T, p);
    1362          42 :   e->S = pol_x(v);
    1363          42 :   e->p  = p; *E = (void*)e;
    1364          42 :   return &FpXQX_algebra;
    1365             : }
    1366             : 
    1367             : /* x over Fq, return lift(x^n) mod S */
    1368             : GEN
    1369        1431 : FpXQXQ_pow(GEN x, GEN n, GEN S, GEN T, GEN p)
    1370             : {
    1371        1431 :   pari_sp ltop = avma;
    1372             :   GEN y;
    1373             :   struct _FpXQXQ D;
    1374        1431 :   long s = signe(n);
    1375        1431 :   if (!s) return pol_1(varn(x));
    1376        1431 :   if (is_pm1(n)) /* +/- 1 */
    1377           0 :     return (s < 0)? FpXQXQ_inv(x,S,T,p): ZXX_copy(x);
    1378        1431 :   if (lgefint(p) == 3)
    1379             :   {
    1380           0 :     ulong pp = to_FlxqX(x, S, T, p, &x, &S, &T);
    1381           0 :     GEN z = FlxqXQ_pow(x, n, S, T, pp);
    1382           0 :     y = FlxX_to_ZXX(z);
    1383             :   }
    1384             :   else
    1385             :   {
    1386        1431 :     T = FpX_get_red(T, p);
    1387        1431 :     S = FpXQX_get_red(S, T, p);
    1388        1431 :     D.S = S; D.T = T; D.p = p;
    1389        1431 :     if (s < 0) x = FpXQXQ_inv(x,S,T,p);
    1390        1431 :     y = gen_pow(x, n, (void*)&D,&_FpXQXQ_sqr,&_FpXQXQ_mul);
    1391             :   }
    1392        1431 :   return gerepileupto(ltop, y);
    1393             : }
    1394             : 
    1395             : /* generates the list of powers of x of degree 0,1,2,...,l*/
    1396             : GEN
    1397        1403 : FpXQXQ_powers(GEN x, long l, GEN S, GEN T, GEN p)
    1398             : {
    1399             :   struct _FpXQXQ D;
    1400        1403 :   int use_sqr = 2*degpol(x) >= get_FpXQX_degree(S);
    1401        1403 :   T = FpX_get_red(T, p);
    1402        1403 :   S = FpXQX_get_red(S, T, p);
    1403        1403 :   D.S = S; D.T = T; D.p = p;
    1404        1403 :   return gen_powers(x, l, use_sqr, (void*)&D, &_FpXQXQ_sqr, &_FpXQXQ_mul,&_FpXQXQ_one);
    1405             : }
    1406             : 
    1407             : /* Let v a linear form, return the linear form z->v(tau*z)
    1408             :    that is, v*(M_tau) */
    1409             : 
    1410             : INLINE GEN
    1411          48 : FpXQX_recipspec(GEN x, long l, long n)
    1412             : {
    1413          48 :   return RgX_recipspec_shallow(x, l, n);
    1414             : }
    1415             : 
    1416             : static GEN
    1417          16 : FpXQXQ_transmul_init(GEN tau, GEN S, GEN T, GEN p)
    1418             : {
    1419             :   GEN bht;
    1420          16 :   GEN h, Sp = get_FpXQX_red(S, &h);
    1421          16 :   long n = degpol(Sp), vT = varn(Sp);
    1422          16 :   GEN ft = FpXQX_recipspec(Sp+2, n+1, n+1);
    1423          16 :   GEN bt = FpXQX_recipspec(tau+2, lgpol(tau), n);
    1424          16 :   setvarn(ft, vT); setvarn(bt, vT);
    1425          16 :   if (h)
    1426           0 :     bht = FpXQXn_mul(bt, h, n-1, T, p);
    1427             :   else
    1428             :   {
    1429          16 :     GEN bh = FpXQX_div(RgX_shift_shallow(tau, n-1), S, T, p);
    1430          16 :     bht = FpXQX_recipspec(bh+2, lgpol(bh), n-1);
    1431          16 :     setvarn(bht, vT);
    1432             :   }
    1433          16 :   return mkvec3(bt, bht, ft);
    1434             : }
    1435             : 
    1436             : static GEN
    1437          40 : FpXQXQ_transmul(GEN tau, GEN a, long n, GEN T, GEN p)
    1438             : {
    1439          40 :   pari_sp ltop = avma;
    1440             :   GEN t1, t2, t3, vec;
    1441          40 :   GEN bt = gel(tau, 1), bht = gel(tau, 2), ft = gel(tau, 3);
    1442          40 :   if (signe(a)==0) return pol_0(varn(a));
    1443          40 :   t2 = RgX_shift_shallow(FpXQX_mul(bt, a, T, p),1-n);
    1444          40 :   if (signe(bht)==0) return gerepilecopy(ltop, t2);
    1445          32 :   t1 = RgX_shift_shallow(FpXQX_mul(ft, a, T, p),-n);
    1446          32 :   t3 = FpXQXn_mul(t1, bht, n-1, T, p);
    1447          32 :   vec = FpXX_sub(t2, RgX_shift_shallow(t3, 1), p);
    1448          32 :   return gerepileupto(ltop, vec);
    1449             : }
    1450             : 
    1451             : static GEN
    1452           8 : polxn_FpXX(long n, long v, long vT)
    1453             : {
    1454           8 :   long i, a = n+2;
    1455           8 :   GEN p = cgetg(a+1, t_POL);
    1456           8 :   p[1] = evalsigne(1)|evalvarn(v);
    1457           8 :   for (i = 2; i < a; i++) gel(p,i) = pol_0(vT);
    1458           8 :   gel(p,a) = pol_1(vT); return p;
    1459             : }
    1460             : 
    1461             : GEN
    1462           8 : FpXQXQ_minpoly(GEN x, GEN S, GEN T, GEN p)
    1463             : {
    1464           8 :   pari_sp ltop = avma;
    1465             :   long vS, vT, n;
    1466             :   GEN v_x, g, tau;
    1467           8 :   vS = get_FpXQX_var(S);
    1468           8 :   vT = get_FpX_var(T);
    1469           8 :   n = get_FpXQX_degree(S);
    1470           8 :   g = pol_1(vS);
    1471           8 :   tau = pol_1(vS);
    1472           8 :   S = FpXQX_get_red(S, T, p);
    1473           8 :   v_x = FpXQXQ_powers(x, usqrt(2*n), S, T, p);
    1474          24 :   while(signe(tau) != 0)
    1475             :   {
    1476             :     long i, j, m, k1;
    1477             :     GEN M, v, tr;
    1478             :     GEN g_prime, c;
    1479           8 :     if (degpol(g) == n) { tau = pol_1(vS); g = pol_1(vS); }
    1480           8 :     v = random_FpXQX(n, vS, T, p);
    1481           8 :     tr = FpXQXQ_transmul_init(tau, S, T, p);
    1482           8 :     v = FpXQXQ_transmul(tr, v, n, T, p);
    1483           8 :     m = 2*(n-degpol(g));
    1484           8 :     k1 = usqrt(m);
    1485           8 :     tr = FpXQXQ_transmul_init(gel(v_x,k1+1), S, T, p);
    1486           8 :     c = cgetg(m+2,t_POL);
    1487           8 :     c[1] = evalsigne(1)|evalvarn(vS);
    1488          40 :     for (i=0; i<m; i+=k1)
    1489             :     {
    1490          32 :       long mj = minss(m-i, k1);
    1491         100 :       for (j=0; j<mj; j++)
    1492          68 :         gel(c,m+1-(i+j)) = FpXQX_dotproduct(v, gel(v_x,j+1), T, p);
    1493          32 :       v = FpXQXQ_transmul(tr, v, n, T, p);
    1494             :     }
    1495           8 :     c = FpXX_renormalize(c, m+2);
    1496             :     /* now c contains <v,x^i> , i = 0..m-1  */
    1497           8 :     M = FpXQX_halfgcd(polxn_FpXX(m, vS, vT), c, T, p);
    1498           8 :     g_prime = gmael(M, 2, 2);
    1499           8 :     if (degpol(g_prime) < 1) continue;
    1500           8 :     g = FpXQX_mul(g, g_prime, T, p);
    1501           8 :     tau = FpXQXQ_mul(tau, FpXQX_FpXQXQV_eval(g_prime, v_x, S, T, p), S, T, p);
    1502             :   }
    1503           8 :   g = FpXQX_normalize(g,T, p);
    1504           8 :   return gerepilecopy(ltop,g);
    1505             : }
    1506             : 
    1507             : GEN
    1508           0 : FpXQXQ_matrix_pow(GEN y, long n, long m, GEN S, GEN T, GEN p)
    1509             : {
    1510           0 :   return RgXV_to_RgM(FpXQXQ_powers(y,m-1,S,T,p),n);
    1511             : }
    1512             : 
    1513             : GEN
    1514        2629 : FpXQX_FpXQXQV_eval(GEN P, GEN V, GEN S, GEN T, GEN p)
    1515             : {
    1516             :   struct _FpXQXQ D;
    1517        2629 :   T = FpX_get_red(T, p);
    1518        2629 :   S = FpXQX_get_red(S, T, p);
    1519        2629 :   D.S=S; D.T=T; D.p=p;
    1520        2629 :   return gen_bkeval_powers(P, degpol(P), V, (void*)&D, &FpXQXQ_algebra,
    1521             :                                                    _FpXQXQ_cmul);
    1522             : }
    1523             : 
    1524             : GEN
    1525         353 : FpXQX_FpXQXQ_eval(GEN Q, GEN x, GEN S, GEN T, GEN p)
    1526             : {
    1527             :   struct _FpXQXQ D;
    1528         353 :   int use_sqr = 2*degpol(x) >= get_FpXQX_degree(S);
    1529         353 :   T = FpX_get_red(T, p);
    1530         353 :   S = FpXQX_get_red(S, T, p);
    1531         353 :   D.S=S; D.T=T; D.p=p;
    1532         353 :   return gen_bkeval(Q, degpol(Q), x, use_sqr, (void*)&D, &FpXQXQ_algebra,
    1533             :       _FpXQXQ_cmul);
    1534             : }
    1535             : 
    1536             : static GEN
    1537         245 : FpXQXQ_autpow_sqr(void * E, GEN x)
    1538             : {
    1539         245 :   struct _FpXQXQ *D = (struct _FpXQXQ *)E;
    1540         245 :   GEN S = D->S, T = D->T, p = D->p;
    1541         245 :   GEN phi = gel(x,1), S1 = gel(x,2);
    1542         245 :   long n = brent_kung_optpow(get_FpX_degree(T)-1,lgpol(S1)+1,1);
    1543         245 :   GEN V = FpXQ_powers(phi, n, T, p);
    1544         245 :   GEN phi2 = FpX_FpXQV_eval(phi, V, T, p);
    1545         245 :   GEN Sphi = FpXY_FpXQV_evalx(S1, V, T, p);
    1546         245 :   GEN S2 = FpXQX_FpXQXQ_eval(Sphi, S1, S, T, p);
    1547         245 :   return mkvec2(phi2, S2);
    1548             : }
    1549             : 
    1550             : static GEN
    1551          85 : FpXQXQ_autpow_mul(void * E, GEN x, GEN y)
    1552             : {
    1553          85 :   struct _FpXQXQ *D = (struct _FpXQXQ *)E;
    1554          85 :   GEN S = D->S, T = D->T, p = D->p;
    1555          85 :   GEN phi1 = gel(x,1), S1 = gel(x,2);
    1556          85 :   GEN phi2 = gel(y,1), S2 = gel(y,2);
    1557          85 :   long n = brent_kung_optpow(get_FpX_degree(T)-1, lgpol(S1)+1, 1);
    1558          85 :   GEN V = FpXQ_powers(phi2, n, T, p);
    1559          85 :   GEN phi3 = FpX_FpXQV_eval(phi1, V, T, p);
    1560          85 :   GEN Sphi = FpXY_FpXQV_evalx(S1, V, T, p);
    1561          85 :   GEN S3 = FpXQX_FpXQXQ_eval(Sphi, S2, S, T, p);
    1562          85 :   return mkvec2(phi3, S3);
    1563             : }
    1564             : 
    1565             : GEN
    1566         231 : FpXQXQ_autpow(GEN aut, long n, GEN S, GEN T, GEN p)
    1567             : {
    1568             :   struct _FpXQXQ D;
    1569         231 :   T = FpX_get_red(T, p);
    1570         231 :   S = FpXQX_get_red(S, T, p);
    1571         231 :   D.S=S; D.T=T; D.p=p;
    1572         231 :   return gen_powu(aut,n,&D,FpXQXQ_autpow_sqr,FpXQXQ_autpow_mul);
    1573             : }
    1574             : 
    1575             : static GEN
    1576           1 : FpXQXQ_auttrace_mul(void *E, GEN x, GEN y)
    1577             : {
    1578           1 :   struct _FpXQXQ *D = (struct _FpXQXQ *)E;
    1579           1 :   GEN S = D->S, T = D->T;
    1580           1 :   GEN p = D->p;
    1581           1 :   GEN S1 = gel(x,1), a1 = gel(x,2);
    1582           1 :   GEN S2 = gel(y,1), a2 = gel(y,2);
    1583           1 :   long n = brent_kung_optpow(maxss(degpol(S1),degpol(a1)),2,1);
    1584           1 :   GEN V = FpXQXQ_powers(S2, n, S, T, p);
    1585           1 :   GEN S3 = FpXQX_FpXQXQV_eval(S1, V, S, T, p);
    1586           1 :   GEN aS = FpXQX_FpXQXQV_eval(a1, V, S, T, p);
    1587           1 :   GEN a3 = FpXX_add(aS, a2, p);
    1588           1 :   return mkvec2(S3, a3);
    1589             : }
    1590             : 
    1591             : static GEN
    1592           1 : FpXQXQ_auttrace_sqr(void *E, GEN x)
    1593           1 : { return FpXQXQ_auttrace_mul(E, x, x); }
    1594             : 
    1595             : GEN
    1596           8 : FpXQXQ_auttrace(GEN aut, long n, GEN S, GEN T, GEN p)
    1597             : {
    1598             :   struct _FpXQXQ D;
    1599           8 :   T = FpX_get_red(T, p);
    1600           8 :   S = FpXQX_get_red(S, T, p);
    1601           8 :   D.S=S; D.T=T; D.p=p;
    1602           8 :   return gen_powu(aut,n,&D,FpXQXQ_auttrace_sqr,FpXQXQ_auttrace_mul);
    1603             : }
    1604             : 
    1605             : static GEN
    1606        1174 : FpXQXQ_autsum_mul(void *E, GEN x, GEN y)
    1607             : {
    1608        1174 :   struct _FpXQXQ *D = (struct _FpXQXQ *) E;
    1609        1174 :   GEN S = D->S, T = D->T, p = D->p;
    1610        1174 :   GEN phi1 = gel(x,1), S1 = gel(x,2), a1 = gel(x,3);
    1611        1174 :   GEN phi2 = gel(y,1), S2 = gel(y,2), a2 = gel(y,3);
    1612        1174 :   long n2 = brent_kung_optpow(get_FpX_degree(T)-1, lgpol(S1)+lgpol(a1)+1, 1);
    1613        1174 :   GEN V2 = FpXQ_powers(phi2, n2, T, p);
    1614        1174 :   GEN phi3 = FpX_FpXQV_eval(phi1, V2, T, p);
    1615        1174 :   GEN Sphi = FpXY_FpXQV_evalx(S1, V2, T, p);
    1616        1174 :   GEN aphi = FpXY_FpXQV_evalx(a1, V2, T, p);
    1617        1174 :   long n = brent_kung_optpow(maxss(degpol(Sphi),degpol(aphi)),2,1);
    1618        1174 :   GEN V = FpXQXQ_powers(S2, n, S, T, p);
    1619        1174 :   GEN S3 = FpXQX_FpXQXQV_eval(Sphi, V, S, T, p);
    1620        1174 :   GEN aS = FpXQX_FpXQXQV_eval(aphi, V, S, T, p);
    1621        1174 :   GEN a3 = FpXQXQ_mul(aS, a2, S, T, p);
    1622        1174 :   return mkvec3(phi3, S3, a3);
    1623             : }
    1624             : 
    1625             : static GEN
    1626        1134 : FpXQXQ_autsum_sqr(void * T, GEN x)
    1627        1134 : { return FpXQXQ_autsum_mul(T,x,x); }
    1628             : 
    1629             : GEN
    1630        1120 : FpXQXQ_autsum(GEN aut, long n, GEN S, GEN T, GEN p)
    1631             : {
    1632             :   struct _FpXQXQ D;
    1633        1120 :   T = FpX_get_red(T, p);
    1634        1120 :   S = FpXQX_get_red(S, T, p);
    1635        1120 :   D.S=S; D.T=T; D.p=p;
    1636        1120 :   return gen_powu(aut,n,&D,FpXQXQ_autsum_sqr,FpXQXQ_autsum_mul);
    1637             : }
    1638             : 
    1639             : GEN
    1640          32 : FpXQXn_mul(GEN x, GEN y, long n, GEN T, GEN p)
    1641             : {
    1642          32 :   pari_sp av = avma;
    1643             :   GEN z, kx, ky;
    1644             :   long d;
    1645          32 :   if (ZXX_is_ZX(y) && ZXX_is_ZX(x))
    1646           0 :     return FpXn_mul(x,y,n,p);
    1647          32 :   d = get_FpX_degree(T);
    1648          32 :   kx = ZXX_to_Kronecker(x, d);
    1649          32 :   ky = ZXX_to_Kronecker(y, d);
    1650          32 :   z = Kronecker_to_FpXQX(ZXn_mul(ky,kx,(2*d-1)*n), T, p);
    1651          32 :   return gerepileupto(av, z);
    1652             : }
    1653             : 
    1654             : GEN
    1655           0 : FpXQXn_sqr(GEN x, long n, GEN T, GEN p)
    1656             : {
    1657           0 :   pari_sp av = avma;
    1658             :   GEN z, kx;
    1659             :   long d;
    1660           0 :   if (ZXX_is_ZX(x)) return ZXn_sqr(x, n);
    1661           0 :   d = get_FpX_degree(T);
    1662           0 :   kx= ZXX_to_Kronecker(x, d);
    1663           0 :   z = Kronecker_to_FpXQX(ZXn_sqr(kx, (2*d-1)*n), T, p);
    1664           0 :   return gerepileupto(av, z);
    1665             : }
    1666             : 
    1667             : INLINE GEN
    1668           0 : FpXXn_red(GEN a, long n)
    1669           0 : { return RgXn_red_shallow(a, n); }
    1670             : 
    1671             : GEN
    1672           0 : FpXQXn_exp(GEN h, long e, GEN T, GEN p)
    1673             : {
    1674           0 :   pari_sp av = avma, av2;
    1675           0 :   long v = varn(h), n=1;
    1676           0 :   GEN f = pol_1(v), g = pol_1(v);
    1677           0 :   ulong mask = quadratic_prec_mask(e);
    1678           0 :   av2 = avma;
    1679           0 :   if (signe(h)==0 || degpol(h)<1 || !gequal0(gel(h,2)))
    1680           0 :     pari_err_DOMAIN("FpXQXn_exp","valuation", "<", gen_1, h);
    1681           0 :   for (;mask>1;)
    1682             :   {
    1683             :     GEN q, w;
    1684           0 :     long n2 = n;
    1685           0 :     n<<=1; if (mask & 1) n--;
    1686           0 :     mask >>= 1;
    1687           0 :     g = FpXX_sub(FpXX_mulu(g,2,p), FpXQXn_mul(f, FpXQXn_sqr(g, n2, T, p), n2, T, p), p);
    1688           0 :     q = FpXX_deriv(FpXXn_red(h,n2), p);
    1689           0 :     w = FpXX_add(q, FpXQXn_mul(g, FpXX_sub(FpXX_deriv(f, p), FpXQXn_mul(f,q,n-1, T, p), p),n-1, T, p), p);
    1690           0 :     f = FpXX_add(f, FpXQXn_mul(f, FpXX_sub(FpXXn_red(h, n), FpXX_integ(w, p), p), n, T, p), p);
    1691           0 :     if (gc_needed(av2,2))
    1692             :     {
    1693           0 :       if (DEBUGMEM>1) pari_warn(warnmem,"FpXQXn_exp, e = %ld", n);
    1694           0 :       gerepileall(av2, 2, &f, &g);
    1695             :     }
    1696             :   }
    1697           0 :   return gerepileupto(av, f);
    1698             : }
    1699             : 
    1700             : /* (f*g) \/ x^n */
    1701             : static GEN
    1702           0 : FpXQX_mulhigh_i(GEN f, GEN g, long n, GEN T, GEN p)
    1703             : {
    1704           0 :   return RgX_shift_shallow(FpXQX_mul(f,g,T, p),-n);
    1705             : }
    1706             : 
    1707             : static GEN
    1708           0 : FpXQXn_mulhigh(GEN f, GEN g, long n2, long n, GEN T, GEN p)
    1709             : {
    1710           0 :   GEN F = RgX_blocks(f, n2, 2), fl = gel(F,1), fh = gel(F,2);
    1711           0 :   return FpXX_add(FpXQX_mulhigh_i(fl, g, n2, T, p), FpXQXn_mul(fh, g, n - n2, T, p), p);
    1712             : }
    1713             : 
    1714             : GEN
    1715           0 : FpXQXn_inv(GEN f, long e, GEN T, GEN p)
    1716             : {
    1717           0 :   pari_sp av = avma, av2;
    1718             :   ulong mask;
    1719             :   GEN W, a;
    1720           0 :   long v = varn(f), n = 1;
    1721             : 
    1722           0 :   if (!signe(f)) pari_err_INV("FpXXn_inv",f);
    1723           0 :   a = Fq_inv(gel(f,2), T, p);
    1724           0 :   if (e == 1) return scalarpol(a, v);
    1725           0 :   else if (e == 2)
    1726             :   {
    1727             :     GEN b;
    1728           0 :     if (degpol(f) <= 0) return scalarpol(a, v);
    1729           0 :     b = Fq_neg(gel(f,3),T,p);
    1730           0 :     if (signe(b)==0) return scalarpol(a, v);
    1731           0 :     if (!is_pm1(a)) b = Fq_mul(b, Fq_sqr(a, T, p), T, p);
    1732           0 :     W = deg1pol_shallow(b, a, v);
    1733           0 :     return gerepilecopy(av, W);
    1734             :   }
    1735           0 :   W = scalarpol_shallow(Fq_inv(gel(f,2), T, p),v);
    1736           0 :   mask = quadratic_prec_mask(e);
    1737           0 :   av2 = avma;
    1738           0 :   for (;mask>1;)
    1739             :   {
    1740             :     GEN u, fr;
    1741           0 :     long n2 = n;
    1742           0 :     n<<=1; if (mask & 1) n--;
    1743           0 :     mask >>= 1;
    1744           0 :     fr = FpXXn_red(f, n);
    1745           0 :     u = FpXQXn_mul(W, FpXQXn_mulhigh(fr, W, n2, n, T, p), n-n2, T, p);
    1746           0 :     W = FpXX_sub(W, RgX_shift_shallow(u, n2), p);
    1747           0 :     if (gc_needed(av2,2))
    1748             :     {
    1749           0 :       if(DEBUGMEM>1) pari_warn(warnmem,"FpXQXn_inv, e = %ld", n);
    1750           0 :       W = gerepileupto(av2, W);
    1751             :     }
    1752             :   }
    1753           0 :   return gerepileupto(av, W);
    1754             : }

Generated by: LCOV version 1.13