Code coverage tests

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

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

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

LCOV - code coverage report
Current view: top level - basemath - Hensel.c (source / functions) Hit Total Coverage
Test: PARI/GP v2.12.0 lcov report (development 23332-367b47754) Lines: 622 652 95.4 %
Date: 2018-12-10 05:41:52 Functions: 62 64 96.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       54772 : BuildTree(GEN link, GEN V, GEN W, GEN a, GEN T, GEN p)
      36             : {
      37       54772 :   long k = lg(a)-1;
      38             :   long i, j, s, minp, mind;
      39             : 
      40       54772 :   for (i=1; i<=k; i++) { gel(V,i) = gel(a,i); link[i] = -i; }
      41             : 
      42       88128 :   for (j=1; j <= 2*k-5; j+=2,i++)
      43             :   {
      44       33356 :     minp = j;
      45       33356 :     mind = degpol(gel(V,j));
      46      432059 :     for (s=j+1; s<i; s++)
      47      398703 :       if (degpol(gel(V,s)) < mind) { minp = s; mind = degpol(gel(V,s)); }
      48             : 
      49       33356 :     swap(gel(V,j), gel(V,minp));
      50       33356 :     lswap(link[j], link[minp]);
      51             : 
      52       33356 :     minp = j+1;
      53       33356 :     mind = degpol(gel(V,j+1));
      54      398703 :     for (s=j+2; s<i; s++)
      55      365347 :       if (degpol(gel(V,s)) < mind) { minp = s; mind = degpol(gel(V,s)); }
      56             : 
      57       33356 :     swap(gel(V,j+1), gel(V,minp));
      58       33356 :     lswap(link[j+1], link[minp]);
      59             : 
      60       33356 :     gel(V,i) = FqX_mul(gel(V,j), gel(V,j+1), T, p);
      61       33356 :     link[i] = j;
      62             :   }
      63             : 
      64      142893 :   for (j=1; j <= 2*k-3; j+=2)
      65             :   {
      66             :     GEN d, u, v;
      67       88128 :     d = FqX_extgcd(gel(V,j), gel(V,j+1), T, p, &u, &v);
      68       88128 :     if (degpol(d) > 0) pari_err_COPRIME("BuildTree", gel(V,j), gel(V,j+1));
      69       88121 :     d = gel(d,2);
      70       88121 :     if (!gequal1(d))
      71             :     {
      72       80492 :       if (typ(d)==t_POL)
      73             :       {
      74        5635 :         d = FpXQ_inv(d, T, p);
      75        5635 :         u = FqX_Fq_mul(u, d, T, p);
      76        5635 :         v = FqX_Fq_mul(v, d, T, p);
      77             :       }
      78             :       else
      79             :       {
      80       74857 :         d = Fp_inv(d, p);
      81       74857 :         u = FqX_Fp_mul(u, d, T,p);
      82       74857 :         v = FqX_Fp_mul(v, d, T,p);
      83             :       }
      84             :     }
      85       88121 :     gel(W,j) = u;
      86       88121 :     gel(W,j+1) = v;
      87             :   }
      88       54765 : }
      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      377379 : ZpX_HenselLift(GEN V, GEN W, long j, GEN f, GEN pd, GEN p0, GEN p1, int noinv)
      94             : {
      95      377379 :   pari_sp av = avma;
      96      377379 :   long space = lg(f) * lgefint(p1);
      97             :   GEN a2, b2, g, z, s, t;
      98      377379 :   GEN a = gel(V,j), b = gel(V,j+1);
      99      377379 :   GEN u = gel(W,j), v = gel(W,j+1);
     100             : 
     101      377379 :   (void)new_chunk(space); /* HACK */
     102      377379 :   g = ZX_sub(f, ZX_mul(a,b));
     103      377379 :   g = ZX_Z_divexact(g, p0);
     104      377379 :   g = FpX_red(g, pd);
     105      377379 :   z = FpX_mul(v,g, pd);
     106      377379 :   t = FpX_divrem(z,a, pd, &s);
     107      377379 :   t = ZX_add(ZX_mul(u,g), ZX_mul(t,b));
     108      377379 :   t = FpX_red(t, pd);
     109      377379 :   t = ZX_Z_mul(t,p0);
     110      377379 :   s = ZX_Z_mul(s,p0);
     111      377379 :   set_avma(av);
     112      377379 :   a2 = ZX_add(a,s);
     113      377379 :   b2 = ZX_add(b,t);
     114             : 
     115             :   /* already reduced mod p1 = pd p0 */
     116      377379 :   gel(V,j)   = a2;
     117      377379 :   gel(V,j+1) = b2;
     118      377379 :   if (noinv) return;
     119             : 
     120      297868 :   av = avma;
     121      297868 :   (void)new_chunk(space); /* HACK */
     122      297868 :   g = ZX_add(ZX_mul(u,a2), ZX_mul(v,b2));
     123      297868 :   g = Z_ZX_sub(gen_1, g);
     124      297868 :   g = ZX_Z_divexact(g, p0);
     125      297868 :   g = FpX_red(g, pd);
     126      297868 :   z = FpX_mul(v,g, pd);
     127      297868 :   t = FpX_divrem(z,a, pd, &s);
     128      297868 :   t = ZX_add(ZX_mul(u,g), ZX_mul(t,b));
     129      297868 :   t = FpX_red(t, pd);
     130      297868 :   t = ZX_Z_mul(t,p0);
     131      297868 :   s = ZX_Z_mul(s,p0);
     132      297868 :   set_avma(av);
     133      297868 :   gel(W,j)   = ZX_add(u,t);
     134      297868 :   gel(W,j+1) = ZX_add(v,s);
     135             : }
     136             : 
     137             : static void
     138       24892 : ZpXQ_HenselLift(GEN V, GEN W, long j, GEN f, GEN Td, GEN T1, GEN pd, GEN p0, GEN p1, int noinv)
     139             : {
     140       24892 :   pari_sp av = avma;
     141       24892 :   const long n = degpol(T1), vT = varn(T1);
     142       24892 :   long space = lg(f) * lgefint(p1) * lg(T1);
     143             :   GEN a2, b2, g, z, s, t;
     144       24892 :   GEN a = gel(V,j), b = gel(V,j+1);
     145       24892 :   GEN u = gel(W,j), v = gel(W,j+1);
     146             : 
     147       24892 :   (void)new_chunk(space); /* HACK */
     148       24892 :   g = RgX_sub(f, Kronecker_to_ZXX(ZXX_mul_Kronecker(a,b,n), n, vT));
     149       24892 :   g = FpXQX_red(g, T1, p1);
     150       24892 :   g = RgX_Rg_divexact(g, p0);
     151       24892 :   z = FpXQX_mul(v,g, Td,pd);
     152       24892 :   t = FpXQX_divrem(z,a, Td,pd, &s);
     153       24892 :   t = ZX_add(ZXX_mul_Kronecker(u,g,n), ZXX_mul_Kronecker(t,b,n));
     154       24892 :   t = Kronecker_to_ZXX(t, n, vT);
     155       24892 :   t = FpXQX_red(t, Td, pd);
     156       24892 :   t = RgX_Rg_mul(t,p0);
     157       24892 :   s = RgX_Rg_mul(s,p0);
     158       24892 :   set_avma(av);
     159             : 
     160       24892 :   a2 = RgX_add(a,s);
     161       24892 :   b2 = RgX_add(b,t);
     162             :   /* already reduced mod p1 = pd p0 */
     163       24892 :   gel(V,j)   = a2;
     164       24892 :   gel(V,j+1) = b2;
     165       24892 :   if (noinv) return;
     166             : 
     167       18802 :   av = avma;
     168       18802 :   (void)new_chunk(space); /* HACK */
     169       18802 :   g = ZX_add(ZXX_mul_Kronecker(u,a2,n), ZXX_mul_Kronecker(v,b2,n));
     170       18802 :   g = Kronecker_to_ZXX(g, n, vT);
     171       18802 :   g = Rg_RgX_sub(gen_1, g);
     172       18802 :   g = FpXQX_red(g, T1, p1);
     173       18802 :   g = RgX_Rg_divexact(g, p0);
     174       18802 :   z = FpXQX_mul(v,g, Td,pd);
     175       18802 :   t = FpXQX_divrem(z,a, Td,pd, &s);
     176       18802 :   t = ZX_add(ZXX_mul_Kronecker(u,g,n), ZXX_mul_Kronecker(t,b,n));
     177       18802 :   t = Kronecker_to_ZXX(t, n, vT);
     178       18802 :   t = FpXQX_red(t, Td, pd);
     179       18802 :   t = RgX_Rg_mul(t,p0);
     180       18802 :   s = RgX_Rg_mul(s,p0);
     181       18802 :   set_avma(av);
     182       18802 :   gel(W,j)   = RgX_add(u,t);
     183       18802 :   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      970321 : ZpX_RecTreeLift(GEN link, GEN v, GEN w, GEN pd, GEN p0, GEN p1,
     190             :                 GEN f, long j, int noinv)
     191             : {
     192      970321 :   if (j < 0) return;
     193      377379 :   ZpX_HenselLift(v, w, j, f, pd, p0,p1, noinv);
     194      377379 :   ZpX_RecTreeLift(link, v, w, pd, p0,p1, gel(v,j)  , link[j  ], noinv);
     195      377379 :   ZpX_RecTreeLift(link, v, w, pd, p0,p1, gel(v,j+1), link[j+1], noinv);
     196             : }
     197             : static void
     198       53200 : 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       53200 :   if (j < 0) return;
     202       24892 :   ZpXQ_HenselLift(v, w, j, f, Td,T1, pd, p0,p1, noinv);
     203       24892 :   ZpXQ_RecTreeLift(link, v, w, Td,T1, pd, p0,p1, gel(v,j)  , link[j  ], noinv);
     204       24892 :   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      995590 : quadratic_prec_mask(long n)
     225             : {
     226      995590 :   long a = n, i;
     227      995590 :   ulong mask = 0;
     228     3841700 :   for(i = 1;; i++, mask <<= 1)
     229             :   {
     230     6687810 :     mask |= (a&1); a = (a+1)>>1;
     231     4837290 :     if (a==1) return mask | (1UL << i);
     232             :   }
     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       54800 : MultiLift(GEN f, GEN a, GEN T, GEN p, long e0, long flag)
     242             : {
     243       54800 :   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       54800 :   if (k < 2) pari_err_DOMAIN("MultiLift", "#(modular factors)", "<", gen_2,a);
     249       54800 :   if (e0 < 1) pari_err_DOMAIN("MultiLift", "precision", "<", gen_1,stoi(e0));
     250       54800 :   if (e0 == 1) return a;
     251             : 
     252       54772 :   if (DEBUGLEVEL > 3) timer_start(&Ti);
     253       54772 :   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       54772 :     e = 1;
     263       54772 :     v = cgetg(2*k-2 + 1, t_VEC);
     264       54772 :     w = cgetg(2*k-2 + 1, t_VEC);
     265       54772 :     link=cgetg(2*k-2 + 1, t_VECSMALL);
     266       54772 :     BuildTree(link, v, w, a, T? FpX_red(T,p): NULL, p);
     267       54765 :     if (DEBUGLEVEL > 3) timer_printf(&Ti, "building tree");
     268             :   }
     269       54765 :   mask = quadratic_prec_mask(e0);
     270       54765 :   eold = 1;
     271       54765 :   penew = NULL;
     272       54765 :   Tnew = NULL;
     273       54765 :   if (DEBUGLEVEL > 3) err_printf("lifting to prec %ld\n", e0);
     274      328509 :   while (mask > 1)
     275             :   {
     276      218979 :     long enew = eold << 1;
     277      218979 :     if (mask & 1) enew--;
     278      218979 :     mask >>= 1;
     279      218979 :     if (enew >= e) { /* mask == 1: last iteration */
     280      218979 :       GEN peold = penew? penew: powiu(p, eold);
     281      218979 :       GEN Td = NULL, pd;
     282      218979 :       long d = enew - eold; /* = eold or eold-1 */
     283             :       /* lift from p^eold to p^enew */
     284      218979 :       pd = (d == eold)? peold: diviiexact(peold, p); /* p^d */
     285      218979 :       penew = mulii(peold,pd);
     286      218979 :       if (T) {
     287        3416 :         if (Tnew)
     288        2450 :           Td = (d == eold)? Tnew: FpX_red(Tnew,pd);
     289             :         else
     290         966 :           Td = FpX_red(T, peold);
     291        3416 :         Tnew = FpX_red(T, penew);
     292        3416 :         ZpXQ_RecTreeLift(link, v, w, Td, Tnew, pd, peold, penew, f, lg(v)-2,
     293             :                          (flag == 0 && mask == 1));
     294             :       }
     295             :       else
     296      215563 :         ZpX_RecTreeLift(link, v, w, pd, peold, penew, f, lg(v)-2,
     297             :                         (flag == 0 && mask == 1));
     298      218979 :       if (DEBUGLEVEL > 3) timer_printf(&Ti, "reaching prec %ld", enew);
     299             :     }
     300      218979 :     eold = enew;
     301             :   }
     302             : 
     303       54765 :   if (flag)
     304         805 :     E = mkvec4(utoipos(e0), link, v, w);
     305             :   else
     306             :   {
     307       53960 :     E = cgetg(k+1, t_VEC);
     308      225162 :     for (i = 1; i <= 2*k-2; i++)
     309             :     {
     310      171202 :       long t = link[i];
     311      171202 :       if (t < 0) gel(E,-t) = gel(v,i);
     312             :     }
     313             :   }
     314       54765 :   return E;
     315             : }
     316             : 
     317             : /* Q list of (coprime, monic) factors of pol mod (T,p). Lift mod p^e = pe.
     318             :  * T may be NULL */
     319             : GEN
     320       93279 : ZpX_liftfact(GEN pol, GEN Q, GEN pe, GEN p, long e)
     321             : {
     322       93279 :   pari_sp av = avma;
     323       93279 :   pol = FpX_normalize(pol, pe);
     324       93279 :   if (lg(Q) == 2) return mkvec(pol);
     325       53022 :   return gerepilecopy(av, MultiLift(pol, Q, NULL, p, e, 0));
     326             : }
     327             : 
     328             : GEN
     329         973 : ZpXQX_liftfact(GEN pol, GEN Q, GEN T, GEN pe, GEN p, long e)
     330             : {
     331         973 :   pari_sp av = avma;
     332         973 :   pol = FpXQX_normalize(pol, T, pe);
     333         973 :   if (lg(Q) == 2) return mkvec(pol);
     334         973 :   return gerepilecopy(av, MultiLift(pol, Q, T, p, e, 0));
     335             : }
     336             : 
     337             : GEN
     338        1946 : ZqX_liftfact(GEN f, GEN a, GEN T, GEN pe, GEN p, long e)
     339        1946 : { return T ? ZpXQX_liftfact(f, a, T, pe, p, e): ZpX_liftfact(f, a, pe, p, e); }
     340             : GEN
     341          42 : ZqX_liftroot(GEN f, GEN a, GEN T, GEN p, long e)
     342          42 : { return T ? ZpXQX_liftroot(f, a,T , p, e): ZpX_liftroot(f, a, p, e); }
     343             : 
     344             : 
     345             : /* U = NULL treated as 1 */
     346             : static void
     347        5845 : BezoutPropagate(GEN link, GEN v, GEN w, GEN pe, GEN U, GEN f, long j)
     348             : {
     349             :   GEN Q, R;
     350        5845 :   if (j < 0) return;
     351             : 
     352        2520 :   Q = FpX_mul(gel(v,j), gel(w,j), pe);
     353        2520 :   if (U)
     354             :   {
     355        1715 :     Q = FpXQ_mul(Q, U, f, pe);
     356        1715 :     R = FpX_sub(U, Q, pe);
     357             :   }
     358             :   else
     359         805 :     R = Fp_FpX_sub(gen_1, Q, pe);
     360        2520 :   gel(w,j+1) = Q; /*  0 mod U v[j],  1 mod (1-U) v[j+1] */
     361        2520 :   gel(w,j) = R; /*  1 mod U v[j],  0 mod (1-U) v[j+1] */
     362        2520 :   BezoutPropagate(link, v, w, pe, R, f, link[j  ]);
     363        2520 :   BezoutPropagate(link, v, w, pe, Q, f, link[j+1]);
     364             : }
     365             : 
     366             : /* as above, but return the Bezout coefficients for the lifted modular factors
     367             :  *   U[i] = 1 mod Qlift[i]
     368             :  *          0 mod Qlift[j], j != i */
     369             : GEN
     370         812 : bezout_lift_fact(GEN pol, GEN Q, GEN p, long e)
     371             : {
     372         812 :   pari_sp av = avma;
     373             :   GEN E, link, v, w, pe;
     374         812 :   long i, k = lg(Q)-1;
     375         812 :   if (k == 1) retmkvec(pol_1(varn(pol)));
     376         805 :   pe = powiu(p, e);
     377         805 :   pol = FpX_normalize(pol, pe);
     378         805 :   E = MultiLift(pol, Q, NULL, p, e, 1);
     379         805 :   link = gel(E,2);
     380         805 :   v    = gel(E,3);
     381         805 :   w    = gel(E,4);
     382         805 :   BezoutPropagate(link, v, w, pe, NULL, pol, lg(v)-2);
     383         805 :   E = cgetg(k+1, t_VEC);
     384        5845 :   for (i = 1; i <= 2*k-2; i++)
     385             :   {
     386        5040 :     long t = link[i];
     387        5040 :     if (t < 0) E[-t] = w[i];
     388             :   }
     389         805 :   return gerepilecopy(av, E);
     390             : }
     391             : 
     392             : /* Front-end for ZpX_liftfact:
     393             :    lift the factorization of pol mod p given by L to p^N (if possible) */
     394             : GEN
     395          28 : polhensellift(GEN pol, GEN L, GEN p, long N)
     396             : {
     397          28 :   GEN T = NULL;
     398             :   long i, l, t;
     399          28 :   pari_sp av = avma;
     400             :   void (*chk)(GEN, const char*);
     401             : 
     402          28 :   if (typ(pol) != t_POL) pari_err_TYPE("polhensellift",pol);
     403          28 :   RgX_check_ZXX(pol, "polhensellift");
     404          28 :   if (!is_vec_t(typ(L)) || lg(L) < 3) pari_err_TYPE("polhensellift",L);
     405          28 :   t = typ(p);
     406          28 :   if (t == t_VEC) /* [p, T] */
     407             :   {
     408           7 :     T = gel(p,2);
     409           7 :     if (typ(T) != t_POL) pari_err_TYPE("polhensellift",pol);
     410           7 :     RgX_check_ZX(T, "polhensellift");
     411           7 :     p = gel(p,1); t = typ(p);
     412             :   }
     413          28 :   chk = T? RgX_check_ZXX: RgX_check_ZX;
     414          28 :   if (t != t_INT) pari_err_TYPE("polhensellift",p);
     415          28 :   if (N < 1) pari_err_DOMAIN("polhensellift", "precision", "<", gen_1,stoi(N));
     416             : 
     417          28 :   l = lg(L); L = leafcopy(L);
     418          70 :   for (i = 1; i < l; i++)
     419             :   {
     420          49 :     GEN q = gel(L,i);
     421          49 :     if (typ(q) != t_POL) gel(L,i) = scalar_ZX_shallow(q, varn(pol));
     422          49 :     else chk(q, "polhensellift");
     423             :   }
     424          21 :   return gerepilecopy(av, ZqX_liftfact(pol, L, T, powiu(p,N), p, N));
     425             : }
     426             : 
     427             : static GEN
     428       43878 : FqV_roots_from_deg1(GEN x, GEN T, GEN p)
     429             : {
     430       43878 :   long i,l = lg(x);
     431       43878 :   GEN r = cgetg(l,t_VEC);
     432       43878 :   for (i=1; i<l; i++) { GEN P = gel(x,i); gel(r,i) = Fq_neg(gel(P,2), T, p); }
     433       43878 :   return r;
     434             : }
     435             : 
     436             : static GEN
     437       43843 : ZpX_liftroots_full(GEN f, GEN S, GEN q, GEN p, long e)
     438             : {
     439       43843 :   pari_sp av = avma;
     440       43843 :   GEN y = ZpX_liftfact(f, deg1_from_roots(S, varn(f)), q, p, e);
     441       43843 :   return gerepileupto(av, FqV_roots_from_deg1(y, NULL, q));
     442             : }
     443             : 
     444             : GEN
     445       43815 : ZpX_roots(GEN F, GEN p, long e)
     446             : {
     447       43815 :   pari_sp av = avma;
     448       43815 :   GEN q = powiu(p, e);
     449       43815 :   GEN f = FpX_normalize(F, p);
     450       43815 :   GEN g = FpX_normalize(FpX_split_part(f, p), p);
     451             :   GEN S;
     452       43815 :   if (degpol(g) < degpol(f))
     453             :   {
     454       39999 :     GEN h = FpX_div(f, g, p);
     455       39999 :     F = gel(ZpX_liftfact(F, mkvec2(g, h), q, p, e), 1);
     456             :   }
     457       43815 :   S = FpX_roots(g, p);
     458       43815 :   return gerepileupto(av, ZpX_liftroots_full(F, S, q, p, e));
     459             : }
     460             : 
     461             : static GEN
     462          35 : ZpXQX_liftroots_full(GEN f, GEN S, GEN T, GEN q, GEN p, long e)
     463             : {
     464          35 :   pari_sp av = avma;
     465          35 :   GEN y = ZpXQX_liftfact(f, deg1_from_roots(S, varn(f)), T, q, p, e);
     466          35 :   return gerepileupto(av, FqV_roots_from_deg1(y, T, q));
     467             : }
     468             : 
     469             : GEN
     470          35 : ZpXQX_roots(GEN F, GEN T, GEN p, long e)
     471             : {
     472          35 :   pari_sp av = avma;
     473          35 :   GEN q = powiu(p, e);
     474          35 :   GEN f = FpXQX_normalize(F, T, p);
     475          35 :   GEN g = FpXQX_normalize(FpXQX_split_part(f, T, p), T, p);
     476             :   GEN S;
     477          35 :   if (degpol(g) < degpol(f))
     478             :   {
     479           7 :     GEN h = FpXQX_div(f, g, T, p);
     480           7 :     F = gel(ZpXQX_liftfact(F, mkvec2(g, h), T, q, p, e), 1);
     481             :   }
     482          35 :   S = FpXQX_roots(g, T, p);
     483          35 :   return gerepileupto(av, ZpXQX_liftroots_full(F, S, T, q, p, e));
     484             : }
     485             : 
     486             : GEN
     487        2353 : ZqX_roots(GEN F, GEN T, GEN p, long e)
     488             : {
     489        2353 :   return T ? ZpXQX_roots(F, T, p, e): ZpX_roots(F, p, e);
     490             : }
     491             : 
     492             : /* SPEC:
     493             :  * p is a t_INT > 1, e >= 1
     494             :  * f is a ZX with leading term prime to p.
     495             :  * a is a simple root mod l for all l|p.
     496             :  * Return roots of f mod p^e, as integers (implicitly mod p^e)
     497             :  * STANDARD USE: p is a prime power */
     498             : GEN
     499         406 : ZpX_liftroot(GEN f, GEN a, GEN p, long e)
     500             : {
     501         406 :   pari_sp av = avma;
     502         406 :   GEN q = p, fr, W;
     503             :   ulong mask;
     504             : 
     505         406 :   a = modii(a,q);
     506         406 :   if (e == 1) return a;
     507         406 :   mask = quadratic_prec_mask(e);
     508         406 :   fr = FpX_red(f,q);
     509         406 :   W = Fp_inv(FpX_eval(ZX_deriv(fr), a, q), q); /* 1/f'(a) mod p */
     510             :   for(;;)
     511             :   {
     512        2632 :     q = sqri(q);
     513        1519 :     if (mask & 1) q = diviiexact(q, p);
     514        1519 :     mask >>= 1;
     515        1519 :     fr = FpX_red(f,q);
     516        1519 :     a = Fp_sub(a, Fp_mul(W, FpX_eval(fr, a,q), q), q);
     517        1519 :     if (mask == 1) return gerepileuptoint(av, a);
     518        1113 :     W = Fp_sub(shifti(W,1), Fp_mul(Fp_sqr(W,q), FpX_eval(ZX_deriv(fr),a,q), q), q);
     519             :   }
     520             : }
     521             : 
     522             : GEN
     523          28 : ZpX_liftroots(GEN f, GEN S, GEN p, long e)
     524             : {
     525          28 :   long i, n = lg(S)-1, d = degpol(f);
     526             :   GEN r;
     527          28 :   if (n == d) return ZpX_liftroots_full(f, S, powiu(p, e), p, e);
     528           0 :   r = cgetg(n+1, typ(S));
     529           0 :   for (i=1; i <= n; i++)
     530           0 :     gel(r,i) = ZpX_liftroot(f, gel(S,i), p, e);
     531           0 :   return r;
     532             : }
     533             : 
     534             : GEN
     535          70 : ZpXQX_liftroot_vald(GEN f, GEN a, long v, GEN T, GEN p, long e)
     536             : {
     537          70 :   pari_sp av = avma, av2;
     538          70 :   GEN pv = p, q, W, df, Tq, fr, dfr;
     539             :   ulong mask;
     540          70 :   a = Fq_red(a, T, p);
     541          70 :   if (e <= v+1) return a;
     542          70 :   df = RgX_deriv(f);
     543          70 :   if (v) { pv = powiu(p,v); df = ZXX_Z_divexact(df, pv); }
     544          70 :   mask = quadratic_prec_mask(e-v);
     545          70 :   Tq = FpXT_red(T, p); dfr = FpXQX_red(df, Tq, p);
     546          70 :   W = Fq_inv(FqX_eval(dfr, a, Tq, p), Tq, p); /* 1/f'(a) mod (T,p) */
     547          70 :   q = p;
     548          70 :   av2 = avma;
     549             :   for (;;)
     550          91 :   {
     551             :     GEN u, fa, qv, q2v, q2, Tq2;
     552         161 :     q2 = q; q = sqri(q);
     553         161 :     if (mask & 1) q = diviiexact(q,p);
     554         161 :     mask >>= 1;
     555         161 :     if (v) { qv = mulii(q, pv); q2v = mulii(q2, pv); }
     556         161 :     else { qv = q; q2v = q2; }
     557         161 :     Tq2 = FpXT_red(T, q2v); Tq = FpXT_red(T, qv);
     558         161 :     fr = FpXQX_red(f, Tq, qv);
     559         161 :     fa = FqX_eval(fr, a, Tq, qv);
     560         161 :     fa = typ(fa)==t_INT? diviiexact(fa,q2v): ZX_Z_divexact(fa, q2v);
     561         161 :     a = Fq_sub(a, Fq_mul(Fq_mul(W,fa,Tq2,q2v),q2, Tq,qv), Tq, qv);
     562         161 :     if (mask == 1) return gerepileupto(av, a);
     563          91 :     dfr = FpXQX_red(df, Tq, q);
     564          91 :     u = Fq_sub(Fq_mul(W,FqX_eval(dfr,a,Tq,q),Tq,q),gen_1,Tq,q);
     565          91 :     u = typ(u)==t_INT? diviiexact(u,q2): ZX_Z_divexact(u,q2);
     566          91 :     W = Fq_sub(W, Fq_mul(Fq_mul(u,W,Tq2,q2),q2, Tq,q),Tq,q);
     567          91 :     if (gc_needed(av2,2))
     568             :     {
     569           0 :       if(DEBUGMEM>1) pari_warn(warnmem,"ZpXQX_liftroot, e = %ld", e);
     570           0 :       gerepileall(av2, 3, &a, &W, &q);
     571             :     }
     572             :   }
     573             : }
     574             : 
     575             : GEN
     576          70 : ZpXQX_liftroot(GEN f, GEN a, GEN T, GEN p, long e) { return ZpXQX_liftroot_vald(f,a,0,T,p,e); }
     577             : 
     578             : GEN
     579           0 : ZpXQX_liftroots(GEN f, GEN S, GEN T, GEN p, long e)
     580             : {
     581           0 :   long i, n = lg(S)-1, d = degpol(f);
     582             :   GEN r;
     583           0 :   if (n == d) return ZpXQX_liftroots_full(f, S, T, powiu(p, e), p, e);
     584           0 :   r = cgetg(n+1, typ(S));
     585           0 :   for (i=1; i <= n; i++)
     586           0 :     gel(r,i) = ZpXQX_liftroot(f, gel(S,i), T, p, e);
     587           0 :   return r;
     588             : }
     589             : 
     590             : /* Same as ZpX_liftroot for the polynomial X^n-b*/
     591             : GEN
     592      176556 : Zp_sqrtnlift(GEN b, GEN n, GEN a, GEN p, long e)
     593             : {
     594      176556 :   pari_sp ltop=avma;
     595             :   GEN q, w, n_1;
     596             :   ulong mask;
     597      176556 :   long pis2 = equalii(n, gen_2)? 1: 0;
     598      176556 :   if (e == 1) return icopy(a);
     599      171649 :   n_1 = subiu(n,1);
     600      171649 :   mask = quadratic_prec_mask(e);
     601      171649 :   w = Fp_inv(pis2 ? shifti(a,1): Fp_mul(n,Fp_pow(a,n_1,p), p), p);
     602      171649 :   q = p;
     603             :   for(;;)
     604             :   {
     605      363671 :     q = sqri(q);
     606      267660 :     if (mask & 1) q = diviiexact(q, p);
     607      267660 :     mask >>= 1;
     608      267660 :     if (lgefint(q) == 3 && lgefint(n) == 3)
     609       95859 :     {
     610      267142 :       ulong Q = uel(q,2), N = uel(n,2);
     611      267142 :       ulong A = umodiu(a, Q);
     612      267142 :       ulong B = umodiu(b, Q);
     613      267142 :       ulong W = umodiu(w, Q);
     614             : 
     615      267142 :       A = Fl_sub(A, Fl_mul(W, Fl_sub(Fl_powu(A,N,Q), B, Q), Q), Q);
     616      267142 :       a = utoi(A);
     617      267142 :       if (mask == 1) break;
     618       95859 :       W = Fl_sub(Fl_add(W,W,Q),
     619             :                  Fl_mul(Fl_sqr(W,Q), Fl_mul(N,Fl_powu(A, N-1, Q), Q), Q), Q);
     620       95859 :       w = utoi(W);
     621             :     }
     622             :     else
     623             :     {
     624             :       /* a -= w (a^n - b) */
     625         518 :       a = modii(subii(a, mulii(w, subii(Fp_pow(a,n,q),b))), q);
     626         518 :       if (mask == 1) break;
     627             :       /* w += w - w^2 n a^(n-1)*/
     628         289 :       w = subii(shifti(w,1), Fp_mul(Fp_sqr(w,q),
     629         137 :                            pis2? shifti(a,1): mulii(n,Fp_pow(a,n_1,q)), q));
     630             :     }
     631             :   }
     632      171649 :   return gerepileuptoint(ltop,a);
     633             : }
     634             : 
     635             : 
     636             : /* Same as ZpX_liftroot for the polynomial X^2-b */
     637             : GEN
     638      139363 : Zp_sqrtlift(GEN b, GEN a, GEN p, long e)
     639             : {
     640      139363 :   return Zp_sqrtnlift(b, gen_2, a, p, e);
     641             : }
     642             : 
     643             : GEN
     644     1354794 : Zp_sqrt(GEN x, GEN p, long e)
     645             : {
     646             :   pari_sp av;
     647             :   GEN z;
     648     1354794 :   if (absequaliu(p,2)) return Z2_sqrt(x,e);
     649     1353317 :   av = avma;
     650     1353317 :   z = Fp_sqrt(Fp_red(x, p), p);
     651     1353317 :   if (!z) return NULL;
     652     1353282 :   if (e > 1) z = Zp_sqrtlift(x, z, p, e);
     653     1353282 :   return gerepileuptoint(av, z);
     654             : }
     655             : 
     656             : /* Compute (x-1)/(x+1)/p^k */
     657             : static GEN
     658       20875 : ZpXQ_log_to_ath(GEN x, long k, GEN T, GEN p, long e, GEN pe)
     659             : {
     660       20875 :   pari_sp av = avma;
     661       20875 :   long vT = get_FpX_var(T);
     662             :   GEN bn, bdi;
     663       20875 :   GEN bd = ZX_Z_add(x, gen_1);
     664       20875 :   if (absequaliu(p,2)) /*For p=2, we need to simplify by 2*/
     665             :   {
     666        6945 :     bn = ZX_shifti(x,-(k+1));
     667        6945 :     bdi= ZpXQ_invlift(ZX_shifti(bd ,-1), pol_1(vT), T, p, e);
     668             :   }
     669             :   else
     670             :   {
     671       13930 :     bn = ZX_Z_divexact(ZX_Z_sub(x, gen_1),powiu(p,k));
     672       13930 :     bdi= ZpXQ_invlift(bd, scalarpol(Fp_inv(gen_2,p),vT), T, p, e);
     673             :   }
     674       20875 :   return gerepileupto(av, FpXQ_mul(bn, bdi, T, pe));
     675             : }
     676             : 
     677             : /* Assume p odd, a = 1 [p], return log(a) */
     678             : GEN
     679       20875 : ZpXQ_log(GEN a, GEN T, GEN p, long N)
     680             : {
     681       20875 :   pari_sp av = avma;
     682             :   pari_timer ti;
     683       20875 :   long is2 = absequaliu(p,2);
     684       20875 :   ulong pp = is2 ? 0: itou_or_0(p);
     685       20875 :   double lp = is2 ? 1: pp ? log2(pp): expi(p);
     686       20875 :   long k = maxss(1 , (long) .5+pow((double)(N>>1)/(lp*lp), 1./3));
     687             :   GEN ak, s, b, pol;
     688       20875 :   long e = is2 ? N-1: N;
     689       20875 :   long i, l = (e-2)/(2*(k+is2));
     690       20875 :   GEN pe = powiu(p,e);
     691       20875 :   GEN TNk, pNk = powiu(p,N+k);
     692       20875 :   if( DEBUGLEVEL>=3) timer_start(&ti);
     693       20875 :   TNk = FpX_get_red(get_FpX_mod(T), pNk);
     694       20875 :   ak = FpXQ_pow(a, powiu(p,k), TNk, pNk);
     695       20875 :   if( DEBUGLEVEL>=3) timer_printf(&ti,"FpXQ_pow(%ld)",k);
     696       20875 :   b = ZpXQ_log_to_ath(ak, k, T, p, e, pe);
     697       20875 :   if( DEBUGLEVEL>=3) timer_printf(&ti,"ZpXQ_log_to_ath");
     698       20875 :   pol= cgetg(l+3,t_POL);
     699       20875 :   pol[1] = evalsigne(1)|evalvarn(0);
     700       56330 :   for(i=0; i<=l; i++)
     701             :   {
     702             :     GEN g;
     703       35455 :     ulong z = 2*i+1;
     704       35455 :     if (pp)
     705             :     {
     706       26222 :       long w = u_lvalrem(z, pp, &z);
     707       26222 :       g = powuu(pp,2*i*k-w);
     708             :     }
     709        9233 :     else g = powiu(p,2*i*k);
     710       35452 :     gel(pol,i+2) = Fp_divu(g, z,pe);
     711             :   }
     712       20875 :   if( DEBUGLEVEL>=3) timer_printf(&ti,"pol(%ld)",l);
     713       20875 :   s = FpX_FpXQ_eval(pol, FpXQ_sqr(b, T, pe), T,  pe);
     714       20875 :   if( DEBUGLEVEL>=3) timer_printf(&ti,"FpX_FpXQ_eval");
     715       20875 :   s = ZX_shifti(FpXQ_mul(b, s, T, pe), 1);
     716       20875 :   return gerepileupto(av, is2? s: FpX_red(s, pe));
     717             : }
     718             : 
     719             : /***********************************************************************/
     720             : /**                                                                   **/
     721             : /**                 Generic quadratic hensel lift over Zp[X]          **/
     722             : /**                                                                   **/
     723             : /***********************************************************************/
     724             : /* q = p^N */
     725             : 
     726             : GEN
     727      199829 : gen_ZpX_Dixon(GEN F, GEN V, GEN q, GEN p, long N, void *E,
     728             :                             GEN lin(void *E, GEN F, GEN d, GEN q),
     729             :                             GEN invl(void *E, GEN d))
     730             : {
     731      199829 :   pari_sp av = avma;
     732             :   long N2, M;
     733             :   GEN VN2, V2, VM, bil;
     734             :   GEN q2, qM;
     735      199829 :   V = FpX_red(V, q);
     736      199829 :   if (N == 1) return invl(E, V);
     737       45458 :   N2 = (N + 1)>>1; M = N - N2;
     738       45458 :   F = FpXT_red(F, q);
     739       45458 :   qM = powiu(p, M);
     740       45458 :   q2 = M == N2? qM: mulii(qM, p);
     741             :   /* q2 = p^N2, qM = p^M, q = q2 * qM */
     742       45458 :   VN2 = gen_ZpX_Dixon(F, V, q2, p, N2, E, lin, invl);
     743       45458 :   bil = lin(E, F, VN2, q);
     744       45458 :   V2 = ZX_Z_divexact(ZX_sub(V, bil), q2);
     745       45458 :   VM = gen_ZpX_Dixon(F, V2, qM, p, M, E, lin, invl);
     746       45458 :   return gerepileupto(av, FpX_red(ZX_add(VN2, ZX_Z_mul(VM, q2)), q));
     747             : }
     748             : 
     749             : GEN
     750      261782 : gen_ZpX_Newton(GEN x, GEN p, long n, void *E,
     751             :                       GEN eval(void *E, GEN f, GEN q),
     752             :                       GEN invd(void *E, GEN V, GEN v, GEN q, long M))
     753             : {
     754      261782 :   pari_sp ltop = avma, av;
     755      261782 :   long N = 1, N2, M;
     756             :   long mask;
     757      261782 :   GEN q = p;
     758      261782 :   if (n == 1) return gcopy(x);
     759      257974 :   mask = quadratic_prec_mask(n);
     760      257974 :   av = avma;
     761     1164368 :   while (mask > 1)
     762             :   {
     763             :     GEN qM, q2, v, V;
     764      648420 :     N2 = N; N <<= 1;
     765      648420 :     q2 = q;
     766      648420 :     if (mask&1UL) { /* can never happen when q2 = p */
     767      271530 :       N--; M = N2-1;
     768      271530 :       qM = diviiexact(q2,p); /* > 1 */
     769      271529 :       q = mulii(qM,q2);
     770             :     } else {
     771      376890 :       M = N2;
     772      376890 :       qM = q2;
     773      376890 :       q = sqri(q2);
     774             :     }
     775             :     /* q2 = p^N2, qM = p^M, q = p^N = q2 * qM */
     776      648415 :     mask >>= 1;
     777      648415 :     v = eval(E, x, q);
     778      648421 :     V = ZX_Z_divexact(gel(v,1), q2);
     779      648408 :     x = FpX_sub(x, ZX_Z_mul(invd(E, V, v, qM, M), q2), q);
     780      648418 :     if (gc_needed(av, 1))
     781             :     {
     782           0 :       if(DEBUGMEM>1) pari_warn(warnmem,"gen_ZpX_Newton");
     783           0 :       gerepileall(av, 2, &x, &q);
     784             :     }
     785             :   }
     786      257974 :   return gerepileupto(ltop, x);
     787             : }
     788             : 
     789             : struct _ZpXQ_inv
     790             : {
     791             :   GEN T, a, p ,n;
     792             : };
     793             : 
     794             : static GEN
     795      506567 : _inv_invd(void *E, GEN V, GEN v, GEN q, long M/*unused*/)
     796             : {
     797      506567 :   struct _ZpXQ_inv *d = (struct _ZpXQ_inv *) E;
     798      506567 :   GEN Tq = FpXT_red(d->T, q);
     799             :   (void)M;
     800      506568 :   return FpXQ_mul(V, gel(v,2), Tq, q);
     801             : }
     802             : 
     803             : static GEN
     804      506568 : _inv_eval(void *E, GEN x, GEN q)
     805             : {
     806      506568 :   struct _ZpXQ_inv *d = (struct _ZpXQ_inv *) E;
     807      506568 :   GEN Tq = FpXT_red(d->T, q);
     808      506568 :   GEN f = FpX_Fp_sub(FpXQ_mul(x, FpX_red(d->a, q), Tq, q), gen_1, q);
     809      506569 :   return mkvec2(f, x);
     810             : }
     811             : 
     812             : GEN
     813      197128 : ZpXQ_invlift(GEN a, GEN x, GEN T, GEN p, long e)
     814             : {
     815             :   struct _ZpXQ_inv d;
     816      197128 :   d.a = a; d.T = T; d.p = p;
     817      197128 :   return gen_ZpX_Newton(x, p, e, &d, _inv_eval, _inv_invd);
     818             : }
     819             : 
     820             : GEN
     821      176253 : ZpXQ_inv(GEN a, GEN T, GEN p, long e)
     822             : {
     823      176253 :   pari_sp av=avma;
     824             :   GEN ai;
     825      176253 :   if (lgefint(p)==3)
     826             :   {
     827      176225 :     ulong pp = p[2];
     828      176225 :     ai = Flx_to_ZX(Flxq_inv(ZX_to_Flx(a,pp), ZXT_to_FlxT(T, pp), pp));
     829             :   } else
     830          28 :     ai = FpXQ_inv(FpX_red(a,p), FpXT_red(T,p),p);
     831      176253 :   return gerepileupto(av, ZpXQ_invlift(a, ai, T, p, e));
     832             : }
     833             : 
     834             : GEN
     835       31584 : ZpXQ_div(GEN a, GEN b, GEN T, GEN q, GEN p, long e)
     836             : {
     837       31584 :   return FpXQ_mul(a, ZpXQ_inv(b, T, p, e), T, q);
     838             : }
     839             : 
     840             : GEN
     841      142422 : ZpXQX_divrem(GEN x, GEN Sp, GEN T, GEN q, GEN p, long e, GEN *pr)
     842             : {
     843      142422 :   pari_sp av = avma;
     844      142422 :   GEN S = get_FpXQX_mod(Sp);
     845      142422 :   GEN b = leading_coeff(S), bi;
     846             :   GEN S2, Q;
     847      142422 :   if (typ(b)==t_INT) return FpXQX_divrem(x, Sp, T, q, pr);
     848      142422 :   bi = ZpXQ_inv(b, T, p, e);
     849      142422 :   S2 = FqX_Fq_mul_to_monic(S, bi, T, q);
     850      142422 :   Q = FpXQX_divrem(x, S2, T, q, pr);
     851      142422 :   if (pr==ONLY_DIVIDES && !Q) { set_avma(av); return NULL; }
     852      142422 :   if (pr==ONLY_REM || pr==ONLY_DIVIDES) return gerepileupto(av, Q);
     853      142422 :   Q = FpXQX_FpXQ_mul(Q, bi, T, q);
     854      142422 :   gerepileall(av, 2, &Q, pr);
     855      142422 :   return Q;
     856             : }
     857             : 
     858             : GEN
     859         588 : ZpXQX_digits(GEN x, GEN B, GEN T, GEN q, GEN p, long e)
     860             : {
     861         588 :   pari_sp av = avma;
     862         588 :   GEN b = leading_coeff(B), bi;
     863             :   GEN B2, P, V, W;
     864             :   long i, lV;
     865         588 :   if (typ(b)==t_INT) return FpXQX_digits(x, B, T, q);
     866         588 :   bi = ZpXQ_inv(b, T, p, e);
     867         588 :   B2 = FqX_Fq_mul_to_monic(B, bi, T, q);
     868         588 :   V = FpXQX_digits(x, B2, T, q);
     869         588 :   lV = lg(V)-1;
     870         588 :   P = FpXQ_powers(bi, lV-1, T, q);
     871         588 :   W = cgetg(lV+1, t_VEC);
     872       11081 :   for(i=1; i<=lV; i++)
     873       10493 :     gel(W, i) = FpXQX_FpXQ_mul(gel(V,i), gel(P, i), T, q);
     874         588 :   return gerepileupto(av, W);
     875             : }
     876             : 
     877             : struct _ZpXQ_sqrtn
     878             : {
     879             :   GEN T, a, n, ai;
     880             : };
     881             : 
     882             : static GEN
     883        2037 : _sqrtn_invd(void *E, GEN V, GEN v, GEN q, long M)
     884             : {
     885        2037 :   struct _ZpXQ_sqrtn *d = (struct _ZpXQ_sqrtn *) E;
     886        2037 :   GEN Tq = FpX_red(d->T, q), aiq = FpX_red(d->ai, q);
     887             :   (void)M;
     888        2037 :   return FpXQ_mul(FpXQ_mul(V, gel(v,2), Tq, q), aiq, Tq, q);
     889             : }
     890             : 
     891             : static GEN
     892        2037 : _sqrtn_eval(void *E, GEN x, GEN q)
     893             : {
     894        2037 :   struct _ZpXQ_sqrtn *d = (struct _ZpXQ_sqrtn *) E;
     895        2037 :   GEN Tq = FpX_red(d->T, q);
     896        2037 :   GEN f = FpX_sub(FpXQ_pow(x, d->n, Tq, q), d->a, q);
     897        2037 :   return mkvec2(f, x);
     898             : }
     899             : 
     900             : GEN
     901        1204 : ZpXQ_sqrtnlift(GEN a, GEN n, GEN x, GEN T, GEN p, long e)
     902             : {
     903             :   struct _ZpXQ_sqrtn d;
     904        1204 :   d.a = a; d.T = T; d.n = n;
     905        1204 :   d.ai = ZpXQ_inv(ZX_Z_mul(a, n),T,p,(e+1)>>1);
     906        1204 :   return gen_ZpX_Newton(x, p, e, &d, _sqrtn_eval, _sqrtn_invd);
     907             : }
     908             : 
     909             : static GEN
     910          98 : to_ZX(GEN a, long v) { return typ(a)==t_INT? scalarpol(a,v): a; }
     911             : 
     912             : GEN
     913         443 : Zq_sqrtnlift(GEN a, GEN n, GEN x, GEN T, GEN p, long e)
     914             : {
     915          49 :   return T? ZpXQ_sqrtnlift(to_ZX(a,varn(T)), n, to_ZX(x,varn(T)), T, p, e)
     916         492 :           : Zp_sqrtnlift(a, n, x, p, e);
     917             : }
     918             : 
     919             : GEN
     920           0 : ZpXQ_sqrt(GEN a, GEN T, GEN p, long e)
     921             : {
     922           0 :   pari_sp av = avma;
     923           0 :   GEN z = FpXQ_sqrt(FpX_red(a, p), T, p);
     924           0 :   if (!z) return NULL;
     925           0 :   if (e <= 1) return gerepileupto(av, z);
     926           0 :   return gerepileupto(av, ZpXQ_sqrtnlift(a, gen_2, z, T, p, e));
     927             : }
     928             : 
     929             : GEN
     930        4187 : ZpX_ZpXQ_liftroot_ea(GEN P, GEN S, GEN T, GEN p, long n, void *E,
     931             :                      GEN early(void *E, GEN x, GEN q))
     932             : {
     933        4187 :   pari_sp ltop = avma, av;
     934             :   long N, r;
     935             :   long mask;
     936             :   GEN q2, q, W, Q, Tq2, Tq, Pq;
     937             :   pari_timer ti;
     938        4187 :   T = FpX_get_red(T, powiu(p, n));
     939        4187 :   if (n == 1) return gcopy(S);
     940        4187 :   mask = quadratic_prec_mask(n);
     941        4187 :   av = avma;
     942        4187 :   q2 = p; q = sqri(p); mask >>= 1; N = 2;
     943        4187 :   if (DEBUGLEVEL > 3) timer_start(&ti);
     944        4187 :   Tq = FpXT_red(T,q);
     945        4187 :   Tq2 = FpXT_red(Tq,q2);
     946        4187 :   Pq = FpX_red(P,q);
     947        4187 :   W = FpXQ_inv(FpX_FpXQ_eval(FpX_deriv(P,q2), S, Tq2, q2), Tq2, q2);
     948        4187 :   Q  = ZX_Z_divexact(FpX_FpXQ_eval(Pq, S, Tq, q), q2);
     949        4187 :   r = brent_kung_optpow(degpol(P), 4, 3);
     950        4187 :   if (DEBUGLEVEL > 3)
     951           0 :     err_printf("ZpX_ZpXQ_liftroot: lifting to prec %ld\n",N);
     952             :   for (;;)
     953        8081 :   {
     954             :     GEN H, Sq, Wq, Spow, dP, qq, Pqq, Tqq;
     955       12268 :     H  = FpXQ_mul(W, Q, Tq2, q2);
     956       12268 :     Sq = FpX_sub(S, ZX_Z_mul(H, q2), q);
     957       12268 :     if (DEBUGLEVEL > 3)
     958           0 :       timer_printf(&ti,"ZpX_ZpXQ_liftroot: reaching prec %ld",N);
     959       12268 :     if (mask==1)
     960        2857 :       return gerepileupto(ltop, Sq);
     961        9411 :     if (early)
     962             :     {
     963        7927 :       GEN Se = early(E, Sq, q);
     964        7927 :       if (Se) return gerepileupto(ltop, Se);
     965             :     }
     966        8081 :     qq = sqri(q); N <<= 1;
     967        8081 :     if (mask&1UL) { qq = diviiexact(qq, p); N--; }
     968        8081 :     mask >>= 1;
     969        8081 :     Pqq  = FpX_red(P, qq);
     970        8081 :     Tqq  = FpXT_red(T, qq);
     971        8081 :     Spow = FpXQ_powers(Sq, r, Tqq, qq);
     972        8081 :     Q  = ZX_Z_divexact(FpX_FpXQV_eval(Pqq, Spow, Tqq, qq), q);
     973        8081 :     dP = FpX_FpXQV_eval(FpX_deriv(Pq, q), FpXV_red(Spow, q), Tq, q);
     974        8081 :     Wq = ZX_Z_divexact(FpX_Fp_sub(FpXQ_mul(W, dP, Tq, q), gen_1, q), q2);
     975        8081 :     Wq = ZX_Z_mul(FpXQ_mul(W, Wq, Tq2, q2), q2);
     976        8081 :     Wq = FpX_sub(W, Wq, q);
     977        8081 :     S = Sq; W = Wq; q2 = q; q = qq; Tq2 = Tq; Tq = Tqq; Pq = Pqq;
     978        8081 :     if (gc_needed(av, 1))
     979             :     {
     980           1 :       if(DEBUGMEM>1) pari_warn(warnmem,"ZpX_ZpXQ_Newton");
     981           1 :       gerepileall(av, 8, &S, &W, &Q, &Tq2, &Tq, &Pq, &q, &q2);
     982             :     }
     983             :   }
     984             : }
     985             : 
     986             : GEN
     987        1680 : ZpX_ZpXQ_liftroot(GEN P, GEN S, GEN T, GEN p, long n)
     988             : {
     989        1680 :   return ZpX_ZpXQ_liftroot_ea(P, S, T, p, n, NULL, NULL);
     990             : }
     991             : 
     992             : GEN
     993         287 : ZpX_Frobenius(GEN T, GEN p, long e)
     994             : {
     995         287 :   return ZpX_ZpXQ_liftroot(get_FpX_mod(T), FpX_Frobenius(T, p), T, p, e);
     996             : }
     997             : 
     998             : GEN
     999         140 : ZpXQM_prodFrobenius(GEN M, GEN T, GEN p, long e)
    1000             : {
    1001         140 :   pari_sp av = avma;
    1002         140 :   GEN xp = ZpX_Frobenius(T, p, e);
    1003         140 :   GEN z = FpXQM_autsum(mkvec2(xp, M), get_FpX_degree(T), T, powiu(p,e));
    1004         140 :   return gerepilecopy(av, gel(z,2));
    1005             : }
    1006             : 
    1007             : /* Canonical lift of polynomial */
    1008             : 
    1009       11130 : static GEN _can_invl(void *E, GEN V) {(void) E; return V; }
    1010             : 
    1011        3654 : static GEN _can_lin(void *E, GEN F, GEN V, GEN q)
    1012             : {
    1013        3654 :   GEN v = RgX_splitting(V, 3);
    1014             :   (void) E;
    1015        3654 :   return FpX_sub(V,ZXV_dotproduct(v, F), q);
    1016             : }
    1017             : 
    1018             : static GEN
    1019        7476 : _can_iter(void *E, GEN f, GEN q)
    1020             : {
    1021        7476 :   GEN h = RgX_splitting(f,3);
    1022        7476 :   GEN h1s = ZX_sqr(gel(h,1)), h2s = ZX_sqr(gel(h,2)), h3s = ZX_sqr(gel(h,3));
    1023        7476 :   GEN h12 = ZX_mul(gel(h,1), gel(h,2));
    1024        7476 :   GEN h13 = ZX_mul(gel(h,1), gel(h,3));
    1025        7476 :   GEN h23 = ZX_mul(gel(h,2), gel(h,3));
    1026        7476 :   GEN h1c = ZX_mul(gel(h,1), h1s);
    1027        7476 :   GEN h3c = ZX_mul(gel(h,3), h3s);
    1028        7476 :   GEN th = ZX_mul(ZX_sub(h2s,ZX_mulu(h13,3)),gel(h,2));
    1029        7476 :   GEN y = FpX_sub(f,ZX_add(RgX_shift_shallow(h3c,2),ZX_add(RgX_shift_shallow(th,1),h1c)),q);
    1030             :   (void) E;
    1031        7476 :   return mkvecn(7,y,h1s,h2s,h3s,h12,h13,h23);
    1032             : }
    1033             : 
    1034             : static GEN
    1035        7476 : _can_invd(void *E, GEN V, GEN v, GEN qM, long M)
    1036             : {
    1037        7476 :   GEN h1s=gel(v,2), h2s=gel(v,3), h3s=gel(v,4);
    1038        7476 :   GEN h12=gel(v,5), h13=gel(v,6), h23=gel(v,7);
    1039        7476 :   GEN F = mkvec3(ZX_sub(h1s,RgX_shift_shallow(h23,1)),RgX_shift_shallow(ZX_sub(h2s,h13),1),
    1040             :                  ZX_sub(RgX_shift_shallow(h3s,2),RgX_shift_shallow(h12,1)));
    1041             :   (void)E;
    1042        7476 :   return gen_ZpX_Dixon(ZXV_Z_mul(F, utoi(3)), V, qM, utoi(3), M, NULL,
    1043             :                                                  _can_lin, _can_invl);
    1044             : }
    1045             : 
    1046             : static GEN
    1047        3717 : F3x_frobeniuslift(GEN P, long n)
    1048        3717 : { return gen_ZpX_Newton(Flx_to_ZX(P),utoi(3), n, NULL, _can_iter, _can_invd); }
    1049             : 
    1050       29736 : static GEN _can5_invl(void *E, GEN V) {(void) E; return V; }
    1051             : 
    1052        9107 : static GEN _can5_lin(void *E, GEN F, GEN V, GEN q)
    1053             : {
    1054        9107 :   ulong p = *(ulong*)E;
    1055        9107 :   GEN v = RgX_splitting(V, p);
    1056        9107 :   return FpX_sub(V,ZXV_dotproduct(v, F), q);
    1057             : }
    1058             : 
    1059             : /* P(X,t) -> P(X*t^n,t) mod (t^p-1) */
    1060             : static GEN
    1061       62146 : _shift(GEN P, long n, ulong p, long v)
    1062             : {
    1063       62146 :   long i, l=lg(P);
    1064       62146 :   GEN r = cgetg(l,t_POL); r[1] = P[1];
    1065      483238 :   for(i=2;i<l;i++)
    1066             :   {
    1067      421092 :     long s = n*(i-2)%p;
    1068      421092 :     GEN ci = gel(P,i);
    1069      421092 :     if (typ(ci)==t_INT)
    1070      104979 :       gel(r,i) = monomial(ci, s, v);
    1071             :     else
    1072      316113 :       gel(r,i) = RgX_rotate_shallow(ci, s, p);
    1073             :   }
    1074       62146 :   return FpXX_renormalize(r, l);
    1075             : }
    1076             : 
    1077             : struct _can_mul
    1078             : {
    1079             :   GEN T, q;
    1080             :   ulong p;
    1081             : };
    1082             : 
    1083             : static GEN
    1084       41517 : _can5_mul(void *E, GEN A, GEN B)
    1085             : {
    1086       41517 :   struct _can_mul *d = (struct _can_mul *)E;
    1087       41517 :   GEN a = gel(A,1), b = gel(B,1);
    1088       41517 :   long n = itos(gel(A,2));
    1089       41517 :   GEN bn = _shift(b, n, d->p, get_FpX_var(d->T));
    1090       41517 :   GEN c = FpXQX_mul(a, bn, d->T, d->q);
    1091       41517 :   return mkvec2(c, addii(gel(A,2), gel(B,2)));
    1092             : }
    1093             : 
    1094             : static GEN
    1095       41349 : _can5_sqr(void *E, GEN A)
    1096             : {
    1097       41349 :   return _can5_mul(E,A,A);
    1098             : }
    1099             : 
    1100             : static GEN
    1101       20629 : ZXX_evalx0(GEN y)
    1102             : {
    1103       20629 :   long i, l = lg(y);
    1104       20629 :   GEN z = cgetg(l,t_POL); z[1] = y[1];
    1105      400260 :   for(i=2; i<l; i++)
    1106             :   {
    1107      379631 :     GEN yi = gel(y,i);
    1108      379631 :     gel(z,i) = typ(yi)==t_INT? yi: constant_coeff(yi);
    1109             :   }
    1110       20629 :   return ZX_renormalize(z,l);
    1111             : }
    1112             : 
    1113             : static GEN
    1114       20629 : _can5_iter(void *E, GEN f, GEN q)
    1115             : {
    1116       20629 :   pari_sp av = avma;
    1117             :   struct _can_mul D;
    1118       20629 :   ulong p = *(ulong*)E;
    1119       20629 :   long i, vT = fetch_var();
    1120             :   GEN N, P, d, V, fs;
    1121       20629 :   D.q = q; D.T = ZX_Z_sub(pol_xn(p,vT),gen_1);
    1122       20629 :   D.p = p;
    1123       20629 :   fs = mkvec2(_shift(f, 1, p, vT), gen_1);
    1124       20629 :   N = gel(gen_powu(fs,p-1,(void*)&D,_can5_sqr,_can5_mul),1);
    1125       20629 :   N = ZXX_evalx0(FpXQX_red(N,polcyclo(p,vT),q));
    1126       20629 :   P = FpX_mul(N,f,q);
    1127       20629 :   P = RgX_deflate(P, p);
    1128       20629 :   d = RgX_splitting(N, p);
    1129       20629 :   V = cgetg(p+1,t_VEC);
    1130       20629 :   gel(V,1) = ZX_mulu(gel(d,1), p);
    1131      103915 :   for(i=2; i<= (long)p; i++)
    1132       83286 :     gel(V,i) = ZX_mulu(RgX_shift_shallow(gel(d,p+2-i), 1), p);
    1133       20629 :   (void)delete_var(); return gerepilecopy(av, mkvec2(ZX_sub(f,P),V));
    1134             : }
    1135             : 
    1136             : static GEN
    1137       20629 : _can5_invd(void *E, GEN H, GEN v, GEN qM, long M)
    1138             : {
    1139       20629 :   ulong p = *(long*)E;
    1140       20629 :   return gen_ZpX_Dixon(gel(v,2), H, qM, utoi(p), M, E, _can5_lin, _can5_invl);
    1141             : }
    1142             : 
    1143             : GEN
    1144       13958 : Flx_Teichmuller(GEN P, ulong p, long n)
    1145             : {
    1146       24199 :   return p==3 ? F3x_frobeniuslift(P,n):
    1147       10241 :          gen_ZpX_Newton(Flx_to_ZX(P),utoi(p), n, &p, _can5_iter, _can5_invd);
    1148             : }
    1149             : 
    1150             : GEN
    1151          35 : polteichmuller(GEN P, ulong p, long n)
    1152             : {
    1153          35 :   pari_sp av = avma;
    1154          35 :   GEN q = NULL;
    1155          35 :   if (typ(P)!=t_POL || !RgX_is_FpX(P,&q))
    1156           0 :     pari_err_TYPE("polteichmuller",P);
    1157          35 :   if (q && !equaliu(q,p))
    1158           0 :     pari_err_MODULUS("polteichmuller",q,utoi(p));
    1159          35 :   if (n <= 0)
    1160           0 :     pari_err_DOMAIN("polteichmuller", "precision", "<=",gen_0,stoi(n));
    1161          63 :   return gerepileupto(av, p==2 ? F2x_Teichmuller(RgX_to_F2x(P), n)
    1162          28 :                                : Flx_Teichmuller(RgX_to_Flx(P, p), p, n));
    1163             : }

Generated by: LCOV version 1.13