Code coverage tests

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

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

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

LCOV - code coverage report
Current view: top level - basemath - Hensel.c (source / functions) Hit Total Coverage
Test: PARI/GP v2.8.0 lcov report (development 19217-a6dcf64) Lines: 514 536 95.9 %
Date: 2016-07-27 07:10:32 Functions: 38 39 97.4 %
Legend: Lines: hit not hit

          Line data    Source code
       1             : /* Copyright (C) 2000  The PARI group.
       2             : 
       3             : This file is part of the PARI/GP package.
       4             : 
       5             : PARI/GP is free software; you can redistribute it and/or modify it under the
       6             : terms of the GNU General Public License as published by the Free Software
       7             : Foundation. It is distributed in the hope that it will be useful, but WITHOUT
       8             : ANY WARRANTY WHATSOEVER.
       9             : 
      10             : Check the License for details. You should have received a copy of it, along
      11             : with the package; see the file 'COPYING'. If not, write to the Free Software
      12             : Foundation, Inc., 51 Franklin Street, Fifth Floor, Boston, MA 02110-1301 USA. */
      13             : #include "pari.h"
      14             : #include "paripriv.h"
      15             : 
      16             : /***********************************************************************/
      17             : /**                                                                   **/
      18             : /**       QUADRATIC HENSEL LIFT (adapted from V. Shoup's NTL)         **/
      19             : /**                                                                   **/
      20             : /***********************************************************************/
      21             : /* Setup for divide/conquer quadratic Hensel lift
      22             :  *  a = set of k t_POL in Z[X] = factors over Fp (T=NULL) or Fp[Y]/(T)
      23             :  *  V = set of products of factors built as follows
      24             :  *  1) V[1..k] = initial a
      25             :  *  2) iterate:
      26             :  *      append to V the two smallest factors (minimal degree) in a, remove them
      27             :  *      from a and replace them by their product [net loss for a = 1 factor]
      28             :  *
      29             :  * W = bezout coeffs W[i]V[i] + W[i+1]V[i+1] = 1
      30             :  *
      31             :  * link[i] = -j if V[i] = a[j]
      32             :  *            j if V[i] = V[j] * V[j+1]
      33             :  * Arrays (link, V, W) pre-allocated for 2k - 2 elements */
      34             : static void
      35       10706 : BuildTree(GEN link, GEN V, GEN W, GEN a, GEN T, GEN p)
      36             : {
      37       10706 :   long k = lg(a)-1;
      38             :   long i, j, s, minp, mind;
      39             : 
      40       10706 :   for (i=1; i<=k; i++) { gel(V,i) = gel(a,i); link[i] = -i; }
      41             : 
      42       40357 :   for (j=1; j <= 2*k-5; j+=2,i++)
      43             :   {
      44       29651 :     minp = j;
      45       29651 :     mind = degpol(gel(V,j));
      46      414385 :     for (s=j+1; s<i; s++)
      47      384734 :       if (degpol(gel(V,s)) < mind) { minp = s; mind = degpol(gel(V,s)); }
      48             : 
      49       29651 :     swap(gel(V,j), gel(V,minp));
      50       29651 :     lswap(link[j], link[minp]);
      51             : 
      52       29651 :     minp = j+1;
      53       29651 :     mind = degpol(gel(V,j+1));
      54      384734 :     for (s=j+2; s<i; s++)
      55      355083 :       if (degpol(gel(V,s)) < mind) { minp = s; mind = degpol(gel(V,s)); }
      56             : 
      57       29651 :     swap(gel(V,j+1), gel(V,minp));
      58       29651 :     lswap(link[j+1], link[minp]);
      59             : 
      60       29651 :     gel(V,i) = FqX_mul(gel(V,j), gel(V,j+1), T, p);
      61       29651 :     link[i] = j;
      62             :   }
      63             : 
      64       51056 :   for (j=1; j <= 2*k-3; j+=2)
      65             :   {
      66             :     GEN d, u, v;
      67       40357 :     d = FqX_extgcd(gel(V,j), gel(V,j+1), T, p, &u, &v);
      68       40357 :     if (degpol(d) > 0) pari_err_COPRIME("BuildTree", gel(V,j), gel(V,j+1));
      69       40350 :     d = gel(d,2);
      70       40350 :     if (!gequal1(d))
      71             :     {
      72       35245 :       if (typ(d)==t_POL)
      73             :       {
      74        5740 :         d = FpXQ_inv(d, T, p);
      75        5740 :         u = FqX_Fq_mul(u, d, T, p);
      76        5740 :         v = FqX_Fq_mul(v, d, T, p);
      77             :       }
      78             :       else
      79             :       {
      80       29505 :         d = Fp_inv(d, p);
      81       29505 :         u = FqX_Fp_mul(u, d, T,p);
      82       29505 :         v = FqX_Fp_mul(v, d, T,p);
      83             :       }
      84             :     }
      85       40350 :     gel(W,j) = u;
      86       40350 :     gel(W,j+1) = v;
      87             :   }
      88       10699 : }
      89             : 
      90             : /* au + bv = 1 (p0), ab = f (p0). Lift mod p1 = p0 pd (<= p0^2).
      91             :  * If noinv is set, don't lift the inverses u and v */
      92             : static void
      93      191860 : ZpX_HenselLift(GEN V, GEN W, long j, GEN f, GEN pd, GEN p0, GEN p1, int noinv)
      94             : {
      95      191860 :   pari_sp av = avma;
      96      191860 :   long space = lg(f) * lgefint(p1);
      97             :   GEN a2, b2, g, z, s, t;
      98      191860 :   GEN a = gel(V,j), b = gel(V,j+1);
      99      191860 :   GEN u = gel(W,j), v = gel(W,j+1);
     100             : 
     101      191860 :   (void)new_chunk(space); /* HACK */
     102      191860 :   g = ZX_sub(f, ZX_mul(a,b));
     103      191860 :   g = ZX_Z_divexact(g, p0);
     104      191860 :   g = FpX_red(g, pd);
     105      191860 :   z = FpX_mul(v,g, pd);
     106      191860 :   t = FpX_divrem(z,a, pd, &s);
     107      191860 :   t = ZX_add(ZX_mul(u,g), ZX_mul(t,b));
     108      191860 :   t = FpX_red(t, pd);
     109      191860 :   t = ZX_Z_mul(t,p0);
     110      191860 :   s = ZX_Z_mul(s,p0);
     111      191860 :   avma = av;
     112      191860 :   a2 = ZX_add(a,s);
     113      191860 :   b2 = ZX_add(b,t);
     114             : 
     115             :   /* already reduced mod p1 = pd p0 */
     116      191860 :   gel(V,j)   = a2;
     117      191860 :   gel(V,j+1) = b2;
     118      383720 :   if (noinv) return;
     119             : 
     120      160953 :   av = avma;
     121      160953 :   (void)new_chunk(space); /* HACK */
     122      160953 :   g = ZX_add(ZX_mul(u,a2), ZX_mul(v,b2));
     123      160953 :   g = Z_ZX_sub(gen_1, g);
     124      160953 :   g = ZX_Z_divexact(g, p0);
     125      160953 :   g = FpX_red(g, pd);
     126      160953 :   z = FpX_mul(v,g, pd);
     127      160953 :   t = FpX_divrem(z,a, pd, &s);
     128      160953 :   t = ZX_add(ZX_mul(u,g), ZX_mul(t,b));
     129      160953 :   t = FpX_red(t, pd);
     130      160953 :   t = ZX_Z_mul(t,p0);
     131      160953 :   s = ZX_Z_mul(s,p0);
     132      160953 :   avma = av;
     133      160953 :   gel(W,j)   = ZX_add(u,t);
     134      160953 :   gel(W,j+1) = ZX_add(v,s);
     135             : }
     136             : 
     137             : static void
     138       34706 : ZpXQ_HenselLift(GEN V, GEN W, long j, GEN f, GEN Td, GEN T1, GEN pd, GEN p0, GEN p1, int noinv)
     139             : {
     140       34706 :   pari_sp av = avma;
     141       34706 :   const long n = degpol(T1), vT = varn(T1);
     142       34706 :   long space = lg(f) * lgefint(p1) * lg(T1);
     143             :   GEN a2, b2, g, z, s, t;
     144       34706 :   GEN a = gel(V,j), b = gel(V,j+1);
     145       34706 :   GEN u = gel(W,j), v = gel(W,j+1);
     146             : 
     147       34706 :   (void)new_chunk(space); /* HACK */
     148       34706 :   g = RgX_sub(f, Kronecker_to_ZXX(ZXX_mul_Kronecker(a,b,n), n, vT));
     149       34706 :   g = FpXQX_red(g, T1, p1);
     150       34706 :   g = RgX_Rg_divexact(g, p0);
     151       34706 :   z = FpXQX_mul(v,g, Td,pd);
     152       34706 :   t = FpXQX_divrem(z,a, Td,pd, &s);
     153       34706 :   t = ZX_add(ZXX_mul_Kronecker(u,g,n), ZXX_mul_Kronecker(t,b,n));
     154       34706 :   t = Kronecker_to_ZXX(t, n, vT);
     155       34706 :   t = FpXQX_red(t, Td, pd);
     156       34706 :   t = RgX_Rg_mul(t,p0);
     157       34706 :   s = RgX_Rg_mul(s,p0);
     158       34706 :   avma = av;
     159             : 
     160       34706 :   a2 = RgX_add(a,s);
     161       34706 :   b2 = RgX_add(b,t);
     162             :   /* already reduced mod p1 = pd p0 */
     163       34706 :   gel(V,j)   = a2;
     164       34706 :   gel(V,j+1) = b2;
     165       69412 :   if (noinv) return;
     166             : 
     167       27748 :   av = avma;
     168       27748 :   (void)new_chunk(space); /* HACK */
     169       27748 :   g = ZX_add(ZXX_mul_Kronecker(u,a2,n), ZXX_mul_Kronecker(v,b2,n));
     170       27748 :   g = Kronecker_to_ZXX(g, n, vT);
     171       27748 :   g = Rg_RgX_sub(gen_1, g);
     172       27748 :   g = FpXQX_red(g, T1, p1);
     173       27748 :   g = RgX_Rg_divexact(g, p0);
     174       27748 :   z = FpXQX_mul(v,g, Td,pd);
     175       27748 :   t = FpXQX_divrem(z,a, Td,pd, &s);
     176       27748 :   t = ZX_add(ZXX_mul_Kronecker(u,g,n), ZXX_mul_Kronecker(t,b,n));
     177       27748 :   t = Kronecker_to_ZXX(t, n, vT);
     178       27748 :   t = FpXQX_red(t, Td, pd);
     179       27748 :   t = RgX_Rg_mul(t,p0);
     180       27748 :   s = RgX_Rg_mul(s,p0);
     181       27748 :   avma = av;
     182       27748 :   gel(W,j)   = RgX_add(u,t);
     183       27748 :   gel(W,j+1) = RgX_add(v,s);
     184             : }
     185             : 
     186             : /* v list of factors, w list of inverses.  f = v[j] v[j+1]
     187             :  * Lift v[j] and v[j+1] mod p0 pd (possibly mod T), then all their divisors */
     188             : static void
     189      430927 : ZpX_RecTreeLift(GEN link, GEN v, GEN w, GEN pd, GEN p0, GEN p1,
     190             :                 GEN f, long j, int noinv)
     191             : {
     192      861854 :   if (j < 0) return;
     193      191860 :   ZpX_HenselLift(v, w, j, f, pd, p0,p1, noinv);
     194      191860 :   ZpX_RecTreeLift(link, v, w, pd, p0,p1, gel(v,j)  , link[j  ], noinv);
     195      191860 :   ZpX_RecTreeLift(link, v, w, pd, p0,p1, gel(v,j+1), link[j+1], noinv);
     196             : }
     197             : static void
     198       73367 : ZpXQ_RecTreeLift(GEN link, GEN v, GEN w, GEN Td, GEN T1, GEN pd, GEN p0, GEN p1,
     199             :                  GEN f, long j, int noinv)
     200             : {
     201      146734 :   if (j < 0) return;
     202       34706 :   ZpXQ_HenselLift(v, w, j, f, Td,T1, pd, p0,p1, noinv);
     203       34706 :   ZpXQ_RecTreeLift(link, v, w, Td,T1, pd, p0,p1, gel(v,j)  , link[j  ], noinv);
     204       34706 :   ZpXQ_RecTreeLift(link, v, w, Td,T1, pd, p0,p1, gel(v,j+1), link[j+1], noinv);
     205             : }
     206             : 
     207             : /* Assume n > 0. We want to go to accuracy n, starting from accuracy 1, using
     208             :  * a quadratically convergent algorithm. Goal: 9 -> 1,2,3,5,9 instead of
     209             :  * 1,2,4,8,9 (sequence of accuracies).
     210             :  *
     211             :  * Let a0 = 1, a1 = 2, a2, ... ak = n, the sequence of accuracies. To obtain
     212             :  * it, work backwards:
     213             :  *   a(k) = n, a(i-1) = (a(i) + 1) \ 2,
     214             :  * but we do not want to store a(i) explicitly, even as a t_VECSMALL, since
     215             :  * this would leave an object on the stack. We store a(i) implicitly in a
     216             :  * MASK: let a(0) = 1, if the i-bit of MASK is set, set a(i+1) = 2 a(i) - 1,
     217             :  * and 2a(i) otherwise.
     218             :  *
     219             :  * In fact, we do something a little more complicated to simplify the
     220             :  * function interface and avoid returning k and MASK separately: we return
     221             :  * MASK + 2^(k+1), so the highest bit of the mask indicates the length of the
     222             :  * sequence, and the following ones are as above. */
     223             : ulong
     224      301388 : quadratic_prec_mask(long n)
     225             : {
     226      301388 :   long a = n, i;
     227      301388 :   ulong mask = 0;
     228      870703 :   for(i = 1;; i++, mask <<= 1)
     229             :   {
     230      870703 :     mask |= (a&1); a = (a+1)>>1;
     231     1172091 :     if (a==1) return mask | (1UL << i);
     232      569315 :   }
     233             : }
     234             : 
     235             : /* Lift to precision p^e0.
     236             :  * a = modular factors of f mod (p,T) [possibly T=NULL]
     237             :  *  OR a TreeLift structure [e, link, v, w]: go on lifting
     238             :  * flag = 0: standard.
     239             :  * flag = 1: return TreeLift structure */
     240             : static GEN
     241       10713 : MultiLift(GEN f, GEN a, GEN T, GEN p, long e0, long flag)
     242             : {
     243       10713 :   long i, eold, e, k = lg(a) - 1;
     244             :   GEN E, v, w, link, penew, Tnew;
     245             :   ulong mask;
     246             :   pari_timer Ti;
     247             : 
     248       10713 :   if (k < 2) pari_err_DOMAIN("MultiLift", "#(modular factors)", "<", gen_2,a);
     249       10713 :   if (e0 < 1) pari_err_DOMAIN("MultiLift", "precision", "<", gen_1,stoi(e0));
     250       10713 :   if (e0 == 1) return a;
     251             : 
     252       10706 :   if (DEBUGLEVEL > 3) timer_start(&Ti);
     253       10706 :   if (typ(gel(a,1)) == t_INT)
     254             :   { /* a = TreeLift structure */
     255           0 :     e = itos(gel(a,1));
     256           0 :     link = gel(a,2);
     257           0 :     v    = gel(a,3);
     258           0 :     w    = gel(a,4);
     259             :   }
     260             :   else
     261             :   {
     262       10706 :     e = 1;
     263       10706 :     v = cgetg(2*k-2 + 1, t_VEC);
     264       10706 :     w = cgetg(2*k-2 + 1, t_VEC);
     265       10706 :     link=cgetg(2*k-2 + 1, t_VECSMALL);
     266       10706 :     BuildTree(link, v, w, a, T? FpX_red(T,p): NULL, p);
     267       10699 :     if (DEBUGLEVEL > 3) timer_printf(&Ti, "building tree");
     268             :   }
     269       10699 :   mask = quadratic_prec_mask(e0);
     270       10699 :   eold = 1;
     271       10699 :   penew = NULL;
     272       10699 :   Tnew = NULL;
     273       72560 :   while (mask > 1)
     274             :   {
     275       51162 :     long enew = eold << 1;
     276       51162 :     if (mask & 1) enew--;
     277       51162 :     mask >>= 1;
     278       51162 :     if (enew >= e) { /* mask == 1: last iteration */
     279       51162 :       GEN peold = penew? penew: powiu(p, eold);
     280       51162 :       GEN Td = NULL, pd;
     281       51162 :       long d = enew - eold; /* = eold or eold-1 */
     282             :       /* lift from p^eold to p^enew */
     283       51162 :       pd = (d == eold)? peold: diviiexact(peold, p); /* p^d */
     284       51162 :       penew = mulii(peold,pd);
     285       51162 :       if (T) {
     286        3955 :         if (Tnew)
     287        2961 :           Td = (d == eold)? Tnew: FpX_red(Tnew,pd);
     288             :         else
     289         994 :           Td = FpX_red(T, peold);
     290        3955 :         Tnew = FpX_red(T, penew);
     291        3955 :         ZpXQ_RecTreeLift(link, v, w, Td, Tnew, pd, peold, penew, f, lg(v)-2,
     292             :                          (flag == 0 && mask == 1));
     293             :       }
     294             :       else
     295       47207 :         ZpX_RecTreeLift(link, v, w, pd, peold, penew, f, lg(v)-2,
     296             :                         (flag == 0 && mask == 1));
     297       51162 :       if (DEBUGLEVEL > 3) timer_printf(&Ti, "lifting to prec %ld", enew);
     298             :     }
     299       51162 :     eold = enew;
     300             :   }
     301             : 
     302       10699 :   if (flag)
     303         588 :     E = mkvec4(utoipos(e0), link, v, w);
     304             :   else
     305             :   {
     306       10111 :     E = cgetg(k+1, t_VEC);
     307       85841 :     for (i = 1; i <= 2*k-2; i++)
     308             :     {
     309       75730 :       long t = link[i];
     310       75730 :       if (t < 0) gel(E,-t) = gel(v,i);
     311             :     }
     312             :   }
     313       10699 :   return E;
     314             : }
     315             : 
     316             : /* Q list of (coprime, monic) factors of pol mod (T,p). Lift mod p^e = pe.
     317             :  * T may be NULL */
     318             : GEN
     319       10447 : ZpX_liftfact(GEN pol, GEN Q, GEN T, GEN p, long e, GEN pe)
     320             : {
     321       10447 :   pari_sp av = avma;
     322       10447 :   if (lg(Q) == 2) return mkvec(pol);
     323       10125 :   pol = FqX_normalize(pol, T, pe);
     324       10125 :   return gerepilecopy(av, MultiLift(pol, Q, T, p, e, 0));
     325             : }
     326             : 
     327             : /* U = NULL treated as 1 */
     328             : static void
     329        5558 : BezoutPropagate(GEN link, GEN v, GEN w, GEN pe, GEN U, GEN f, long j)
     330             : {
     331             :   GEN Q, R;
     332       11116 :   if (j < 0) return;
     333             : 
     334        2485 :   Q = FpX_mul(gel(v,j), gel(w,j), pe);
     335        2485 :   if (U)
     336             :   {
     337        1897 :     Q = FpXQ_mul(Q, U, f, pe);
     338        1897 :     R = FpX_sub(U, Q, pe);
     339             :   }
     340             :   else
     341         588 :     R = Fp_FpX_sub(gen_1, Q, pe);
     342        2485 :   gel(w,j+1) = Q; /*  0 mod U v[j],  1 mod (1-U) v[j+1] */
     343        2485 :   gel(w,j) = R; /*  1 mod U v[j],  0 mod (1-U) v[j+1] */
     344        2485 :   BezoutPropagate(link, v, w, pe, R, f, link[j  ]);
     345        2485 :   BezoutPropagate(link, v, w, pe, Q, f, link[j+1]);
     346             : }
     347             : 
     348             : /* as above, but return the Bezout coefficients for the lifted modular factors
     349             :  *   U[i] = 1 mod Qlift[i]
     350             :  *          0 mod Qlift[j], j != i */
     351             : GEN
     352         595 : bezout_lift_fact(GEN pol, GEN Q, GEN p, long e)
     353             : {
     354         595 :   pari_sp av = avma;
     355             :   GEN E, link, v, w, pe;
     356         595 :   long i, k = lg(Q)-1;
     357         595 :   if (k == 1) return mkvec(pol);
     358         588 :   pe = powiu(p, e);
     359         588 :   pol = FpX_normalize(pol, pe);
     360         588 :   E = MultiLift(pol, Q, NULL, p, e, 1);
     361         588 :   link = gel(E,2);
     362         588 :   v    = gel(E,3);
     363         588 :   w    = gel(E,4);
     364         588 :   BezoutPropagate(link, v, w, pe, NULL, pol, lg(v)-2);
     365         588 :   E = cgetg(k+1, t_VEC);
     366        5558 :   for (i = 1; i <= 2*k-2; i++)
     367             :   {
     368        4970 :     long t = link[i];
     369        4970 :     if (t < 0) E[-t] = w[i];
     370             :   }
     371         588 :   return gerepilecopy(av, E);
     372             : }
     373             : 
     374             : /* Front-end for ZpX_liftfact:
     375             :    lift the factorization of pol mod p given by L to p^N (if possible) */
     376             : GEN
     377          21 : polhensellift(GEN pol, GEN L, GEN p, long N)
     378             : {
     379          21 :   GEN T = NULL;
     380             :   long i, l, t;
     381          21 :   pari_sp av = avma;
     382             : 
     383          21 :   if (typ(pol) != t_POL) pari_err_TYPE("polhensellift",pol);
     384          21 :   RgX_check_ZXX(pol, "polhensellift");
     385          21 :   if (!is_vec_t(typ(L)) || lg(L) < 3) pari_err_TYPE("polhensellift",L);
     386          21 :   t = typ(p);
     387          21 :   if (t == t_VEC) /* [p, T] */
     388             :   {
     389           7 :     T = gel(p,2);
     390           7 :     if (typ(T) != t_POL) pari_err_TYPE("polhensellift",pol);
     391           7 :     RgX_check_ZX(T, "polhensellift");
     392           7 :     p = gel(p,1); t = typ(p);
     393             :   }
     394          21 :   if (t != t_INT) pari_err_TYPE("polhensellift",p);
     395          21 :   if (N < 1) pari_err_DOMAIN("polhensellift", "precision", "<", gen_1,stoi(N));
     396             : 
     397          21 :   l = lg(L); L = leafcopy(L);
     398          63 :   for (i = 1; i < l; i++)
     399             :   {
     400          42 :     if (typ(gel(L,i)) != t_POL)
     401           0 :       gel(L,i) = scalar_ZX_shallow(gel(L,i), varn(pol));
     402          42 :     RgX_check_ZXX(gel(L,i), "polhensellift");
     403             :   }
     404          21 :   return gerepilecopy(av, ZpX_liftfact(pol, L, T, p, N, powiu(p,N)));
     405             : }
     406             : 
     407             : static GEN
     408          59 : ZpX_liftroots_full(GEN f, GEN S, GEN p, long e)
     409             : {
     410          59 :   long i, n = lg(S)-1;
     411          59 :   GEN r = cgetg(n+1, typ(S));
     412          59 :   pari_sp av = avma;
     413          59 :   long v = varn(f);
     414             :   GEN y, q;
     415         953 :   for (i=1; i<=n; i++)
     416         894 :     gel(r,i) = deg1pol(gen_1, Fp_neg(gel(S, i), p), v);
     417          59 :   q = powiu(p, e);
     418          59 :   y = ZpX_liftfact(f, r, NULL, p, e, q);
     419          59 :   r = cgetg(n+1 ,t_COL);
     420         953 :   for (i=1; i<=n; i++)
     421         894 :     gel(r,i) = Fp_neg(gmael(y, i, 2), q);
     422          59 :   return gerepileupto(av, r);
     423             : }
     424             : 
     425             : GEN
     426       41930 : ZpX_roots(GEN F, GEN p, long e)
     427             : {
     428       41930 :   pari_sp av = avma;
     429       41930 :   long i, v = varn(F);
     430             :   GEN y, r, q;
     431       41930 :   GEN f = FpX_normalize(F, p);
     432       41930 :   GEN g = FpX_normalize(FpX_split_part(f, p), p);
     433       41930 :   GEN S = FpX_roots(g, p);
     434       41930 :   long l = lg(S)-1, n = l < degpol(f)? l+1: l;
     435       41930 :   if (l==0) return cgetg(1, t_COL);
     436       41930 :   if (l==1) return mkcol(ZpX_liftroot(F,gel(S,1),p,e));
     437        3290 :   r = cgetg(n+1 ,t_COL);
     438       23121 :   for (i=1; i<=l; i++)
     439       19831 :     gel(r,i) = deg1pol_shallow(gen_1, Fp_neg(gel(S,i), p), v);
     440        3290 :   if (l < degpol(f)) gel(r, n) = FpX_div(f, g, p);
     441        3290 :   q = powiu(p, e);
     442        3290 :   y = ZpX_liftfact(F, r, NULL, p, e, q);
     443        3290 :   r = cgetg(l+1 ,t_COL);
     444        3290 :   for (i=1; i<=l; i++) gel(r,i) = Fp_neg(gmael(y,i,2), q);
     445        3290 :   return gerepileupto(av, r);
     446             : }
     447             : 
     448             : /* SPEC:
     449             :  * p is a t_INT > 1, e >= 1
     450             :  * f is a ZX with leading term prime to p.
     451             :  * a is a simple root mod l for all l|p.
     452             :  * Return roots of f mod p^e, as integers (implicitly mod p^e)
     453             :  * STANDARD USE: p is a prime power */
     454             : GEN
     455       38948 : ZpX_liftroot(GEN f, GEN a, GEN p, long e)
     456             : {
     457       38948 :   pari_sp av = avma;
     458       38948 :   GEN q = p, fr, W;
     459             :   ulong mask;
     460             : 
     461       38948 :   a = modii(a,q);
     462       38948 :   if (e == 1) return a;
     463       38934 :   mask = quadratic_prec_mask(e);
     464       38934 :   fr = FpX_red(f,q);
     465       38934 :   W = Fp_inv(FpX_eval(ZX_deriv(fr), a, q), q); /* 1/f'(a) mod p */
     466             :   for(;;)
     467             :   {
     468      148581 :     q = sqri(q);
     469      148581 :     if (mask & 1) q = diviiexact(q, p);
     470      148581 :     mask >>= 1;
     471      148581 :     fr = FpX_red(f,q);
     472      148581 :     a = Fp_sub(a, Fp_mul(W, FpX_eval(fr, a,q), q), q);
     473      148581 :     if (mask == 1) return gerepileuptoint(av, a);
     474      109647 :     W = Fp_sub(shifti(W,1), Fp_mul(Fp_sqr(W,q), FpX_eval(ZX_deriv(fr),a,q), q), q);
     475      109647 :   }
     476             : }
     477             : 
     478             : GEN
     479          59 : ZpX_liftroots(GEN f, GEN S, GEN p, long e)
     480             : {
     481          59 :   long i, n = lg(S)-1, d = degpol(f);
     482             :   GEN r;
     483          59 :   if (n == d) return ZpX_liftroots_full(f, S, p, e);
     484           0 :   r = cgetg(n+1, typ(S));
     485           0 :   for (i=1; i <= n; i++)
     486           0 :     gel(r,i) = ZpX_liftroot(f, gel(S,i), p, e);
     487           0 :   return r;
     488             : }
     489             : 
     490             : GEN
     491          91 : ZpXQX_liftroot_vald(GEN f, GEN a, long v, GEN T, GEN p, long e)
     492             : {
     493          91 :   pari_sp av = avma, av2;
     494          91 :   GEN pv = p, q, W, df, Tq, fr, dfr;
     495             :   ulong mask;
     496          91 :   a = Fq_red(a, T, p);
     497          91 :   if (e <= v+1) return a;
     498          91 :   df = RgX_deriv(f);
     499          91 :   if (v) { pv = powiu(p,v); df = ZXX_Z_divexact(df, pv); }
     500          91 :   mask = quadratic_prec_mask(e-v);
     501          91 :   Tq = FpXT_red(T, p); dfr = FpXQX_red(df, Tq, p);
     502          91 :   W = Fq_inv(FqX_eval(dfr, a, Tq, p), Tq, p); /* 1/f'(a) mod (T,p) */
     503          91 :   q = p;
     504          91 :   av2 = avma;
     505             :   for (;;)
     506             :   {
     507             :     GEN u, fa, qv, q2v, q2, Tq2;
     508         161 :     q2 = q; q = sqri(q);
     509         161 :     if (mask & 1) q = diviiexact(q,p);
     510         161 :     mask >>= 1;
     511         161 :     if (v) { qv = mulii(q, pv); q2v = mulii(q2, pv); }
     512         161 :     else { qv = q; q2v = q2; }
     513         161 :     Tq2 = FpXT_red(T, q2v); Tq = FpXT_red(T, qv);
     514         161 :     fr = FpXQX_red(f, Tq, qv);
     515         161 :     fa = FqX_eval(fr, a, Tq, qv);
     516         161 :     fa = typ(fa)==t_INT? diviiexact(fa,q2v): ZX_Z_divexact(fa, q2v);
     517         161 :     a = Fq_sub(a, ZX_Z_mul(Fq_mul(W, fa, Tq2, q2v), q2), Tq, qv);
     518         161 :     if (mask == 1) return gerepileupto(av, a);
     519          70 :     dfr = FpXQX_red(df, Tq, q);
     520          70 :     u = ZX_Z_divexact(FpX_Fp_sub(Fq_mul(W,FqX_eval(dfr,a,Tq,q),Tq,q),gen_1,q),q2);
     521          70 :     W = Fq_sub(W,ZX_Z_mul(Fq_mul(u,W,Tq2,q2),q2),Tq,q);
     522          70 :     if (gc_needed(av2,2))
     523             :     {
     524           0 :       if(DEBUGMEM>1) pari_warn(warnmem,"ZpXQX_liftroot, e = %ld", e);
     525           0 :       gerepileall(av2, 3, &a, &W, &q);
     526             :     }
     527          70 :   }
     528             : }
     529             : 
     530             : GEN
     531          91 : ZpXQX_liftroot(GEN f, GEN a, GEN T, GEN p, long e) { return ZpXQX_liftroot_vald(f,a,0,T,p,e); }
     532             : 
     533             : /* Same as ZpX_liftroot for the polynomial X^n-b*/
     534             : GEN
     535       46879 : Zp_sqrtnlift(GEN b, GEN n, GEN a, GEN p, long e)
     536             : {
     537       46879 :   pari_sp ltop=avma;
     538             :   GEN q, w, n_1;
     539             :   ulong mask;
     540       46879 :   long pis2 = equalii(n, gen_2)? 1: 0;
     541       46879 :   if (e == 1) return icopy(a);
     542       44380 :   n_1 = subis(n,1);
     543       44380 :   mask = quadratic_prec_mask(e);
     544       44380 :   w = Fp_inv(pis2 ? shifti(a,1): Fp_mul(n,Fp_pow(a,n_1,p), p), p);
     545       44380 :   q = p;
     546             :   for(;;)
     547             :   {
     548      104776 :     q = sqri(q);
     549      104776 :     if (mask & 1) q = diviiexact(q, p);
     550      104776 :     mask >>= 1;
     551      104776 :     if (lgefint(q) == 3 && lgefint(n) == 3)
     552       60259 :     {
     553      104377 :       ulong Q = uel(q,2), N = uel(n,2);
     554      104377 :       ulong A = umodiu(a, Q);
     555      104377 :       ulong B = umodiu(b, Q);
     556      104377 :       ulong W = umodiu(w, Q);
     557             : 
     558      104377 :       A = Fl_sub(A, Fl_mul(W, Fl_sub(Fl_powu(A,N,Q), B, Q), Q), Q);
     559      104377 :       a = utoi(A);
     560      104377 :       if (mask == 1) break;
     561       60259 :       W = Fl_sub(Fl_add(W,W,Q),
     562             :                  Fl_mul(Fl_sqr(W,Q), Fl_mul(N,Fl_powu(A, N-1, Q), Q), Q), Q);
     563       60259 :       w = utoi(W);
     564             :     }
     565             :     else
     566             :     {
     567             :       /* a -= w (a^n - b) */
     568         399 :       a = modii(subii(a, mulii(w, subii(Fp_pow(a,n,q),b))), q);
     569         399 :       if (mask == 1) break;
     570             :       /* w += w - w^2 n a^(n-1)*/
     571         259 :       w = subii(shifti(w,1), Fp_mul(Fp_sqr(w,q),
     572         122 :                            pis2? shifti(a,1): mulii(n,Fp_pow(a,n_1,q)), q));
     573             :     }
     574       60396 :   }
     575       44380 :   return gerepileuptoint(ltop,a);
     576             : }
     577             : 
     578             : 
     579             : /* Same as ZpX_liftroot for the polynomial X^2-b */
     580             : GEN
     581       12831 : Zp_sqrtlift(GEN b, GEN a, GEN p, long e)
     582             : {
     583       12831 :   return Zp_sqrtnlift(b, gen_2, a, p, e);
     584             : }
     585             : 
     586             : GEN
     587        2408 : Zp_sqrt(GEN x, GEN p, long e)
     588             : {
     589             :   pari_sp av;
     590             :   GEN z;
     591        2408 :   if (equaliu(p,2)) return Z2_sqrt(x,e);
     592        1603 :   av = avma;
     593        1603 :   z = Fp_sqrt(Fp_red(x, p), p);
     594        1603 :   if (!z) return NULL;
     595        1575 :   if (e > 1) z = Zp_sqrtlift(x, z, p, e);
     596        1575 :   return gerepileuptoint(av, z);
     597             : }
     598             : 
     599             : /* Compute (x-1)/(x+1)/p^k */
     600             : static GEN
     601       20750 : ZpXQ_log_to_ath(GEN x, long k, GEN T, GEN p, long e, GEN pe)
     602             : {
     603       20750 :   pari_sp av = avma;
     604       20750 :   long vT = get_FpX_var(T);
     605             :   GEN bn, bdi;
     606       20750 :   GEN bd = ZX_Z_add(x, gen_1);
     607       20750 :   if (equaliu(p,2)) /*For p=2, we need to simplify by 2*/
     608             :   {
     609        6848 :     bn = ZX_shifti(x,-(k+1));
     610        6848 :     bdi= ZpXQ_invlift(ZX_shifti(bd ,-1), pol_1(vT), T, p, e);
     611             :   }
     612             :   else
     613             :   {
     614       13902 :     bn = ZX_Z_divexact(ZX_Z_sub(x, gen_1),powiu(p,k));
     615       13902 :     bdi= ZpXQ_invlift(bd, scalarpol(Fp_inv(gen_2,p),vT), T, p, e);
     616             :   }
     617       20750 :   return gerepileupto(av, FpXQ_mul(bn, bdi, T, pe));
     618             : }
     619             : 
     620             : /* Assume p odd, a = 1 [p], return log(a) */
     621             : GEN
     622       20750 : ZpXQ_log(GEN a, GEN T, GEN p, long N)
     623             : {
     624       20750 :   pari_sp av = avma;
     625             :   pari_timer ti;
     626       20750 :   long is2 = equaliu(p,2);
     627       20750 :   ulong pp = is2 ? 0: itou_or_0(p);
     628       20750 :   double lp = is2 ? 1: pp ? log2(pp): expi(p);
     629       20750 :   long k = maxss(1 , (long) .5+pow((double)(N>>1)/(lp*lp), 1./3));
     630             :   GEN ak, s, b, pol;
     631       20750 :   long e = is2 ? N-1: N;
     632       20750 :   long i, l = (e-2)/(2*(k+is2));
     633       20750 :   GEN pe = powiu(p,e);
     634       20750 :   GEN TNk, pNk = powiu(p,N+k);
     635       20750 :   if( DEBUGLEVEL>=3) timer_start(&ti);
     636       20750 :   TNk = FpX_get_red(get_FpX_mod(T), pNk);
     637       20750 :   ak = FpXQ_pow(a, powiu(p,k), TNk, pNk);
     638       20750 :   if( DEBUGLEVEL>=3) timer_printf(&ti,"FpXQ_pow(%ld)",k);
     639       20750 :   b = ZpXQ_log_to_ath(ak, k, T, p, e, pe);
     640       20750 :   if( DEBUGLEVEL>=3) timer_printf(&ti,"ZpXQ_log_to_ath");
     641       20750 :   pol= cgetg(l+3,t_POL);
     642       20750 :   pol[1] = evalsigne(1)|evalvarn(0);
     643       55317 :   for(i=0; i<=l; i++)
     644             :   {
     645             :     GEN g;
     646       34567 :     ulong z = 2*i+1;
     647       34567 :     if (pp)
     648             :     {
     649       26110 :       long w = u_lvalrem(z, pp, &z);
     650       26110 :       g = powuu(pp,2*i*k-w);
     651             :     }
     652        8457 :     else g = powiu(p,2*i*k);
     653       34568 :     gel(pol,i+2) = Fp_div(g, utoi(z),pe);
     654             :   }
     655       20750 :   if( DEBUGLEVEL>=3) timer_printf(&ti,"pol(%ld)",l);
     656       20750 :   s = FpX_FpXQ_eval(pol, FpXQ_sqr(b, T, pe), T,  pe);
     657       20750 :   if( DEBUGLEVEL>=3) timer_printf(&ti,"FpX_FpXQ_eval");
     658       20750 :   s = ZX_shifti(FpXQ_mul(b, s, T, pe), 1);
     659       20750 :   return gerepileupto(av, is2? s: FpX_red(s, pe));
     660             : }
     661             : 
     662             : /***********************************************************************/
     663             : /**                                                                   **/
     664             : /**                 Generic quadratic hensel lift over Zp[X]          **/
     665             : /**                                                                   **/
     666             : /***********************************************************************/
     667             : /* q = p^N */
     668             : 
     669             : GEN
     670      198212 : gen_ZpX_Dixon(GEN F, GEN V, GEN q, GEN p, long N, void *E,
     671             :                             GEN lin(void *E, GEN F, GEN d, GEN q),
     672             :                             GEN invl(void *E, GEN d))
     673             : {
     674      198212 :   pari_sp av = avma;
     675             :   long N2, M;
     676             :   GEN VN2, V2, VM, bil;
     677             :   GEN q2, qM;
     678      198212 :   V = FpX_red(V, q);
     679      198212 :   if (N == 1) return invl(E, V);
     680       44842 :   N2 = (N + 1)>>1; M = N - N2;
     681       44842 :   F = FpXT_red(F, q);
     682       44842 :   qM = powiu(p, M);
     683       44842 :   q2 = M == N2? qM: mulii(qM, p);
     684             :   /* q2 = p^N2, qM = p^M, q = q2 * qM */
     685       44842 :   VN2 = gen_ZpX_Dixon(F, V, q2, p, N2, E, lin, invl);
     686       44842 :   bil = lin(E, F, VN2, q);
     687       44842 :   V2 = ZX_Z_divexact(ZX_sub(V, bil), q2);
     688       44842 :   VM = gen_ZpX_Dixon(F, V2, qM, p, M, E, lin, invl);
     689       44842 :   return gerepileupto(av, FpX_red(ZX_add(VN2, ZX_Z_mul(VM, q2)), q));
     690             : }
     691             : 
     692             : GEN
     693      180417 : gen_ZpX_Newton(GEN x, GEN p, long n, void *E,
     694             :                       GEN eval(void *E, GEN f, GEN q),
     695             :                       GEN invd(void *E, GEN V, GEN v, GEN q, long M))
     696             : {
     697      180417 :   pari_sp ltop = avma, av;
     698      180417 :   long N = 1, N2, M;
     699             :   long mask;
     700      180417 :   GEN q = p;
     701      180417 :   if (n == 1) return gcopy(x);
     702      178443 :   mask = quadratic_prec_mask(n);
     703      178443 :   av = avma;
     704      788109 :   while (mask > 1)
     705             :   {
     706             :     GEN qM, q2, v, V;
     707      431223 :     N2 = N; N <<= 1;
     708      431223 :     q2 = q;
     709      431223 :     if (mask&1UL) { /* can never happen when q2 = p */
     710      162601 :       N--; M = N2-1;
     711      162601 :       qM = diviiexact(q2,p); /* > 1 */
     712      162601 :       q = mulii(qM,q2);
     713             :     } else {
     714      268622 :       M = N2;
     715      268622 :       qM = q2;
     716      268622 :       q = sqri(q2);
     717             :     }
     718             :     /* q2 = p^N2, qM = p^M, q = p^N = q2 * qM */
     719      431222 :     mask >>= 1;
     720      431222 :     v = eval(E, x, q);
     721      431223 :     V = ZX_Z_divexact(gel(v,1), q2);
     722      431224 :     x = FpX_sub(x, ZX_Z_mul(invd(E, V, v, qM, M), q2), q);
     723      431223 :     if (gc_needed(av, 1))
     724             :     {
     725           0 :       if(DEBUGMEM>1) pari_warn(warnmem,"gen_ZpX_Newton");
     726           0 :       gerepileall(av, 2, &x, &q);
     727             :     }
     728             :   }
     729      178443 :   return gerepileupto(ltop, x);
     730             : }
     731             : 
     732             : struct _ZpXQ_inv
     733             : {
     734             :   GEN T, a, p ,n;
     735             : };
     736             : 
     737             : static GEN
     738      291207 : _inv_invd(void *E, GEN V, GEN v, GEN q, long M/*unused*/)
     739             : {
     740      291207 :   struct _ZpXQ_inv *d = (struct _ZpXQ_inv *) E;
     741      291207 :   GEN Tq = FpXT_red(d->T, q);
     742             :   (void)M;
     743      291207 :   return FpXQ_mul(V, gel(v,2), Tq, q);
     744             : }
     745             : 
     746             : static GEN
     747      291206 : _inv_eval(void *E, GEN x, GEN q)
     748             : {
     749      291206 :   struct _ZpXQ_inv *d = (struct _ZpXQ_inv *) E;
     750      291206 :   GEN Tq = FpXT_red(d->T, q);
     751      291207 :   GEN f = FpX_Fp_sub(FpXQ_mul(x, FpX_red(d->a, q), Tq, q), gen_1, q);
     752      291207 :   return mkvec2(f, x);
     753             : }
     754             : 
     755             : GEN
     756      116153 : ZpXQ_invlift(GEN a, GEN x, GEN T, GEN p, long e)
     757             : {
     758             :   struct _ZpXQ_inv d;
     759      116153 :   d.a = a; d.T = T; d.p = p;
     760      116153 :   return gen_ZpX_Newton(x, p, e, &d, _inv_eval, _inv_invd);
     761             : }
     762             : 
     763             : GEN
     764       95403 : ZpXQ_inv(GEN a, GEN T, GEN p, long e)
     765             : {
     766       95403 :   pari_sp av=avma;
     767             :   GEN ai;
     768       95403 :   if (lgefint(p)==3)
     769             :   {
     770       95375 :     ulong pp = p[2];
     771       95375 :     ai = Flx_to_ZX(Flxq_inv(ZX_to_Flx(a,pp), ZXT_to_FlxT(T, pp), pp));
     772             :   } else
     773          28 :     ai = FpXQ_inv(FpX_red(a,p), FpXT_red(T,p),p);
     774       95403 :   return gerepileupto(av, ZpXQ_invlift(a, ai, T, p, e));
     775             : }
     776             : 
     777             : GEN
     778       31549 : ZpXQ_div(GEN a, GEN b, GEN T, GEN q, GEN p, long e)
     779             : {
     780       31549 :   return FpXQ_mul(a, ZpXQ_inv(b, T, p, e), T, q);
     781             : }
     782             : 
     783             : GEN
     784      142422 : ZpXQX_divrem(GEN x, GEN Sp, GEN T, GEN q, GEN p, long e, GEN *pr)
     785             : {
     786      142422 :   pari_sp av = avma;
     787      142422 :   GEN S = get_FpXQX_mod(Sp);
     788      142422 :   GEN b = leading_coeff(S), bi;
     789             :   GEN S2, Q;
     790      142422 :   if (typ(b)==t_INT) return FpXQX_divrem(x, Sp, T, q, pr);
     791       61796 :   bi = ZpXQ_inv(b, T, p, e);
     792       61796 :   S2 = FqX_Fq_mul_to_monic(S, bi, T, q);
     793       61796 :   Q = FpXQX_divrem(x, S2, T, q, pr);
     794       61796 :   if (pr==ONLY_DIVIDES && !Q) { avma = av; return NULL; }
     795       61796 :   if (pr==ONLY_REM || pr==ONLY_DIVIDES) return gerepileupto(av, Q);
     796       61796 :   Q = FpXQX_FpXQ_mul(Q, bi, T, q);
     797       61796 :   gerepileall(av, 2, &Q, pr);
     798       61796 :   return Q;
     799             : }
     800             : 
     801             : GEN
     802         588 : ZpXQX_digits(GEN x, GEN B, GEN T, GEN q, GEN p, long e)
     803             : {
     804         588 :   pari_sp av = avma;
     805         588 :   GEN b = leading_coeff(B), bi;
     806             :   GEN B2, P, V, W;
     807             :   long i, lV;
     808         588 :   if (typ(b)==t_INT) return FpXQX_digits(x, B, T, q);
     809         294 :   bi = ZpXQ_inv(b, T, p, e);
     810         294 :   B2 = FqX_Fq_mul_to_monic(B, bi, T, q);
     811         294 :   V = FpXQX_digits(x, B2, T, q);
     812         294 :   lV = lg(V)-1;
     813         294 :   P = FpXQ_powers(bi, lV-1, T, q);
     814         294 :   W = cgetg(lV+1, t_VEC);
     815        4816 :   for(i=1; i<=lV; i++)
     816        4522 :     gel(W, i) = FpXQX_FpXQ_mul(gel(V,i), gel(P, i), T, q);
     817         294 :   return gerepileupto(av, W);
     818             : }
     819             : 
     820             : struct _ZpXQ_sqrtn
     821             : {
     822             :   GEN T, a, n, ai;
     823             : };
     824             : 
     825             : static GEN
     826        1778 : _sqrtn_invd(void *E, GEN V, GEN v, GEN q, long M)
     827             : {
     828        1778 :   struct _ZpXQ_sqrtn *d = (struct _ZpXQ_sqrtn *) E;
     829        1778 :   GEN Tq = FpX_red(d->T, q), aiq = FpX_red(d->ai, q);
     830             :   (void)M;
     831        1778 :   return FpXQ_mul(FpXQ_mul(V, gel(v,2), Tq, q), aiq, Tq, q);
     832             : }
     833             : 
     834             : static GEN
     835        1778 : _sqrtn_eval(void *E, GEN x, GEN q)
     836             : {
     837        1778 :   struct _ZpXQ_sqrtn *d = (struct _ZpXQ_sqrtn *) E;
     838        1778 :   GEN Tq = FpX_red(d->T, q);
     839        1778 :   GEN f = FpX_sub(FpXQ_pow(x, d->n, Tq, q), d->a, q);
     840        1778 :   return mkvec2(f, x);
     841             : }
     842             : 
     843             : GEN
     844        1141 : ZpXQ_sqrtnlift(GEN a, GEN n, GEN x, GEN T, GEN p, long e)
     845             : {
     846             :   struct _ZpXQ_sqrtn d;
     847        1141 :   d.a = a; d.T = T; d.n = n;
     848        1141 :   d.ai = ZpXQ_inv(ZX_Z_mul(a, n),T,p,(e+1)>>1);
     849        1141 :   return gen_ZpX_Newton(x, p, e, &d, _sqrtn_eval, _sqrtn_invd);
     850             : }
     851             : 
     852             : GEN
     853           0 : ZpXQ_sqrt(GEN a, GEN T, GEN p, long e)
     854             : {
     855           0 :   pari_sp av = avma;
     856           0 :   GEN z = FpXQ_sqrt(FpX_red(a, p), T, p);
     857           0 :   if (!z) return NULL;
     858           0 :   if (e <= 1) return gerepileupto(av, z);
     859           0 :   return gerepileupto(av, ZpXQ_sqrtnlift(a, gen_2, z, T, p, e));
     860             : }
     861             : 
     862             : GEN
     863        2618 : ZpX_ZpXQ_liftroot_ea(GEN P, GEN S, GEN T, GEN p, long n, void *E,
     864             :                      int early(void *E, GEN x, GEN q))
     865             : {
     866        2618 :   pari_sp ltop = avma, av;
     867             :   long N, r;
     868             :   long mask;
     869             :   GEN q2, q, W, Q, Tq2, Tq, Pq;
     870             :   pari_timer ti;
     871        2618 :   T = FpX_get_red(T, powiu(p, n));
     872        2618 :   if (n == 1) return gcopy(S);
     873        2618 :   mask = quadratic_prec_mask(n);
     874        2618 :   av = avma;
     875        2618 :   q2 = p; q = sqri(p); mask >>= 1; N = 2;
     876        2618 :   if (DEBUGLEVEL > 3) timer_start(&ti);
     877        2618 :   Tq = FpXT_red(T,q);
     878        2618 :   Tq2 = FpXT_red(Tq,q2);
     879        2618 :   Pq = FpX_red(P,q);
     880        2618 :   W = FpXQ_inv(FpX_FpXQ_eval(FpX_deriv(P,q2), S, Tq2, q2), Tq2, q2);
     881        2618 :   Q  = ZX_Z_divexact(FpX_FpXQ_eval(Pq, S, Tq, q), q2);
     882        2618 :   r = brent_kung_optpow(degpol(P), 4, 3);
     883             :   for (;;)
     884             :   {
     885             :     GEN H, Sq, Wq, Spow, dP, qq, Pqq, Tqq;
     886        7992 :     H  = FpXQ_mul(W, Q, Tq2, q2);
     887        7992 :     Sq = FpX_sub(S, ZX_Z_mul(H, q2), q);
     888        7992 :     if (DEBUGLEVEL > 3)
     889           0 :       timer_printf(&ti,"ZpX_ZpXQ_liftroot: lift to prec %ld",N);
     890        7992 :     if (mask==1 || (early && early(E, Sq, q)))
     891        2618 :       return gerepileupto(ltop, Sq);
     892        5374 :     qq = sqri(q); N <<= 1;
     893        5374 :     if (mask&1UL) { qq = diviiexact(qq, p); N--; }
     894        5374 :     mask >>= 1;
     895        5374 :     Pqq  = FpX_red(P, qq);
     896        5374 :     Tqq  = FpXT_red(T, qq);
     897        5374 :     Spow = FpXQ_powers(Sq, r, Tqq, qq);
     898        5374 :     Q  = ZX_Z_divexact(FpX_FpXQV_eval(Pqq, Spow, Tqq, qq), q);
     899        5374 :     dP = FpX_FpXQV_eval(FpX_deriv(Pq, q), FpXV_red(Spow, q), Tq, q);
     900        5374 :     Wq = ZX_Z_divexact(FpX_Fp_sub(FpXQ_mul(W, dP, Tq, q), gen_1, q), q2);
     901        5374 :     Wq = ZX_Z_mul(FpXQ_mul(W, Wq, Tq2, q2), q2);
     902        5374 :     Wq = FpX_sub(W, Wq, q);
     903        5374 :     S = Sq; W = Wq; q2 = q; q = qq; Tq2 = Tq; Tq = Tqq; Pq = Pqq;
     904        5374 :     if (gc_needed(av, 1))
     905             :     {
     906           0 :       if(DEBUGMEM>1) pari_warn(warnmem,"ZpX_ZpXQ_Newton");
     907           0 :       gerepileall(av, 8, &S, &W, &Q, &Tq2, &Tq, &Pq, &q, &q2);
     908             :     }
     909        5374 :   }
     910             : }
     911             : 
     912             : GEN
     913         441 : ZpX_ZpXQ_liftroot(GEN P, GEN S, GEN T, GEN p, long n)
     914             : {
     915         441 :   return ZpX_ZpXQ_liftroot_ea(P, S, T, p, n, NULL, NULL);
     916             : }
     917             : 
     918             : GEN
     919         287 : ZpX_Frobenius(GEN T, GEN p, long e)
     920             : {
     921         287 :   return ZpX_ZpXQ_liftroot(get_FpX_mod(T), FpX_Frobenius(T, p), T, p, e);
     922             : }
     923             : 
     924             : GEN
     925         140 : ZpXQM_prodFrobenius(GEN M, GEN T, GEN p, long e)
     926             : {
     927         140 :   pari_sp av = avma;
     928         140 :   GEN xp = ZpX_Frobenius(T, p, e);
     929         140 :   GEN z = FpXQM_autsum(mkvec2(xp, M), get_FpX_degree(T), T, powiu(p,e));
     930         140 :   return gerepilecopy(av, gel(z,2));
     931             : }

Generated by: LCOV version 1.11