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 19608-b7b365e) Lines: 530 552 96.0 %
Date: 2016-09-26 05:54:43 Functions: 46 47 97.9 %
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       49132 : BuildTree(GEN link, GEN V, GEN W, GEN a, GEN T, GEN p)
      36             : {
      37       49132 :   long k = lg(a)-1;
      38             :   long i, j, s, minp, mind;
      39             : 
      40       49132 :   for (i=1; i<=k; i++) { gel(V,i) = gel(a,i); link[i] = -i; }
      41             : 
      42       78379 :   for (j=1; j <= 2*k-5; j+=2,i++)
      43             :   {
      44       29247 :     minp = j;
      45       29247 :     mind = degpol(gel(V,j));
      46      413055 :     for (s=j+1; s<i; s++)
      47      383808 :       if (degpol(gel(V,s)) < mind) { minp = s; mind = degpol(gel(V,s)); }
      48             : 
      49       29247 :     swap(gel(V,j), gel(V,minp));
      50       29247 :     lswap(link[j], link[minp]);
      51             : 
      52       29247 :     minp = j+1;
      53       29247 :     mind = degpol(gel(V,j+1));
      54      383808 :     for (s=j+2; s<i; s++)
      55      354561 :       if (degpol(gel(V,s)) < mind) { minp = s; mind = degpol(gel(V,s)); }
      56             : 
      57       29247 :     swap(gel(V,j+1), gel(V,minp));
      58       29247 :     lswap(link[j+1], link[minp]);
      59             : 
      60       29247 :     gel(V,i) = FqX_mul(gel(V,j), gel(V,j+1), T, p);
      61       29247 :     link[i] = j;
      62             :   }
      63             : 
      64      127504 :   for (j=1; j <= 2*k-3; j+=2)
      65             :   {
      66             :     GEN d, u, v;
      67       78379 :     d = FqX_extgcd(gel(V,j), gel(V,j+1), T, p, &u, &v);
      68       78379 :     if (degpol(d) > 0) pari_err_COPRIME("BuildTree", gel(V,j), gel(V,j+1));
      69       78372 :     d = gel(d,2);
      70       78372 :     if (!gequal1(d))
      71             :     {
      72       71734 :       if (typ(d)==t_POL)
      73             :       {
      74        5936 :         d = FpXQ_inv(d, T, p);
      75        5936 :         u = FqX_Fq_mul(u, d, T, p);
      76        5936 :         v = FqX_Fq_mul(v, d, T, p);
      77             :       }
      78             :       else
      79             :       {
      80       65798 :         d = Fp_inv(d, p);
      81       65798 :         u = FqX_Fp_mul(u, d, T,p);
      82       65798 :         v = FqX_Fp_mul(v, d, T,p);
      83             :       }
      84             :     }
      85       78372 :     gel(W,j) = u;
      86       78372 :     gel(W,j+1) = v;
      87             :   }
      88       49125 : }
      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      336746 : ZpX_HenselLift(GEN V, GEN W, long j, GEN f, GEN pd, GEN p0, GEN p1, int noinv)
      94             : {
      95      336746 :   pari_sp av = avma;
      96      336746 :   long space = lg(f) * lgefint(p1);
      97             :   GEN a2, b2, g, z, s, t;
      98      336746 :   GEN a = gel(V,j), b = gel(V,j+1);
      99      336746 :   GEN u = gel(W,j), v = gel(W,j+1);
     100             : 
     101      336746 :   (void)new_chunk(space); /* HACK */
     102      336746 :   g = ZX_sub(f, ZX_mul(a,b));
     103      336746 :   g = ZX_Z_divexact(g, p0);
     104      336746 :   g = FpX_red(g, pd);
     105      336746 :   z = FpX_mul(v,g, pd);
     106      336746 :   t = FpX_divrem(z,a, pd, &s);
     107      336746 :   t = ZX_add(ZX_mul(u,g), ZX_mul(t,b));
     108      336746 :   t = FpX_red(t, pd);
     109      336746 :   t = ZX_Z_mul(t,p0);
     110      336746 :   s = ZX_Z_mul(s,p0);
     111      336746 :   avma = av;
     112      336746 :   a2 = ZX_add(a,s);
     113      336746 :   b2 = ZX_add(b,t);
     114             : 
     115             :   /* already reduced mod p1 = pd p0 */
     116      336746 :   gel(V,j)   = a2;
     117      336746 :   gel(V,j+1) = b2;
     118      673492 :   if (noinv) return;
     119             : 
     120      267733 :   av = avma;
     121      267733 :   (void)new_chunk(space); /* HACK */
     122      267733 :   g = ZX_add(ZX_mul(u,a2), ZX_mul(v,b2));
     123      267733 :   g = Z_ZX_sub(gen_1, g);
     124      267733 :   g = ZX_Z_divexact(g, p0);
     125      267733 :   g = FpX_red(g, pd);
     126      267733 :   z = FpX_mul(v,g, pd);
     127      267733 :   t = FpX_divrem(z,a, pd, &s);
     128      267733 :   t = ZX_add(ZX_mul(u,g), ZX_mul(t,b));
     129      267733 :   t = FpX_red(t, pd);
     130      267733 :   t = ZX_Z_mul(t,p0);
     131      267733 :   s = ZX_Z_mul(s,p0);
     132      267733 :   avma = av;
     133      267733 :   gel(W,j)   = ZX_add(u,t);
     134      267733 :   gel(W,j+1) = ZX_add(v,s);
     135             : }
     136             : 
     137             : static void
     138       32949 : ZpXQ_HenselLift(GEN V, GEN W, long j, GEN f, GEN Td, GEN T1, GEN pd, GEN p0, GEN p1, int noinv)
     139             : {
     140       32949 :   pari_sp av = avma;
     141       32949 :   const long n = degpol(T1), vT = varn(T1);
     142       32949 :   long space = lg(f) * lgefint(p1) * lg(T1);
     143             :   GEN a2, b2, g, z, s, t;
     144       32949 :   GEN a = gel(V,j), b = gel(V,j+1);
     145       32949 :   GEN u = gel(W,j), v = gel(W,j+1);
     146             : 
     147       32949 :   (void)new_chunk(space); /* HACK */
     148       32949 :   g = RgX_sub(f, Kronecker_to_ZXX(ZXX_mul_Kronecker(a,b,n), n, vT));
     149       32949 :   g = FpXQX_red(g, T1, p1);
     150       32949 :   g = RgX_Rg_divexact(g, p0);
     151       32949 :   z = FpXQX_mul(v,g, Td,pd);
     152       32949 :   t = FpXQX_divrem(z,a, Td,pd, &s);
     153       32949 :   t = ZX_add(ZXX_mul_Kronecker(u,g,n), ZXX_mul_Kronecker(t,b,n));
     154       32949 :   t = Kronecker_to_ZXX(t, n, vT);
     155       32949 :   t = FpXQX_red(t, Td, pd);
     156       32949 :   t = RgX_Rg_mul(t,p0);
     157       32949 :   s = RgX_Rg_mul(s,p0);
     158       32949 :   avma = av;
     159             : 
     160       32949 :   a2 = RgX_add(a,s);
     161       32949 :   b2 = RgX_add(b,t);
     162             :   /* already reduced mod p1 = pd p0 */
     163       32949 :   gel(V,j)   = a2;
     164       32949 :   gel(V,j+1) = b2;
     165       65898 :   if (noinv) return;
     166             : 
     167       26082 :   av = avma;
     168       26082 :   (void)new_chunk(space); /* HACK */
     169       26082 :   g = ZX_add(ZXX_mul_Kronecker(u,a2,n), ZXX_mul_Kronecker(v,b2,n));
     170       26082 :   g = Kronecker_to_ZXX(g, n, vT);
     171       26082 :   g = Rg_RgX_sub(gen_1, g);
     172       26082 :   g = FpXQX_red(g, T1, p1);
     173       26082 :   g = RgX_Rg_divexact(g, p0);
     174       26082 :   z = FpXQX_mul(v,g, Td,pd);
     175       26082 :   t = FpXQX_divrem(z,a, Td,pd, &s);
     176       26082 :   t = ZX_add(ZXX_mul_Kronecker(u,g,n), ZXX_mul_Kronecker(t,b,n));
     177       26082 :   t = Kronecker_to_ZXX(t, n, vT);
     178       26082 :   t = FpXQX_red(t, Td, pd);
     179       26082 :   t = RgX_Rg_mul(t,p0);
     180       26082 :   s = RgX_Rg_mul(s,p0);
     181       26082 :   avma = av;
     182       26082 :   gel(W,j)   = RgX_add(u,t);
     183       26082 :   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      866030 : ZpX_RecTreeLift(GEN link, GEN v, GEN w, GEN pd, GEN p0, GEN p1,
     190             :                 GEN f, long j, int noinv)
     191             : {
     192     1732060 :   if (j < 0) return;
     193      336746 :   ZpX_HenselLift(v, w, j, f, pd, p0,p1, noinv);
     194      336746 :   ZpX_RecTreeLift(link, v, w, pd, p0,p1, gel(v,j)  , link[j  ], noinv);
     195      336746 :   ZpX_RecTreeLift(link, v, w, pd, p0,p1, gel(v,j+1), link[j+1], noinv);
     196             : }
     197             : static void
     198       69580 : 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      139160 :   if (j < 0) return;
     202       32949 :   ZpXQ_HenselLift(v, w, j, f, Td,T1, pd, p0,p1, noinv);
     203       32949 :   ZpXQ_RecTreeLift(link, v, w, Td,T1, pd, p0,p1, gel(v,j)  , link[j  ], noinv);
     204       32949 :   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      301693 : quadratic_prec_mask(long n)
     225             : {
     226      301693 :   long a = n, i;
     227      301693 :   ulong mask = 0;
     228      864456 :   for(i = 1;; i++, mask <<= 1)
     229             :   {
     230      864456 :     mask |= (a&1); a = (a+1)>>1;
     231     1166149 :     if (a==1) return mask | (1UL << i);
     232      562763 :   }
     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       49153 : MultiLift(GEN f, GEN a, GEN T, GEN p, long e0, long flag)
     242             : {
     243       49153 :   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       49153 :   if (k < 2) pari_err_DOMAIN("MultiLift", "#(modular factors)", "<", gen_2,a);
     249       49153 :   if (e0 < 1) pari_err_DOMAIN("MultiLift", "precision", "<", gen_1,stoi(e0));
     250       49153 :   if (e0 == 1) return a;
     251             : 
     252       49132 :   if (DEBUGLEVEL > 3) timer_start(&Ti);
     253       49132 :   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       49132 :     e = 1;
     263       49132 :     v = cgetg(2*k-2 + 1, t_VEC);
     264       49132 :     w = cgetg(2*k-2 + 1, t_VEC);
     265       49132 :     link=cgetg(2*k-2 + 1, t_VECSMALL);
     266       49132 :     BuildTree(link, v, w, a, T? FpX_red(T,p): NULL, p);
     267       49125 :     if (DEBUGLEVEL > 3) timer_printf(&Ti, "building tree");
     268             :   }
     269       49125 :   mask = quadratic_prec_mask(e0);
     270       49125 :   eold = 1;
     271       49125 :   penew = NULL;
     272       49125 :   Tnew = NULL;
     273      294470 :   while (mask > 1)
     274             :   {
     275      196220 :     long enew = eold << 1;
     276      196220 :     if (mask & 1) enew--;
     277      196220 :     mask >>= 1;
     278      196220 :     if (enew >= e) { /* mask == 1: last iteration */
     279      196220 :       GEN peold = penew? penew: powiu(p, eold);
     280      196220 :       GEN Td = NULL, pd;
     281      196220 :       long d = enew - eold; /* = eold or eold-1 */
     282             :       /* lift from p^eold to p^enew */
     283      196220 :       pd = (d == eold)? peold: diviiexact(peold, p); /* p^d */
     284      196220 :       penew = mulii(peold,pd);
     285      196220 :       if (T) {
     286        3682 :         if (Tnew)
     287        2716 :           Td = (d == eold)? Tnew: FpX_red(Tnew,pd);
     288             :         else
     289         966 :           Td = FpX_red(T, peold);
     290        3682 :         Tnew = FpX_red(T, penew);
     291        3682 :         ZpXQ_RecTreeLift(link, v, w, Td, Tnew, pd, peold, penew, f, lg(v)-2,
     292             :                          (flag == 0 && mask == 1));
     293             :       }
     294             :       else
     295      192538 :         ZpX_RecTreeLift(link, v, w, pd, peold, penew, f, lg(v)-2,
     296             :                         (flag == 0 && mask == 1));
     297      196220 :       if (DEBUGLEVEL > 3) timer_printf(&Ti, "lifting to prec %ld", enew);
     298             :     }
     299      196220 :     eold = enew;
     300             :   }
     301             : 
     302       49125 :   if (flag)
     303         595 :     E = mkvec4(utoipos(e0), link, v, w);
     304             :   else
     305             :   {
     306       48530 :     E = cgetg(k+1, t_VEC);
     307      200290 :     for (i = 1; i <= 2*k-2; i++)
     308             :     {
     309      151760 :       long t = link[i];
     310      151760 :       if (t < 0) gel(E,-t) = gel(v,i);
     311             :     }
     312             :   }
     313       49125 :   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       86610 : ZpX_liftfact(GEN pol, GEN Q, GEN pe, GEN p, long e)
     320             : {
     321       86610 :   pari_sp av = avma;
     322       86610 :   pol = FpX_normalize(pol, pe);
     323       86610 :   if (lg(Q) == 2) return mkvec(pol);
     324       47585 :   return gerepilecopy(av, MultiLift(pol, Q, NULL, p, e, 0));
     325             : }
     326             : 
     327             : GEN
     328         973 : ZpXQX_liftfact(GEN pol, GEN Q, GEN T, GEN pe, GEN p, long e)
     329             : {
     330         973 :   pari_sp av = avma;
     331         973 :   pol = FpXQX_normalize(pol, T, pe);
     332         973 :   if (lg(Q) == 2) return mkvec(pol);
     333         973 :   return gerepilecopy(av, MultiLift(pol, Q, T, p, e, 0));
     334             : }
     335             : 
     336             : GEN
     337        1582 : ZqX_liftfact(GEN f, GEN a, GEN T, GEN pe, GEN p, long e)
     338        1582 : { return T ? ZpXQX_liftfact(f, a, T, pe, p, e): ZpX_liftfact(f, a, pe, p, e); }
     339             : GEN
     340          70 : ZqX_liftroot(GEN f, GEN a, GEN T, GEN p, long e)
     341          70 : { return T ? ZpXQX_liftroot(f, a,T , p, e): ZpX_liftroot(f, a, p, e); }
     342             : 
     343             : 
     344             : /* U = NULL treated as 1 */
     345             : static void
     346        5579 : BezoutPropagate(GEN link, GEN v, GEN w, GEN pe, GEN U, GEN f, long j)
     347             : {
     348             :   GEN Q, R;
     349       11158 :   if (j < 0) return;
     350             : 
     351        2492 :   Q = FpX_mul(gel(v,j), gel(w,j), pe);
     352        2492 :   if (U)
     353             :   {
     354        1897 :     Q = FpXQ_mul(Q, U, f, pe);
     355        1897 :     R = FpX_sub(U, Q, pe);
     356             :   }
     357             :   else
     358         595 :     R = Fp_FpX_sub(gen_1, Q, pe);
     359        2492 :   gel(w,j+1) = Q; /*  0 mod U v[j],  1 mod (1-U) v[j+1] */
     360        2492 :   gel(w,j) = R; /*  1 mod U v[j],  0 mod (1-U) v[j+1] */
     361        2492 :   BezoutPropagate(link, v, w, pe, R, f, link[j  ]);
     362        2492 :   BezoutPropagate(link, v, w, pe, Q, f, link[j+1]);
     363             : }
     364             : 
     365             : /* as above, but return the Bezout coefficients for the lifted modular factors
     366             :  *   U[i] = 1 mod Qlift[i]
     367             :  *          0 mod Qlift[j], j != i */
     368             : GEN
     369         602 : bezout_lift_fact(GEN pol, GEN Q, GEN p, long e)
     370             : {
     371         602 :   pari_sp av = avma;
     372             :   GEN E, link, v, w, pe;
     373         602 :   long i, k = lg(Q)-1;
     374         602 :   if (k == 1) return mkvec(pol);
     375         595 :   pe = powiu(p, e);
     376         595 :   pol = FpX_normalize(pol, pe);
     377         595 :   E = MultiLift(pol, Q, NULL, p, e, 1);
     378         595 :   link = gel(E,2);
     379         595 :   v    = gel(E,3);
     380         595 :   w    = gel(E,4);
     381         595 :   BezoutPropagate(link, v, w, pe, NULL, pol, lg(v)-2);
     382         595 :   E = cgetg(k+1, t_VEC);
     383        5579 :   for (i = 1; i <= 2*k-2; i++)
     384             :   {
     385        4984 :     long t = link[i];
     386        4984 :     if (t < 0) E[-t] = w[i];
     387             :   }
     388         595 :   return gerepilecopy(av, E);
     389             : }
     390             : 
     391             : /* Front-end for ZpX_liftfact:
     392             :    lift the factorization of pol mod p given by L to p^N (if possible) */
     393             : GEN
     394          21 : polhensellift(GEN pol, GEN L, GEN p, long N)
     395             : {
     396          21 :   GEN T = NULL;
     397             :   long i, l, t;
     398          21 :   pari_sp av = avma;
     399             : 
     400          21 :   if (typ(pol) != t_POL) pari_err_TYPE("polhensellift",pol);
     401          21 :   RgX_check_ZXX(pol, "polhensellift");
     402          21 :   if (!is_vec_t(typ(L)) || lg(L) < 3) pari_err_TYPE("polhensellift",L);
     403          21 :   t = typ(p);
     404          21 :   if (t == t_VEC) /* [p, T] */
     405             :   {
     406           7 :     T = gel(p,2);
     407           7 :     if (typ(T) != t_POL) pari_err_TYPE("polhensellift",pol);
     408           7 :     RgX_check_ZX(T, "polhensellift");
     409           7 :     p = gel(p,1); t = typ(p);
     410             :   }
     411          21 :   if (t != t_INT) pari_err_TYPE("polhensellift",p);
     412          21 :   if (N < 1) pari_err_DOMAIN("polhensellift", "precision", "<", gen_1,stoi(N));
     413             : 
     414          21 :   l = lg(L); L = leafcopy(L);
     415          63 :   for (i = 1; i < l; i++)
     416             :   {
     417          42 :     if (typ(gel(L,i)) != t_POL)
     418           0 :       gel(L,i) = scalar_ZX_shallow(gel(L,i), varn(pol));
     419          42 :     RgX_check_ZXX(gel(L,i), "polhensellift");
     420             :   }
     421          21 :   return gerepilecopy(av, ZqX_liftfact(pol, L, T, powiu(p,N), p, N));
     422             : }
     423             : 
     424             : static GEN
     425       41143 : FqV_roots_from_deg1(GEN x, GEN T, GEN p)
     426             : {
     427       41143 :   long i,l = lg(x);
     428       41143 :   GEN r = cgetg(l,t_VEC);
     429       41143 :   for (i=1; i<l; i++) { GEN P = gel(x,i); gel(r,i) = Fq_neg(gel(P,2), T, p); }
     430       41143 :   return r;
     431             : }
     432             : 
     433             : static GEN
     434       41115 : ZpX_liftroots_full(GEN f, GEN S, GEN q, GEN p, long e)
     435             : {
     436       41115 :   pari_sp av = avma;
     437       41115 :   GEN y = ZpX_liftfact(f, deg1_from_roots(S, varn(f)), q, p, e);
     438       41115 :   return gerepileupto(av, FqV_roots_from_deg1(y, NULL, q));
     439             : }
     440             : 
     441             : GEN
     442       41049 : ZpX_roots(GEN F, GEN p, long e)
     443             : {
     444       41049 :   pari_sp av = avma;
     445       41049 :   GEN q = powiu(p, e);
     446       41049 :   GEN f = FpX_normalize(F, p);
     447       41049 :   GEN g = FpX_normalize(FpX_split_part(f, p), p);
     448             :   GEN S;
     449       41049 :   if (degpol(g) < degpol(f))
     450             :   {
     451       39152 :     GEN h = FpX_div(f, g, p);
     452       39152 :     F = gel(ZpX_liftfact(F, mkvec2(g, h), q, p, e), 1);
     453             :   }
     454       41049 :   S = FpX_roots(g, p);
     455       41049 :   return gerepileupto(av, ZpX_liftroots_full(F, S, q, p, e));
     456             : }
     457             : 
     458             : static GEN
     459          28 : ZpXQX_liftroots_full(GEN f, GEN S, GEN T, GEN q, GEN p, long e)
     460             : {
     461          28 :   pari_sp av = avma;
     462          28 :   GEN y = ZpXQX_liftfact(f, deg1_from_roots(S, varn(f)), T, q, p, e);
     463          28 :   return gerepileupto(av, FqV_roots_from_deg1(y, T, q));
     464             : }
     465             : 
     466             : GEN
     467          28 : ZpXQX_roots(GEN F, GEN T, GEN p, long e)
     468             : {
     469          28 :   pari_sp av = avma;
     470          28 :   GEN q = powiu(p, e);
     471          28 :   GEN f = FpXQX_normalize(F, T, p);
     472          28 :   GEN g = FpXQX_normalize(FpXQX_split_part(f, T, p), T, p);
     473             :   GEN S;
     474          28 :   if (degpol(g) < degpol(f))
     475             :   {
     476           7 :     GEN h = FpXQX_div(f, g, T, p);
     477           7 :     F = gel(ZpXQX_liftfact(F, mkvec2(g, h), T, q, p, e), 1);
     478             :   }
     479          28 :   S = FpXQX_roots(g, T, p);
     480          28 :   return gerepileupto(av, ZpXQX_liftroots_full(F, S, T, q, p, e));
     481             : }
     482             : 
     483             : GEN
     484         589 : ZqX_roots(GEN F, GEN T, GEN p, long e)
     485             : {
     486         589 :   return T ? ZpXQX_roots(F, T, p, e): ZpX_roots(F, p, e);
     487             : }
     488             : 
     489             : /* SPEC:
     490             :  * p is a t_INT > 1, e >= 1
     491             :  * f is a ZX with leading term prime to p.
     492             :  * a is a simple root mod l for all l|p.
     493             :  * Return roots of f mod p^e, as integers (implicitly mod p^e)
     494             :  * STANDARD USE: p is a prime power */
     495             : GEN
     496         364 : ZpX_liftroot(GEN f, GEN a, GEN p, long e)
     497             : {
     498         364 :   pari_sp av = avma;
     499         364 :   GEN q = p, fr, W;
     500             :   ulong mask;
     501             : 
     502         364 :   a = modii(a,q);
     503         364 :   if (e == 1) return a;
     504         364 :   mask = quadratic_prec_mask(e);
     505         364 :   fr = FpX_red(f,q);
     506         364 :   W = Fp_inv(FpX_eval(ZX_deriv(fr), a, q), q); /* 1/f'(a) mod p */
     507             :   for(;;)
     508             :   {
     509        1372 :     q = sqri(q);
     510        1372 :     if (mask & 1) q = diviiexact(q, p);
     511        1372 :     mask >>= 1;
     512        1372 :     fr = FpX_red(f,q);
     513        1372 :     a = Fp_sub(a, Fp_mul(W, FpX_eval(fr, a,q), q), q);
     514        1372 :     if (mask == 1) return gerepileuptoint(av, a);
     515        1008 :     W = Fp_sub(shifti(W,1), Fp_mul(Fp_sqr(W,q), FpX_eval(ZX_deriv(fr),a,q), q), q);
     516        1008 :   }
     517             : }
     518             : 
     519             : GEN
     520          66 : ZpX_liftroots(GEN f, GEN S, GEN p, long e)
     521             : {
     522          66 :   long i, n = lg(S)-1, d = degpol(f);
     523             :   GEN r;
     524          66 :   if (n == d) return ZpX_liftroots_full(f, S, powiu(p, e), p, e);
     525           0 :   r = cgetg(n+1, typ(S));
     526           0 :   for (i=1; i <= n; i++)
     527           0 :     gel(r,i) = ZpX_liftroot(f, gel(S,i), p, e);
     528           0 :   return r;
     529             : }
     530             : 
     531             : GEN
     532          91 : ZpXQX_liftroot_vald(GEN f, GEN a, long v, GEN T, GEN p, long e)
     533             : {
     534          91 :   pari_sp av = avma, av2;
     535          91 :   GEN pv = p, q, W, df, Tq, fr, dfr;
     536             :   ulong mask;
     537          91 :   a = Fq_red(a, T, p);
     538          91 :   if (e <= v+1) return a;
     539          91 :   df = RgX_deriv(f);
     540          91 :   if (v) { pv = powiu(p,v); df = ZXX_Z_divexact(df, pv); }
     541          91 :   mask = quadratic_prec_mask(e-v);
     542          91 :   Tq = FpXT_red(T, p); dfr = FpXQX_red(df, Tq, p);
     543          91 :   W = Fq_inv(FqX_eval(dfr, a, Tq, p), Tq, p); /* 1/f'(a) mod (T,p) */
     544          91 :   q = p;
     545          91 :   av2 = avma;
     546             :   for (;;)
     547             :   {
     548             :     GEN u, fa, qv, q2v, q2, Tq2;
     549         161 :     q2 = q; q = sqri(q);
     550         161 :     if (mask & 1) q = diviiexact(q,p);
     551         161 :     mask >>= 1;
     552         161 :     if (v) { qv = mulii(q, pv); q2v = mulii(q2, pv); }
     553         161 :     else { qv = q; q2v = q2; }
     554         161 :     Tq2 = FpXT_red(T, q2v); Tq = FpXT_red(T, qv);
     555         161 :     fr = FpXQX_red(f, Tq, qv);
     556         161 :     fa = FqX_eval(fr, a, Tq, qv);
     557         161 :     fa = typ(fa)==t_INT? diviiexact(fa,q2v): ZX_Z_divexact(fa, q2v);
     558         161 :     a = Fq_sub(a, ZX_Z_mul(Fq_mul(W, fa, Tq2, q2v), q2), Tq, qv);
     559         161 :     if (mask == 1) return gerepileupto(av, a);
     560          70 :     dfr = FpXQX_red(df, Tq, q);
     561          70 :     u = ZX_Z_divexact(FpX_Fp_sub(Fq_mul(W,FqX_eval(dfr,a,Tq,q),Tq,q),gen_1,q),q2);
     562          70 :     W = Fq_sub(W,ZX_Z_mul(Fq_mul(u,W,Tq2,q2),q2),Tq,q);
     563          70 :     if (gc_needed(av2,2))
     564             :     {
     565           0 :       if(DEBUGMEM>1) pari_warn(warnmem,"ZpXQX_liftroot, e = %ld", e);
     566           0 :       gerepileall(av2, 3, &a, &W, &q);
     567             :     }
     568          70 :   }
     569             : }
     570             : 
     571             : GEN
     572          91 : ZpXQX_liftroot(GEN f, GEN a, GEN T, GEN p, long e) { return ZpXQX_liftroot_vald(f,a,0,T,p,e); }
     573             : 
     574             : /* Same as ZpX_liftroot for the polynomial X^n-b*/
     575             : GEN
     576       47419 : Zp_sqrtnlift(GEN b, GEN n, GEN a, GEN p, long e)
     577             : {
     578       47419 :   pari_sp ltop=avma;
     579             :   GEN q, w, n_1;
     580             :   ulong mask;
     581       47419 :   long pis2 = equalii(n, gen_2)? 1: 0;
     582       47419 :   if (e == 1) return icopy(a);
     583       44899 :   n_1 = subis(n,1);
     584       44899 :   mask = quadratic_prec_mask(e);
     585       44899 :   w = Fp_inv(pis2 ? shifti(a,1): Fp_mul(n,Fp_pow(a,n_1,p), p), p);
     586       44899 :   q = p;
     587             :   for(;;)
     588             :   {
     589      106600 :     q = sqri(q);
     590      106600 :     if (mask & 1) q = diviiexact(q, p);
     591      106600 :     mask >>= 1;
     592      106600 :     if (lgefint(q) == 3 && lgefint(n) == 3)
     593       61551 :     {
     594      106090 :       ulong Q = uel(q,2), N = uel(n,2);
     595      106090 :       ulong A = umodiu(a, Q);
     596      106090 :       ulong B = umodiu(b, Q);
     597      106090 :       ulong W = umodiu(w, Q);
     598             : 
     599      106090 :       A = Fl_sub(A, Fl_mul(W, Fl_sub(Fl_powu(A,N,Q), B, Q), Q), Q);
     600      106090 :       a = utoi(A);
     601      106090 :       if (mask == 1) break;
     602       61551 :       W = Fl_sub(Fl_add(W,W,Q),
     603             :                  Fl_mul(Fl_sqr(W,Q), Fl_mul(N,Fl_powu(A, N-1, Q), Q), Q), Q);
     604       61551 :       w = utoi(W);
     605             :     }
     606             :     else
     607             :     {
     608             :       /* a -= w (a^n - b) */
     609         510 :       a = modii(subii(a, mulii(w, subii(Fp_pow(a,n,q),b))), q);
     610         510 :       if (mask == 1) break;
     611             :       /* w += w - w^2 n a^(n-1)*/
     612         285 :       w = subii(shifti(w,1), Fp_mul(Fp_sqr(w,q),
     613         135 :                            pis2? shifti(a,1): mulii(n,Fp_pow(a,n_1,q)), q));
     614             :     }
     615       61701 :   }
     616       44899 :   return gerepileuptoint(ltop,a);
     617             : }
     618             : 
     619             : 
     620             : /* Same as ZpX_liftroot for the polynomial X^2-b */
     621             : GEN
     622       12929 : Zp_sqrtlift(GEN b, GEN a, GEN p, long e)
     623             : {
     624       12929 :   return Zp_sqrtnlift(b, gen_2, a, p, e);
     625             : }
     626             : 
     627             : GEN
     628        3185 : Zp_sqrt(GEN x, GEN p, long e)
     629             : {
     630             :   pari_sp av;
     631             :   GEN z;
     632        3185 :   if (absequaliu(p,2)) return Z2_sqrt(x,e);
     633        1708 :   av = avma;
     634        1708 :   z = Fp_sqrt(Fp_red(x, p), p);
     635        1708 :   if (!z) return NULL;
     636        1673 :   if (e > 1) z = Zp_sqrtlift(x, z, p, e);
     637        1673 :   return gerepileuptoint(av, z);
     638             : }
     639             : 
     640             : /* Compute (x-1)/(x+1)/p^k */
     641             : static GEN
     642       20781 : ZpXQ_log_to_ath(GEN x, long k, GEN T, GEN p, long e, GEN pe)
     643             : {
     644       20781 :   pari_sp av = avma;
     645       20781 :   long vT = get_FpX_var(T);
     646             :   GEN bn, bdi;
     647       20781 :   GEN bd = ZX_Z_add(x, gen_1);
     648       20781 :   if (absequaliu(p,2)) /*For p=2, we need to simplify by 2*/
     649             :   {
     650        6851 :     bn = ZX_shifti(x,-(k+1));
     651        6851 :     bdi= ZpXQ_invlift(ZX_shifti(bd ,-1), pol_1(vT), T, p, e);
     652             :   }
     653             :   else
     654             :   {
     655       13930 :     bn = ZX_Z_divexact(ZX_Z_sub(x, gen_1),powiu(p,k));
     656       13930 :     bdi= ZpXQ_invlift(bd, scalarpol(Fp_inv(gen_2,p),vT), T, p, e);
     657             :   }
     658       20781 :   return gerepileupto(av, FpXQ_mul(bn, bdi, T, pe));
     659             : }
     660             : 
     661             : /* Assume p odd, a = 1 [p], return log(a) */
     662             : GEN
     663       20781 : ZpXQ_log(GEN a, GEN T, GEN p, long N)
     664             : {
     665       20781 :   pari_sp av = avma;
     666             :   pari_timer ti;
     667       20781 :   long is2 = absequaliu(p,2);
     668       20781 :   ulong pp = is2 ? 0: itou_or_0(p);
     669       20781 :   double lp = is2 ? 1: pp ? log2(pp): expi(p);
     670       20781 :   long k = maxss(1 , (long) .5+pow((double)(N>>1)/(lp*lp), 1./3));
     671             :   GEN ak, s, b, pol;
     672       20781 :   long e = is2 ? N-1: N;
     673       20781 :   long i, l = (e-2)/(2*(k+is2));
     674       20781 :   GEN pe = powiu(p,e);
     675       20781 :   GEN TNk, pNk = powiu(p,N+k);
     676       20781 :   if( DEBUGLEVEL>=3) timer_start(&ti);
     677       20781 :   TNk = FpX_get_red(get_FpX_mod(T), pNk);
     678       20781 :   ak = FpXQ_pow(a, powiu(p,k), TNk, pNk);
     679       20781 :   if( DEBUGLEVEL>=3) timer_printf(&ti,"FpXQ_pow(%ld)",k);
     680       20781 :   b = ZpXQ_log_to_ath(ak, k, T, p, e, pe);
     681       20781 :   if( DEBUGLEVEL>=3) timer_printf(&ti,"ZpXQ_log_to_ath");
     682       20781 :   pol= cgetg(l+3,t_POL);
     683       20781 :   pol[1] = evalsigne(1)|evalvarn(0);
     684       55400 :   for(i=0; i<=l; i++)
     685             :   {
     686             :     GEN g;
     687       34619 :     ulong z = 2*i+1;
     688       34619 :     if (pp)
     689             :     {
     690       26138 :       long w = u_lvalrem(z, pp, &z);
     691       26138 :       g = powuu(pp,2*i*k-w);
     692             :     }
     693        8481 :     else g = powiu(p,2*i*k);
     694       34620 :     gel(pol,i+2) = Fp_div(g, utoi(z),pe);
     695             :   }
     696       20781 :   if( DEBUGLEVEL>=3) timer_printf(&ti,"pol(%ld)",l);
     697       20781 :   s = FpX_FpXQ_eval(pol, FpXQ_sqr(b, T, pe), T,  pe);
     698       20781 :   if( DEBUGLEVEL>=3) timer_printf(&ti,"FpX_FpXQ_eval");
     699       20781 :   s = ZX_shifti(FpXQ_mul(b, s, T, pe), 1);
     700       20781 :   return gerepileupto(av, is2? s: FpX_red(s, pe));
     701             : }
     702             : 
     703             : /***********************************************************************/
     704             : /**                                                                   **/
     705             : /**                 Generic quadratic hensel lift over Zp[X]          **/
     706             : /**                                                                   **/
     707             : /***********************************************************************/
     708             : /* q = p^N */
     709             : 
     710             : GEN
     711      198380 : gen_ZpX_Dixon(GEN F, GEN V, GEN q, GEN p, long N, void *E,
     712             :                             GEN lin(void *E, GEN F, GEN d, GEN q),
     713             :                             GEN invl(void *E, GEN d))
     714             : {
     715      198380 :   pari_sp av = avma;
     716             :   long N2, M;
     717             :   GEN VN2, V2, VM, bil;
     718             :   GEN q2, qM;
     719      198380 :   V = FpX_red(V, q);
     720      198380 :   if (N == 1) return invl(E, V);
     721       44842 :   N2 = (N + 1)>>1; M = N - N2;
     722       44842 :   F = FpXT_red(F, q);
     723       44842 :   qM = powiu(p, M);
     724       44842 :   q2 = M == N2? qM: mulii(qM, p);
     725             :   /* q2 = p^N2, qM = p^M, q = q2 * qM */
     726       44842 :   VN2 = gen_ZpX_Dixon(F, V, q2, p, N2, E, lin, invl);
     727       44842 :   bil = lin(E, F, VN2, q);
     728       44842 :   V2 = ZX_Z_divexact(ZX_sub(V, bil), q2);
     729       44842 :   VM = gen_ZpX_Dixon(F, V2, qM, p, M, E, lin, invl);
     730       44842 :   return gerepileupto(av, FpX_red(ZX_add(VN2, ZX_Z_mul(VM, q2)), q));
     731             : }
     732             : 
     733             : GEN
     734      180650 : gen_ZpX_Newton(GEN x, GEN p, long n, void *E,
     735             :                       GEN eval(void *E, GEN f, GEN q),
     736             :                       GEN invd(void *E, GEN V, GEN v, GEN q, long M))
     737             : {
     738      180650 :   pari_sp ltop = avma, av;
     739      180650 :   long N = 1, N2, M;
     740             :   long mask;
     741      180650 :   GEN q = p;
     742      180650 :   if (n == 1) return gcopy(x);
     743      178676 :   mask = quadratic_prec_mask(n);
     744      178676 :   av = avma;
     745      789263 :   while (mask > 1)
     746             :   {
     747             :     GEN qM, q2, v, V;
     748      431911 :     N2 = N; N <<= 1;
     749      431911 :     q2 = q;
     750      431911 :     if (mask&1UL) { /* can never happen when q2 = p */
     751      162907 :       N--; M = N2-1;
     752      162907 :       qM = diviiexact(q2,p); /* > 1 */
     753      162906 :       q = mulii(qM,q2);
     754             :     } else {
     755      269004 :       M = N2;
     756      269004 :       qM = q2;
     757      269004 :       q = sqri(q2);
     758             :     }
     759             :     /* q2 = p^N2, qM = p^M, q = p^N = q2 * qM */
     760      431910 :     mask >>= 1;
     761      431910 :     v = eval(E, x, q);
     762      431910 :     V = ZX_Z_divexact(gel(v,1), q2);
     763      431911 :     x = FpX_sub(x, ZX_Z_mul(invd(E, V, v, qM, M), q2), q);
     764      431910 :     if (gc_needed(av, 1))
     765             :     {
     766           0 :       if(DEBUGMEM>1) pari_warn(warnmem,"gen_ZpX_Newton");
     767           0 :       gerepileall(av, 2, &x, &q);
     768             :     }
     769             :   }
     770      178676 :   return gerepileupto(ltop, x);
     771             : }
     772             : 
     773             : struct _ZpXQ_inv
     774             : {
     775             :   GEN T, a, p ,n;
     776             : };
     777             : 
     778             : static GEN
     779      291494 : _inv_invd(void *E, GEN V, GEN v, GEN q, long M/*unused*/)
     780             : {
     781      291494 :   struct _ZpXQ_inv *d = (struct _ZpXQ_inv *) E;
     782      291494 :   GEN Tq = FpXT_red(d->T, q);
     783             :   (void)M;
     784      291494 :   return FpXQ_mul(V, gel(v,2), Tq, q);
     785             : }
     786             : 
     787             : static GEN
     788      291494 : _inv_eval(void *E, GEN x, GEN q)
     789             : {
     790      291494 :   struct _ZpXQ_inv *d = (struct _ZpXQ_inv *) E;
     791      291494 :   GEN Tq = FpXT_red(d->T, q);
     792      291494 :   GEN f = FpX_Fp_sub(FpXQ_mul(x, FpX_red(d->a, q), Tq, q), gen_1, q);
     793      291494 :   return mkvec2(f, x);
     794             : }
     795             : 
     796             : GEN
     797      116254 : ZpXQ_invlift(GEN a, GEN x, GEN T, GEN p, long e)
     798             : {
     799             :   struct _ZpXQ_inv d;
     800      116254 :   d.a = a; d.T = T; d.p = p;
     801      116254 :   return gen_ZpX_Newton(x, p, e, &d, _inv_eval, _inv_invd);
     802             : }
     803             : 
     804             : GEN
     805       95473 : ZpXQ_inv(GEN a, GEN T, GEN p, long e)
     806             : {
     807       95473 :   pari_sp av=avma;
     808             :   GEN ai;
     809       95473 :   if (lgefint(p)==3)
     810             :   {
     811       95445 :     ulong pp = p[2];
     812       95445 :     ai = Flx_to_ZX(Flxq_inv(ZX_to_Flx(a,pp), ZXT_to_FlxT(T, pp), pp));
     813             :   } else
     814          28 :     ai = FpXQ_inv(FpX_red(a,p), FpXT_red(T,p),p);
     815       95473 :   return gerepileupto(av, ZpXQ_invlift(a, ai, T, p, e));
     816             : }
     817             : 
     818             : GEN
     819       31577 : ZpXQ_div(GEN a, GEN b, GEN T, GEN q, GEN p, long e)
     820             : {
     821       31577 :   return FpXQ_mul(a, ZpXQ_inv(b, T, p, e), T, q);
     822             : }
     823             : 
     824             : GEN
     825      142422 : ZpXQX_divrem(GEN x, GEN Sp, GEN T, GEN q, GEN p, long e, GEN *pr)
     826             : {
     827      142422 :   pari_sp av = avma;
     828      142422 :   GEN S = get_FpXQX_mod(Sp);
     829      142422 :   GEN b = leading_coeff(S), bi;
     830             :   GEN S2, Q;
     831      142422 :   if (typ(b)==t_INT) return FpXQX_divrem(x, Sp, T, q, pr);
     832       61796 :   bi = ZpXQ_inv(b, T, p, e);
     833       61796 :   S2 = FqX_Fq_mul_to_monic(S, bi, T, q);
     834       61796 :   Q = FpXQX_divrem(x, S2, T, q, pr);
     835       61796 :   if (pr==ONLY_DIVIDES && !Q) { avma = av; return NULL; }
     836       61796 :   if (pr==ONLY_REM || pr==ONLY_DIVIDES) return gerepileupto(av, Q);
     837       61796 :   Q = FpXQX_FpXQ_mul(Q, bi, T, q);
     838       61796 :   gerepileall(av, 2, &Q, pr);
     839       61796 :   return Q;
     840             : }
     841             : 
     842             : GEN
     843         588 : ZpXQX_digits(GEN x, GEN B, GEN T, GEN q, GEN p, long e)
     844             : {
     845         588 :   pari_sp av = avma;
     846         588 :   GEN b = leading_coeff(B), bi;
     847             :   GEN B2, P, V, W;
     848             :   long i, lV;
     849         588 :   if (typ(b)==t_INT) return FpXQX_digits(x, B, T, q);
     850         294 :   bi = ZpXQ_inv(b, T, p, e);
     851         294 :   B2 = FqX_Fq_mul_to_monic(B, bi, T, q);
     852         294 :   V = FpXQX_digits(x, B2, T, q);
     853         294 :   lV = lg(V)-1;
     854         294 :   P = FpXQ_powers(bi, lV-1, T, q);
     855         294 :   W = cgetg(lV+1, t_VEC);
     856        4816 :   for(i=1; i<=lV; i++)
     857        4522 :     gel(W, i) = FpXQX_FpXQ_mul(gel(V,i), gel(P, i), T, q);
     858         294 :   return gerepileupto(av, W);
     859             : }
     860             : 
     861             : struct _ZpXQ_sqrtn
     862             : {
     863             :   GEN T, a, n, ai;
     864             : };
     865             : 
     866             : static GEN
     867        1974 : _sqrtn_invd(void *E, GEN V, GEN v, GEN q, long M)
     868             : {
     869        1974 :   struct _ZpXQ_sqrtn *d = (struct _ZpXQ_sqrtn *) E;
     870        1974 :   GEN Tq = FpX_red(d->T, q), aiq = FpX_red(d->ai, q);
     871             :   (void)M;
     872        1974 :   return FpXQ_mul(FpXQ_mul(V, gel(v,2), Tq, q), aiq, Tq, q);
     873             : }
     874             : 
     875             : static GEN
     876        1974 : _sqrtn_eval(void *E, GEN x, GEN q)
     877             : {
     878        1974 :   struct _ZpXQ_sqrtn *d = (struct _ZpXQ_sqrtn *) E;
     879        1974 :   GEN Tq = FpX_red(d->T, q);
     880        1974 :   GEN f = FpX_sub(FpXQ_pow(x, d->n, Tq, q), d->a, q);
     881        1974 :   return mkvec2(f, x);
     882             : }
     883             : 
     884             : GEN
     885        1183 : ZpXQ_sqrtnlift(GEN a, GEN n, GEN x, GEN T, GEN p, long e)
     886             : {
     887             :   struct _ZpXQ_sqrtn d;
     888        1183 :   d.a = a; d.T = T; d.n = n;
     889        1183 :   d.ai = ZpXQ_inv(ZX_Z_mul(a, n),T,p,(e+1)>>1);
     890        1183 :   return gen_ZpX_Newton(x, p, e, &d, _sqrtn_eval, _sqrtn_invd);
     891             : }
     892             : GEN
     893         414 : Zq_sqrtnlift(GEN a, GEN n, GEN x, GEN T, GEN p, long e)
     894             : {
     895         456 :   return T? ZpXQ_sqrtnlift(typ(a) == t_INT? scalarpol_shallow(a, varn(T)):a,
     896             :                            n, x, T, p, e)
     897         456 :           : Zp_sqrtnlift(a, n, x, p, e);
     898             : }
     899             : 
     900             : GEN
     901           0 : ZpXQ_sqrt(GEN a, GEN T, GEN p, long e)
     902             : {
     903           0 :   pari_sp av = avma;
     904           0 :   GEN z = FpXQ_sqrt(FpX_red(a, p), T, p);
     905           0 :   if (!z) return NULL;
     906           0 :   if (e <= 1) return gerepileupto(av, z);
     907           0 :   return gerepileupto(av, ZpXQ_sqrtnlift(a, gen_2, z, T, p, e));
     908             : }
     909             : 
     910             : GEN
     911        3136 : ZpX_ZpXQ_liftroot_ea(GEN P, GEN S, GEN T, GEN p, long n, void *E,
     912             :                      int early(void *E, GEN x, GEN q))
     913             : {
     914        3136 :   pari_sp ltop = avma, av;
     915             :   long N, r;
     916             :   long mask;
     917             :   GEN q2, q, W, Q, Tq2, Tq, Pq;
     918             :   pari_timer ti;
     919        3136 :   T = FpX_get_red(T, powiu(p, n));
     920        3136 :   if (n == 1) return gcopy(S);
     921        3136 :   mask = quadratic_prec_mask(n);
     922        3136 :   av = avma;
     923        3136 :   q2 = p; q = sqri(p); mask >>= 1; N = 2;
     924        3136 :   if (DEBUGLEVEL > 3) timer_start(&ti);
     925        3136 :   Tq = FpXT_red(T,q);
     926        3136 :   Tq2 = FpXT_red(Tq,q2);
     927        3136 :   Pq = FpX_red(P,q);
     928        3136 :   W = FpXQ_inv(FpX_FpXQ_eval(FpX_deriv(P,q2), S, Tq2, q2), Tq2, q2);
     929        3136 :   Q  = ZX_Z_divexact(FpX_FpXQ_eval(Pq, S, Tq, q), q2);
     930        3136 :   r = brent_kung_optpow(degpol(P), 4, 3);
     931             :   for (;;)
     932             :   {
     933             :     GEN H, Sq, Wq, Spow, dP, qq, Pqq, Tqq;
     934        9279 :     H  = FpXQ_mul(W, Q, Tq2, q2);
     935        9279 :     Sq = FpX_sub(S, ZX_Z_mul(H, q2), q);
     936        9279 :     if (DEBUGLEVEL > 3)
     937           0 :       timer_printf(&ti,"ZpX_ZpXQ_liftroot: lift to prec %ld",N);
     938        9279 :     if (mask==1 || (early && early(E, Sq, q)))
     939        3136 :       return gerepileupto(ltop, Sq);
     940        6143 :     qq = sqri(q); N <<= 1;
     941        6143 :     if (mask&1UL) { qq = diviiexact(qq, p); N--; }
     942        6143 :     mask >>= 1;
     943        6143 :     Pqq  = FpX_red(P, qq);
     944        6143 :     Tqq  = FpXT_red(T, qq);
     945        6143 :     Spow = FpXQ_powers(Sq, r, Tqq, qq);
     946        6143 :     Q  = ZX_Z_divexact(FpX_FpXQV_eval(Pqq, Spow, Tqq, qq), q);
     947        6143 :     dP = FpX_FpXQV_eval(FpX_deriv(Pq, q), FpXV_red(Spow, q), Tq, q);
     948        6143 :     Wq = ZX_Z_divexact(FpX_Fp_sub(FpXQ_mul(W, dP, Tq, q), gen_1, q), q2);
     949        6143 :     Wq = ZX_Z_mul(FpXQ_mul(W, Wq, Tq2, q2), q2);
     950        6143 :     Wq = FpX_sub(W, Wq, q);
     951        6143 :     S = Sq; W = Wq; q2 = q; q = qq; Tq2 = Tq; Tq = Tqq; Pq = Pqq;
     952        6143 :     if (gc_needed(av, 1))
     953             :     {
     954           0 :       if(DEBUGMEM>1) pari_warn(warnmem,"ZpX_ZpXQ_Newton");
     955           0 :       gerepileall(av, 8, &S, &W, &Q, &Tq2, &Tq, &Pq, &q, &q2);
     956             :     }
     957        6143 :   }
     958             : }
     959             : 
     960             : GEN
     961        1582 : ZpX_ZpXQ_liftroot(GEN P, GEN S, GEN T, GEN p, long n)
     962             : {
     963        1582 :   return ZpX_ZpXQ_liftroot_ea(P, S, T, p, n, NULL, NULL);
     964             : }
     965             : 
     966             : GEN
     967         287 : ZpX_Frobenius(GEN T, GEN p, long e)
     968             : {
     969         287 :   return ZpX_ZpXQ_liftroot(get_FpX_mod(T), FpX_Frobenius(T, p), T, p, e);
     970             : }
     971             : 
     972             : GEN
     973         140 : ZpXQM_prodFrobenius(GEN M, GEN T, GEN p, long e)
     974             : {
     975         140 :   pari_sp av = avma;
     976         140 :   GEN xp = ZpX_Frobenius(T, p, e);
     977         140 :   GEN z = FpXQM_autsum(mkvec2(xp, M), get_FpX_degree(T), T, powiu(p,e));
     978         140 :   return gerepilecopy(av, gel(z,2));
     979             : }

Generated by: LCOV version 1.11