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 - nffactor.c (source / functions) Hit Total Coverage
Test: PARI/GP v2.10.0 lcov report (development 20459-9710128) Lines: 1110 1207 92.0 %
Date: 2017-04-27 05:33:52 Functions: 65 69 94.2 %
Legend: Lines: hit not hit

          Line data    Source code
       1             : /* Copyright (C) 2000-2004  The PARI group.
       2             : 
       3             : This file is part of the PARI/GP package.
       4             : 
       5             : PARI/GP is free software; you can redistribute it and/or modify it under the
       6             : terms of the GNU General Public License as published by the Free Software
       7             : Foundation. It is distributed in the hope that it will be useful, but WITHOUT
       8             : ANY WARRANTY WHATSOEVER.
       9             : 
      10             : Check the License for details. You should have received a copy of it, along
      11             : with the package; see the file 'COPYING'. If not, write to the Free Software
      12             : Foundation, Inc., 51 Franklin Street, Fifth Floor, Boston, MA 02110-1301 USA. */
      13             : 
      14             : /*******************************************************************/
      15             : /*                                                                 */
      16             : /*            POLYNOMIAL FACTORIZATION IN A NUMBER FIELD           */
      17             : /*                                                                 */
      18             : /*******************************************************************/
      19             : #include "pari.h"
      20             : #include "paripriv.h"
      21             : 
      22             : static GEN nfsqff(GEN nf,GEN pol,long fl,GEN den);
      23             : static int nfsqff_use_Trager(long n, long dpol);
      24             : 
      25             : enum { FACTORS = 0, ROOTS, ROOTS_SPLIT };
      26             : 
      27             : /* for nf_bestlift: reconstruction of algebraic integers known mod P^k,
      28             :  * P maximal ideal above p */
      29             : typedef struct {
      30             :   long k;    /* input known mod P^k */
      31             :   GEN p, pk; /* p^k */
      32             :   GEN den;   /* denom(prk^-1) = p^k [ assume pr unramified ] */
      33             :   GEN prk;   /* |.|^2 LLL-reduced basis (b_i) of P^k  (NOT T2-reduced) */
      34             :   GEN prkHNF;/* HNF basis of P^k */
      35             :   GEN iprk;  /* den * prk^-1 */
      36             :   GEN GSmin; /* min |b_i^*|^2 */
      37             : 
      38             :   GEN Tp; /* Tpk mod p */
      39             :   GEN Tpk;
      40             :   GEN ZqProj;/* projector to Zp / P^k = Z/p^k[X] / Tpk */
      41             : 
      42             :   GEN tozk;
      43             :   GEN topow;
      44             :   GEN topowden; /* topow x / topowden = basistoalg(x) */
      45             :   GEN dn; /* NULL (we trust nf.zk) or a t_INT > 1 (an alg. integer has
      46             :              denominator dividing dn, when expressed on nf.zk */
      47             : } nflift_t;
      48             : 
      49             : typedef struct
      50             : {
      51             :   nflift_t *L;
      52             :   GEN nf;
      53             :   GEN pol, polbase; /* leading coeff is a t_INT */
      54             :   GEN fact;
      55             :   GEN pr;
      56             :   GEN Br, bound, ZC, BS_2;
      57             : } nfcmbf_t;
      58             : 
      59             : /*******************************************************************/
      60             : /*              RATIONAL RECONSTRUCTION (use ratlift)              */
      61             : /*******************************************************************/
      62             : /* NOT stack clean. a, b stay on the stack */
      63             : static GEN
      64     2973606 : lift_to_frac(GEN t, GEN mod, GEN amax, GEN bmax, GEN denom)
      65             : {
      66             :   GEN a, b;
      67     2973606 :   if (signe(t) < 0) t = addii(t, mod); /* in case t is a centerlift */
      68     2973606 :   if (!Fp_ratlift(t, mod, amax,bmax, &a,&b)
      69     2965920 :      || (denom && !dvdii(denom,b))
      70     2964720 :      || !is_pm1(gcdii(a,b))) return NULL;
      71     2964624 :   if (is_pm1(b)) { cgiv(b); return a; }
      72     1317006 :   return mkfrac(a, b);
      73             : }
      74             : 
      75             : /* Compute rational lifting for all the components of M modulo mod.
      76             :  * Assume all Fp_ratlift preconditions are met; we allow centerlifts but in
      77             :  * that case are no longer stack clean. If one component fails, return NULL.
      78             :  * If denom != NULL, check that the denominators divide denom.
      79             :  *
      80             :  * We suppose gcd(mod, denom) = 1, then a and b are coprime; so we can use
      81             :  * mkfrac rather than gdiv */
      82             : GEN
      83      100992 : FpM_ratlift(GEN M, GEN mod, GEN amax, GEN bmax, GEN denom)
      84             : {
      85      100992 :   pari_sp av = avma;
      86      100992 :   long i, j, h, l = lg(M);
      87      100992 :   GEN a, N = cgetg_copy(M, &l);
      88      100992 :   if (l == 1) return N;
      89       63786 :   h = lgcols(M);
      90      173190 :   for (j = 1; j < l; ++j)
      91             :   {
      92      114204 :     gel(N,j) = cgetg(h, t_COL);
      93     1969650 :     for (i = 1; i < h; ++i)
      94             :     {
      95     1860246 :       a = lift_to_frac(gcoeff(M,i,j), mod, amax,bmax,denom);
      96     1860246 :       if (!a) { avma = av; return NULL; }
      97     1855446 :       gcoeff(N,i,j) = a;
      98             :     }
      99             :   }
     100       58986 :   return N;
     101             : }
     102             : GEN
     103      362910 : FpC_ratlift(GEN P, GEN mod, GEN amax, GEN bmax, GEN denom)
     104             : {
     105      362910 :   pari_sp ltop = avma;
     106             :   long j, l;
     107      362910 :   GEN a, Q = cgetg_copy(P, &l);
     108     1464858 :   for (j = 1; j < l; ++j)
     109             :   {
     110     1101948 :     a = lift_to_frac(gel(P,j), mod, amax,bmax,denom);
     111     1101948 :     if (!a) { avma = ltop; return NULL; }
     112     1101948 :     gel(Q,j) = a;
     113             :   }
     114      362910 :   return Q;
     115             : }
     116             : GEN
     117        5142 : FpX_ratlift(GEN P, GEN mod, GEN amax, GEN bmax, GEN denom)
     118             : {
     119        5142 :   pari_sp ltop = avma;
     120             :   long j, l;
     121        5142 :   GEN a, Q = cgetg_copy(P, &l);
     122        5142 :   Q[1] = P[1];
     123       12372 :   for (j = 2; j < l; ++j)
     124             :   {
     125       11412 :     a = lift_to_frac(gel(P,j), mod, amax,bmax,denom);
     126       11412 :     if (!a) { avma = ltop; return NULL; }
     127        7230 :     gel(Q,j) = a;
     128             :   }
     129         960 :   return Q;
     130             : }
     131             : 
     132             : /*******************************************************************/
     133             : /*              GCD in K[X], K NUMBER FIELD                        */
     134             : /*******************************************************************/
     135             : /* P,Q in Z[X,Y], T in Z[Y] irreducible. compute GCD in Q[Y]/(T)[X].
     136             :  *
     137             :  * M. Encarnacion "On a modular Algorithm for computing GCDs of polynomials
     138             :  * over number fields" (ISSAC'94).
     139             :  *
     140             :  * We procede as follows
     141             :  *  1:compute the gcd modulo primes discarding bad primes as they are detected.
     142             :  *  2:reconstruct the result via FpM_ratlift, stoping as soon as we get weird
     143             :  *    denominators.
     144             :  *  3:if FpM_ratlift succeeds, try the full division.
     145             :  * Suppose accuracy is insufficient to get the result right: FpM_ratlift will
     146             :  * rarely succeed, and even if it does the polynomial we get has sensible
     147             :  * coefficients, so the full division will not be too costly.
     148             :  *
     149             :  * If not NULL, den must be a multiple of the denominator of the gcd,
     150             :  * for example the discriminant of T.
     151             :  *
     152             :  * NOTE: if T is not irreducible, nfgcd may loop forever, esp. if gcd | T */
     153             : GEN
     154        2550 : nfgcd_all(GEN P, GEN Q, GEN T, GEN den, GEN *Pnew)
     155             : {
     156        2550 :   pari_sp btop, ltop = avma;
     157        2550 :   GEN lP, lQ, M, dsol, R, bo, sol, mod = NULL;
     158        2550 :   long vP = varn(P), vT = varn(T), dT = degpol(T), dM = 0, dR;
     159             :   forprime_t S;
     160             : 
     161        2550 :   if (!signe(P)) { if (Pnew) *Pnew = pol_0(vT); return gcopy(Q); }
     162        2550 :   if (!signe(Q)) { if (Pnew) *Pnew = pol_1(vT);   return gcopy(P); }
     163             :   /*Compute denominators*/
     164        2550 :   if (!den) den = ZX_disc(T);
     165        2550 :   lP = leading_coeff(P);
     166        2550 :   lQ = leading_coeff(Q);
     167        2550 :   if ( !((typ(lP)==t_INT && is_pm1(lP)) || (typ(lQ)==t_INT && is_pm1(lQ))) )
     168         618 :     den = mulii(den, gcdii(ZX_resultant(lP, T), ZX_resultant(lQ, T)));
     169             : 
     170        2550 :   init_modular_small(&S);
     171        2550 :   btop = avma;
     172             :   for(;;)
     173             :   {
     174        3378 :     ulong p = u_forprime_next(&S);
     175        3378 :     if (!p) pari_err_OVERFLOW("nfgcd [ran out of primes]");
     176             :     /*Discard primes dividing disc(T) or lc(PQ) */
     177        3378 :     if (!umodiu(den, p)) continue;
     178        3378 :     if (DEBUGLEVEL>5) err_printf("nfgcd: p=%lu\n",p);
     179             :     /*Discard primes when modular gcd does not exist*/
     180        3378 :     if ((R = FlxqX_safegcd(ZXX_to_FlxX(P,p,vT),
     181             :                            ZXX_to_FlxX(Q,p,vT),
     182           0 :                            ZX_to_Flx(T,p), p)) == NULL) continue;
     183        3378 :     dR = degpol(R);
     184        3378 :     if (dR == 0) { avma = ltop; if (Pnew) *Pnew = P; return pol_1(vP); }
     185        1350 :     if (mod && dR > dM) continue; /* p divides Res(P/gcd, Q/gcd). Discard. */
     186             : 
     187        1350 :     R = FlxX_to_Flm(R, dT);
     188             :     /* previous primes divided Res(P/gcd, Q/gcd)? Discard them. */
     189        1350 :     if (!mod || dR < dM) { M = ZM_init_CRT(R, p); mod = utoipos(p); dM = dR; continue; }
     190         828 :     (void)ZM_incremental_CRT(&M,R, &mod,p);
     191         828 :     if (gc_needed(btop, 1))
     192             :     {
     193           0 :       if (DEBUGMEM>1) pari_warn(warnmem,"nfgcd");
     194           0 :       gerepileall(btop, 2, &M, &mod);
     195             :     }
     196             :     /* I suspect it must be better to take amax > bmax*/
     197         828 :     bo = sqrti(shifti(mod, -1));
     198         828 :     if ((sol = FpM_ratlift(M, mod, bo, bo, den)) == NULL) continue;
     199         522 :     sol = RgM_to_RgXX(sol,vP,vT);
     200         522 :     dsol = Q_primpart(sol);
     201             : 
     202         522 :     if (!ZXQX_dvd(Q, dsol, T)) continue;
     203         522 :     if (Pnew)
     204             :     {
     205         102 :       *Pnew = RgXQX_pseudodivrem(P, dsol, T, &R);
     206         102 :       if (signe(R)) continue;
     207             :     }
     208             :     else
     209             :     {
     210         420 :       if (!ZXQX_dvd(P, dsol, T)) continue;
     211             :     }
     212         522 :     gerepileall(ltop, Pnew? 2: 1, &dsol, Pnew);
     213         522 :     return dsol; /* both remainders are 0 */
     214         828 :   }
     215             : }
     216             : GEN
     217        1110 : nfgcd(GEN P, GEN Q, GEN T, GEN den)
     218        1110 : { return nfgcd_all(P,Q,T,den,NULL); }
     219             : 
     220             : int
     221        1242 : nfissquarefree(GEN nf, GEN x)
     222             : {
     223        1242 :   pari_sp av = avma;
     224        1242 :   GEN g, y = RgX_deriv(x);
     225        1242 :   if (RgX_is_rational(x))
     226         576 :     g = QX_gcd(x, y);
     227             :   else
     228             :   {
     229         666 :     GEN T = get_nfpol(nf,&nf);
     230         666 :     x = Q_primpart( liftpol_shallow(x) );
     231         666 :     y = Q_primpart( liftpol_shallow(y) );
     232         666 :     g = nfgcd(x, y, T, nf? nf_get_index(nf): NULL);
     233             :   }
     234        1242 :   avma = av; return (degpol(g) == 0);
     235             : }
     236             : 
     237             : /*******************************************************************/
     238             : /*             FACTOR OVER (Z_K/pr)[X] --> FqX_factor              */
     239             : /*******************************************************************/
     240             : GEN
     241           6 : nffactormod(GEN nf, GEN x, GEN pr)
     242             : {
     243           6 :   long j, l, vx = varn(x), vn;
     244           6 :   pari_sp av = avma;
     245             :   GEN F, E, rep, xrd, modpr, T, p;
     246             : 
     247           6 :   nf = checknf(nf);
     248           6 :   vn = nf_get_varn(nf);
     249           6 :   if (typ(x)!=t_POL) pari_err_TYPE("nffactormod",x);
     250           6 :   if (varncmp(vx,vn) >= 0) pari_err_PRIORITY("nffactormod", x, ">=", vn);
     251             : 
     252           6 :   modpr = nf_to_Fq_init(nf, &pr, &T, &p);
     253           6 :   xrd = nfX_to_FqX(x, nf, modpr);
     254           6 :   rep = FqX_factor(xrd,T,p);
     255           6 :   settyp(rep, t_MAT);
     256           6 :   F = gel(rep,1); l = lg(F);
     257           6 :   E = gel(rep,2); settyp(E, t_COL);
     258          12 :   for (j = 1; j < l; j++) {
     259           6 :     gel(F,j) = FqX_to_nfX(gel(F,j), modpr);
     260           6 :     gel(E,j) = stoi(E[j]);
     261             :   }
     262           6 :   return gerepilecopy(av, rep);
     263             : }
     264             : 
     265             : /*******************************************************************/
     266             : /*               MAIN ROUTINES nfroots / nffactor                  */
     267             : /*******************************************************************/
     268             : static GEN
     269        1980 : QXQX_normalize(GEN P, GEN T)
     270             : {
     271        1980 :   GEN P0 = leading_coeff(P);
     272        1980 :   long t = typ(P0);
     273        1980 :   if (t == t_POL)
     274             :   {
     275         504 :     if (degpol(P0)) return RgXQX_RgXQ_mul(P, QXQ_inv(P0,T), T);
     276         480 :     P0 = gel(P0,2); t = typ(P0);
     277             :   }
     278             :   /* t = t_INT/t_FRAC */
     279        1956 :   if (t == t_INT && is_pm1(P0) && signe(P0) > 0) return P; /* monic */
     280         726 :   return RgX_Rg_div(P, P0);
     281             : }
     282             : /* assume leading term of P is an integer */
     283             : static GEN
     284        1356 : RgX_int_normalize(GEN P)
     285             : {
     286        1356 :   GEN P0 = leading_coeff(P);
     287             :   /* cater for t_POL */
     288        1356 :   if (typ(P0) == t_POL)
     289             :   {
     290          48 :     P0 = gel(P0,2); /* non-0 constant */
     291          48 :     P = shallowcopy(P);
     292          48 :     gel(P,lg(P)-1) = P0; /* now leading term is a t_INT */
     293             :   }
     294        1356 :   if (typ(P0) != t_INT) pari_err_BUG("RgX_int_normalize");
     295        1356 :   if (is_pm1(P0)) return signe(P0) > 0? P: RgX_neg(P);
     296         810 :   return RgX_Rg_div(P, P0);
     297             : }
     298             : 
     299             : /* discard change of variable if nf is of the form [nf,c] as return by nfinit
     300             :  * for non-monic polynomials */
     301             : static GEN
     302         216 : proper_nf(GEN nf)
     303         216 : { return (lg(nf) == 3)? gel(nf,1): nf; }
     304             : 
     305             : /* if *pnf = NULL replace if by a "quick" K = nfinit(T), ensuring maximality
     306             :  * by small primes only. Return a multiplicative bound for the denominator of
     307             :  * algebraic integers in Z_K in terms of K.zk */
     308             : static GEN
     309        1260 : fix_nf(GEN *pnf, GEN *pT, GEN *pA)
     310             : {
     311        1260 :   GEN nf, NF, fa, P, Q, q, D, T = *pT;
     312             :   nfmaxord_t S;
     313             :   long i, l;
     314             : 
     315        1260 :   if (*pnf) return gen_1;
     316         216 :   nfmaxord(&S, T, nf_PARTIALFACT);
     317         216 :   NF = nfinit_complete(&S, 0, DEFAULTPREC);
     318         216 :   *pnf = nf = proper_nf(NF);
     319         216 :   if (nf != NF) { /* t_POL defining base field changed (not monic) */
     320          30 :     GEN A = *pA, a = cgetg_copy(A, &l);
     321          30 :     GEN rev = gel(NF,2), pow, dpow;
     322             : 
     323          30 :     *pT = T = nf_get_pol(nf); /* need to update T */
     324          30 :     pow = QXQ_powers(lift_shallow(rev), degpol(T)-1, T);
     325          30 :     pow = Q_remove_denom(pow, &dpow);
     326          30 :     a[1] = A[1];
     327         132 :     for (i=2; i<l; i++) {
     328         102 :       GEN c = gel(A,i);
     329         102 :       if (typ(c) == t_POL) c = QX_ZXQV_eval(c, pow, dpow);
     330         102 :       gel(a,i) = c;
     331             :     }
     332          30 :     *pA = Q_primpart(a); /* need to update A */
     333             :   }
     334             : 
     335         216 :   D = nf_get_disc(nf);
     336         216 :   if (is_pm1(D)) return gen_1;
     337         210 :   fa = absZ_factor_limit(D, 0);
     338         210 :   P = gel(fa,1); q = gel(P, lg(P)-1);
     339         210 :   if (BPSW_psp(q)) return gen_1;
     340             :   /* nf_get_disc(nf) may be incorrect */
     341           6 :   P = nf_get_ramified_primes(nf);
     342           6 :   l = lg(P);
     343           6 :   Q = q; q = gen_1;
     344          36 :   for (i = 1; i < l; i++)
     345             :   {
     346          30 :     GEN p = gel(P,i);
     347          30 :     if (Z_pvalrem(Q, p, &Q) && !BPSW_psp(p)) q = mulii(q, p);
     348             :   }
     349           6 :   return q;
     350             : }
     351             : 
     352             : /* lt(A) is an integer; ensure it is not a constant t_POL. In place */
     353             : static void
     354        1386 : ensure_lt_INT(GEN A)
     355             : {
     356        1386 :   long n = lg(A)-1;
     357        1386 :   GEN lt = gel(A,n);
     358        1386 :   while (typ(lt) != t_INT) gel(A,n) = lt = gel(lt,2);
     359        1386 : }
     360             : 
     361             : /* set B = A/gcd(A,A'), squarefree */
     362             : static GEN
     363        1290 : get_nfsqff_data(GEN *pnf, GEN *pT, GEN *pA, GEN *pB, GEN *ptbad)
     364             : {
     365        1290 :   GEN den, bad, D, B, A = *pA, T = *pT;
     366        1290 :   long n = degpol(T);
     367             : 
     368        1290 :   A = Q_primpart( QXQX_normalize(A, T) );
     369        1290 :   if (nfsqff_use_Trager(n, degpol(A)))
     370             :   {
     371         102 :     *pnf = T;
     372         102 :     bad = den = ZX_disc(T);
     373         102 :     if (is_pm1(leading_coeff(T))) den = indexpartial(T, den);
     374             :   }
     375             :   else
     376             :   {
     377        1188 :     den = fix_nf(pnf, &T, &A);
     378        1188 :     bad = nf_get_index(*pnf);
     379        1188 :     if (den != gen_1) bad = mulii(bad, den);
     380             :   }
     381        1290 :   D = nfgcd_all(A, RgX_deriv(A), T, bad, &B);
     382        1290 :   if (degpol(D)) B = Q_primpart( QXQX_normalize(B, T) );
     383        1290 :   if (ptbad) *ptbad = bad;
     384        1290 :   *pA = A;
     385        1290 :   *pB = B; ensure_lt_INT(B);
     386        1290 :   *pT = T; return den;
     387             : }
     388             : 
     389             : /* return the roots of pol in nf */
     390             : GEN
     391        1722 : nfroots(GEN nf,GEN pol)
     392             : {
     393        1722 :   pari_sp av = avma;
     394             :   GEN z, A, B, T, den;
     395             :   long d, dT;
     396             : 
     397        1722 :   if (!nf) return nfrootsQ(pol);
     398         930 :   T = get_nfpol(nf, &nf);
     399         930 :   RgX_check_ZX(T,"nfroots");
     400         930 :   A = RgX_nffix("nfroots", T,pol,1);
     401         930 :   d = degpol(A);
     402         930 :   if (d < 0) pari_err_ROOTS0("nfroots");
     403         930 :   if (d == 0) return cgetg(1,t_VEC);
     404         930 :   if (d == 1)
     405             :   {
     406           6 :     A = QXQX_normalize(A,T);
     407           6 :     A = mkpolmod(gneg_i(gel(A,2)), T);
     408           6 :     return gerepilecopy(av, mkvec(A));
     409             :   }
     410         924 :   dT = degpol(T);
     411         924 :   if (dT == 1) return gerepileupto(av, nfrootsQ(simplify_shallow(A)));
     412             : 
     413         924 :   den = get_nfsqff_data(&nf, &T, &A, &B, NULL);
     414         924 :   if (RgX_is_ZX(B))
     415             :   {
     416         312 :     GEN v = gel(ZX_factor(B), 1);
     417         312 :     long i, l = lg(v), p = mael(factoru(dT),1,1); /* smallest prime divisor */
     418         312 :     z = cgetg(1, t_VEC);
     419         852 :     for (i = 1; i < l; i++)
     420             :     {
     421         540 :       GEN b = gel(v,i); /* irreducible / Q */
     422         540 :       long db = degpol(b);
     423         540 :       if (db != 1 && degpol(b) < p) continue;
     424         540 :       z = shallowconcat(z, nfsqff(nf, b, ROOTS, den));
     425             :     }
     426             :   }
     427             :   else
     428         612 :     z = nfsqff(nf,B, ROOTS, den);
     429         924 :   z = gerepileupto(av, QXQV_to_mod(z, T));
     430         924 :   gen_sort_inplace(z, (void*)&cmp_RgX, &cmp_nodata, NULL);
     431         924 :   return z;
     432             : }
     433             : 
     434             : static GEN
     435      109116 : _norml2(GEN x) { return RgC_fpnorml2(x, DEFAULTPREC); }
     436             : 
     437             : /* return a minimal lift of elt modulo id, as a ZC */
     438             : static GEN
     439       31542 : nf_bestlift(GEN elt, GEN bound, nflift_t *L)
     440             : {
     441             :   GEN u;
     442       31542 :   long i,l = lg(L->prk), t = typ(elt);
     443       31542 :   if (t != t_INT)
     444             :   {
     445        8340 :     if (t == t_POL) elt = ZM_ZX_mul(L->tozk, elt);
     446        8340 :     u = ZM_ZC_mul(L->iprk,elt);
     447        8340 :     for (i=1; i<l; i++) gel(u,i) = diviiround(gel(u,i), L->den);
     448             :   }
     449             :   else
     450             :   {
     451       23202 :     u = ZC_Z_mul(gel(L->iprk,1), elt);
     452       23202 :     for (i=1; i<l; i++) gel(u,i) = diviiround(gel(u,i), L->den);
     453       23202 :     elt = scalarcol(elt, l-1);
     454             :   }
     455       31542 :   u = ZC_sub(elt, ZM_ZC_mul(L->prk, u));
     456       31542 :   if (bound && gcmp(_norml2(u), bound) > 0) u = NULL;
     457       31542 :   return u;
     458             : }
     459             : 
     460             : /* Warning: return L->topowden * (best lift). */
     461             : static GEN
     462       23850 : nf_bestlift_to_pol(GEN elt, GEN bound, nflift_t *L)
     463             : {
     464       23850 :   pari_sp av = avma;
     465       23850 :   GEN u,v = nf_bestlift(elt,bound,L);
     466       23850 :   if (!v) return NULL;
     467       19464 :   if (ZV_isscalar(v))
     468             :   {
     469       13494 :     if (L->topowden)
     470       13494 :       u = mulii(L->topowden, gel(v,1));
     471             :     else
     472           0 :       u = icopy(gel(v,1));
     473       13494 :     u = gerepileuptoint(av, u);
     474             :   }
     475             :   else
     476             :   {
     477        5970 :     v = gclone(v); avma = av;
     478        5970 :     u = RgV_dotproduct(L->topow, v);
     479        5970 :     gunclone(v);
     480             :   }
     481       19464 :   return u;
     482             : }
     483             : 
     484             : /* return the T->powden * (lift of pol with coefficients of T2-norm <= C)
     485             :  * if it exists. */
     486             : static GEN
     487        5400 : nf_pol_lift(GEN pol, GEN bound, nflift_t *L)
     488             : {
     489        5400 :   long i, l = lg(pol);
     490        5400 :   GEN x = cgetg(l,t_POL);
     491             : 
     492        5400 :   x[1] = pol[1];
     493        5400 :   gel(x,l-1) = mul_content(gel(pol,l-1), L->topowden);
     494       22506 :   for (i=l-2; i>1; i--)
     495             :   {
     496       21492 :     GEN t = nf_bestlift_to_pol(gel(pol,i), bound, L);
     497       21492 :     if (!t) return NULL;
     498       17106 :     gel(x,i) = t;
     499             :   }
     500        1014 :   return x;
     501             : }
     502             : 
     503             : static GEN
     504           0 : zerofact(long v)
     505             : {
     506           0 :   GEN z = cgetg(3, t_MAT);
     507           0 :   gel(z,1) = mkcol(pol_0(v));
     508           0 :   gel(z,2) = mkcol(gen_1); return z;
     509             : }
     510             : 
     511             : /* Return the factorization of A in Q[X]/(T) in rep [pre-allocated with
     512             :  * cgetg(3,t_MAT)], reclaiming all memory between avma and rep.
     513             :  * y is the vector of irreducible factors of B = Q_primpart( A/gcd(A,A') ).
     514             :  * Bad primes divide 'bad' */
     515             : static void
     516         462 : fact_from_sqff(GEN rep, GEN A, GEN B, GEN y, GEN T, GEN bad)
     517             : {
     518         462 :   pari_sp av = (pari_sp)rep;
     519         462 :   long n = lg(y)-1;
     520             :   GEN ex;
     521             : 
     522         462 :   if (A != B)
     523             :   { /* not squarefree */
     524          18 :     if (n == 1)
     525             :     { /* perfect power, simple ! */
     526           0 :       long e = degpol(A) / degpol(gel(y,1));
     527           0 :       y = gerepileupto(av, QXQXV_to_mod(y, T));
     528           0 :       ex = mkcol(utoipos(e));
     529             :     }
     530             :     else
     531             :     { /* compute valuations mod a prime of degree 1 (avoid coeff explosion) */
     532          18 :       GEN quo, p, r, Bp, lb = leading_coeff(B), E = cgetalloc(t_VECSMALL,n+1);
     533          18 :       pari_sp av1 = avma;
     534             :       ulong pp;
     535             :       long j;
     536             :       forprime_t S;
     537          18 :       u_forprime_init(&S, degpol(T), ULONG_MAX);
     538          90 :       for (; ; avma = av1)
     539             :       {
     540         108 :         pp = u_forprime_next(&S);
     541         108 :         if (! umodiu(bad,pp) || !umodiu(lb, pp)) continue;
     542          96 :         p = utoipos(pp);
     543          96 :         r = FpX_oneroot(T, p);
     544          96 :         if (!r) continue;
     545          42 :         Bp = FpXY_evalx(B, r, p);
     546          42 :         if (FpX_is_squarefree(Bp, p)) break;
     547          90 :       }
     548             : 
     549          18 :       quo = FpXY_evalx(Q_primpart(A), r, p);
     550          48 :       for (j=n; j>=2; j--)
     551             :       {
     552          30 :         GEN junk, fact = Q_remove_denom(gel(y,j), &junk);
     553          30 :         long e = 0;
     554          30 :         fact = FpXY_evalx(fact, r, p);
     555          84 :         for(;; e++)
     556             :         {
     557         114 :           GEN q = FpX_divrem(quo,fact,p,ONLY_DIVIDES);
     558         114 :           if (!q) break;
     559          84 :           quo = q;
     560          84 :         }
     561          30 :         E[j] = e;
     562             :       }
     563          18 :       E[1] = degpol(quo) / degpol(gel(y,1));
     564          18 :       y = gerepileupto(av, QXQXV_to_mod(y, T));
     565          18 :       ex = zc_to_ZC(E); pari_free((void*)E);
     566             :     }
     567             :   }
     568             :   else
     569             :   {
     570         444 :     y = gerepileupto(av, QXQXV_to_mod(y, T));
     571         444 :     ex = const_col(n, gen_1);
     572             :   }
     573         462 :   gel(rep,1) = y; settyp(y, t_COL);
     574         462 :   gel(rep,2) = ex;
     575         462 : }
     576             : 
     577             : /* return the factorization of x in nf */
     578             : GEN
     579         462 : nffactor(GEN nf,GEN pol)
     580             : {
     581         462 :   GEN bad, A, B, y, T, den, rep = cgetg(3, t_MAT);
     582         462 :   pari_sp av = avma;
     583             :   long dA;
     584             :   pari_timer ti;
     585             : 
     586         462 :   if (DEBUGLEVEL>2) { timer_start(&ti); err_printf("\nEntering nffactor:\n"); }
     587         462 :   T = get_nfpol(nf, &nf);
     588         462 :   RgX_check_ZX(T,"nffactor");
     589         462 :   A = RgX_nffix("nffactor",T,pol,1);
     590         456 :   dA = degpol(A);
     591         456 :   if (dA <= 0) {
     592           0 :     avma = (pari_sp)(rep + 3);
     593           0 :     return (dA == 0)? trivial_fact(): zerofact(varn(pol));
     594             :   }
     595         456 :   if (dA == 1) {
     596             :     GEN c;
     597          36 :     A = Q_primpart( QXQX_normalize(A, T) );
     598          36 :     A = gerepilecopy(av, A); c = gel(A,2);
     599          36 :     if (typ(c) == t_POL && degpol(c) > 0) gel(A,2) = mkpolmod(c, ZX_copy(T));
     600          36 :     gel(rep,1) = mkcol(A);
     601          36 :     gel(rep,2) = mkcol(gen_1); return rep;
     602             :   }
     603         420 :   if (degpol(T) == 1) return gerepileupto(av, QX_factor(simplify_shallow(A)));
     604             : 
     605         366 :   den = get_nfsqff_data(&nf, &T, &A, &B, &bad);
     606         366 :   if (DEBUGLEVEL>2) timer_printf(&ti, "squarefree test");
     607         366 :   if (RgX_is_ZX(B))
     608             :   {
     609         306 :     GEN v = gel(ZX_factor(B), 1);
     610         306 :     long i, l = lg(v);
     611         306 :     y = cgetg(1, t_VEC);
     612         630 :     for (i = 1; i < l; i++)
     613             :     {
     614         324 :       GEN b = gel(v,i); /* irreducible / Q */
     615         324 :       y = shallowconcat(y, nfsqff(nf, b, 0, den));
     616             :     }
     617             :   }
     618             :   else
     619          60 :     y = nfsqff(nf,B, 0, den);
     620         366 :   if (DEBUGLEVEL>3) err_printf("number of factor(s) found: %ld\n", lg(y)-1);
     621             : 
     622         366 :   fact_from_sqff(rep, A, B, y, T, bad);
     623         366 :   return sort_factor_pol(rep, cmp_RgX);
     624             : }
     625             : 
     626             : /* assume x scalar or t_COL, G t_MAT */
     627             : static GEN
     628       12150 : arch_for_T2(GEN G, GEN x)
     629             : {
     630       24300 :   return (typ(x) == t_COL)? RgM_RgC_mul(G,x)
     631       12150 :                           : RgC_Rg_mul(gel(G,1),x);
     632             : }
     633             : 
     634             : /* polbase a zkX with t_INT leading coeff; return a bound for T_2(P),
     635             :  * P | polbase in C[X]. NB: Mignotte bound: A | S ==>
     636             :  *  |a_i| <= binom(d-1, i-1) || S ||_2 + binom(d-1, i) lc(S)
     637             :  *
     638             :  * Apply to sigma(S) for all embeddings sigma, then take the L_2 norm over
     639             :  * sigma, then take the sup over i */
     640             : static GEN
     641         324 : nf_Mignotte_bound(GEN nf, GEN polbase)
     642         324 : { GEN lS = leading_coeff(polbase); /* t_INT */
     643             :   GEN p1, C, N2, binlS, bin;
     644         324 :   long prec = nf_get_prec(nf), n = nf_get_degree(nf), r1 = nf_get_r1(nf);
     645         324 :   long i, j, d = degpol(polbase);
     646             : 
     647         324 :   binlS = bin = vecbinomial(d-1);
     648         324 :   if (!isint1(lS)) binlS = ZC_Z_mul(bin,lS);
     649             : 
     650         324 :   N2 = cgetg(n+1, t_VEC);
     651             :   for (;;)
     652             :   {
     653         324 :     GEN G = nf_get_G(nf), matGS = cgetg(d+2, t_MAT);
     654             : 
     655         324 :     for (j=0; j<=d; j++) gel(matGS,j+1) = arch_for_T2(G, gel(polbase,j+2));
     656         324 :     matGS = shallowtrans(matGS);
     657         792 :     for (j=1; j <= r1; j++) /* N2[j] = || sigma_j(S) ||_2 */
     658             :     {
     659         468 :       GEN c = sqrtr( _norml2(gel(matGS,j)) );
     660         468 :       gel(N2,j) = c; if (!signe(c)) goto PRECPB;
     661             :     }
     662        1242 :     for (   ; j <= n; j+=2)
     663             :     {
     664         918 :       GEN q1 = _norml2(gel(matGS, j));
     665         918 :       GEN q2 = _norml2(gel(matGS, j+1));
     666         918 :       GEN c = sqrtr( gmul2n(addrr(q1, q2), -1) );
     667         918 :       gel(N2,j) = gel(N2,j+1) = c; if (!signe(c)) goto PRECPB;
     668             :     }
     669         324 :     break; /* done */
     670             : PRECPB:
     671           0 :     prec = precdbl(prec);
     672           0 :     nf = nfnewprec_shallow(nf, prec);
     673           0 :     if (DEBUGLEVEL>1) pari_warn(warnprec, "nf_factor_bound", prec);
     674           0 :   }
     675             : 
     676             :   /* Take sup over 0 <= i <= d of
     677             :    * sum_j | binom(d-1, i-1) ||sigma_j(S)||_2 + binom(d-1,i) lc(S) |^2 */
     678             : 
     679             :   /* i = 0: n lc(S)^2 */
     680         324 :   C = mului(n, sqri(lS));
     681             :   /* i = d: sum_sigma ||sigma(S)||_2^2 */
     682         324 :   p1 = gnorml2(N2); if (gcmp(C, p1) < 0) C = p1;
     683        7380 :   for (i = 1; i < d; i++)
     684             :   {
     685        7056 :     GEN B = gel(bin,i), L = gel(binlS,i+1);
     686        7056 :     GEN s = sqrr(addri(mulir(B, gel(N2,1)),  L)); /* j=1 */
     687        7056 :     for (j = 2; j <= n; j++) s = addrr(s, sqrr(addri(mulir(B, gel(N2,j)), L)));
     688        7056 :     if (mpcmp(C, s) < 0) C = s;
     689             :   }
     690         324 :   return C;
     691             : }
     692             : 
     693             : /* return a bound for T_2(P), P | polbase
     694             :  * max |b_i|^2 <= 3^{3/2 + d} / (4 \pi d) [P]_2,
     695             :  * where [P]_2 is Bombieri's 2-norm
     696             :  * Sum over conjugates */
     697             : static GEN
     698         324 : nf_Beauzamy_bound(GEN nf, GEN polbase)
     699             : {
     700             :   GEN lt, C, s, POL, bin;
     701         324 :   long d = degpol(polbase), n = nf_get_degree(nf), prec = nf_get_prec(nf);
     702         324 :   bin = vecbinomial(d);
     703         324 :   POL = polbase + 2;
     704             :   /* compute [POL]_2 */
     705             :   for (;;)
     706             :   {
     707         324 :     GEN G = nf_get_G(nf);
     708             :     long i;
     709             : 
     710         324 :     s = real_0(prec);
     711        8028 :     for (i=0; i<=d; i++)
     712             :     {
     713        7704 :       GEN c = gel(POL,i);
     714        7704 :       if (gequal0(c)) continue;
     715        4446 :       c = _norml2(arch_for_T2(G,c));
     716        4446 :       if (!signe(c)) goto PRECPB;
     717             :       /* s += T2(POL[i]) / binomial(d,i) */
     718        4446 :       s = addrr(s, divri(c, gel(bin,i+1)));
     719             :     }
     720         324 :     break;
     721             : PRECPB:
     722           0 :     prec = precdbl(prec);
     723           0 :     nf = nfnewprec_shallow(nf, prec);
     724           0 :     if (DEBUGLEVEL>1) pari_warn(warnprec, "nf_factor_bound", prec);
     725           0 :   }
     726         324 :   lt = leading_coeff(polbase);
     727         324 :   s = mulri(s, muliu(sqri(lt), n));
     728         324 :   C = powruhalf(stor(3,DEFAULTPREC), 3 + 2*d); /* 3^{3/2 + d} */
     729         324 :   return divrr(mulrr(C, s), mulur(d, mppi(DEFAULTPREC)));
     730             : }
     731             : 
     732             : static GEN
     733         324 : nf_factor_bound(GEN nf, GEN polbase)
     734             : {
     735         324 :   pari_sp av = avma;
     736         324 :   GEN a = nf_Mignotte_bound(nf, polbase);
     737         324 :   GEN b = nf_Beauzamy_bound(nf, polbase);
     738         324 :   if (DEBUGLEVEL>2)
     739             :   {
     740           0 :     err_printf("Mignotte bound: %Ps\n",a);
     741           0 :     err_printf("Beauzamy bound: %Ps\n",b);
     742             :   }
     743         324 :   return gerepileupto(av, gmin(a, b));
     744             : }
     745             : 
     746             : /* return Bs: if r a root of sigma_i(P), |r| < Bs[i] */
     747             : static GEN
     748        1014 : nf_root_bounds(GEN P, GEN T)
     749             : {
     750             :   long lR, i, j, l, prec;
     751             :   GEN Ps, R, V, nf;
     752             : 
     753        1014 :   if (RgX_is_rational(P)) return polrootsbound(P);
     754         540 :   T = get_nfpol(T, &nf);
     755             : 
     756         540 :   P = Q_primpart(P);
     757         540 :   prec = ZXX_max_lg(P) + 1;
     758         540 :   l = lg(P);
     759         540 :   if (nf && nf_get_prec(nf) >= prec)
     760         504 :     R = nf_get_roots(nf);
     761             :   else
     762          36 :     R = QX_complex_roots(T, prec);
     763         540 :   lR = lg(R);
     764         540 :   V = cgetg(lR, t_VEC);
     765         540 :   Ps = cgetg(l, t_POL); /* sigma (P) */
     766         540 :   Ps[1] = P[1];
     767        1584 :   for (j=1; j<lg(R); j++)
     768             :   {
     769        1044 :     GEN r = gel(R,j);
     770        1044 :     for (i=2; i<l; i++) gel(Ps,i) = poleval(gel(P,i), r);
     771        1044 :     gel(V,j) = polrootsbound(Ps);
     772             :   }
     773         540 :   return V;
     774             : }
     775             : 
     776             : /* return B such that if x in O_K, K = Z[X]/(T), then the L2-norm of the
     777             :  * coordinates of the numerator of x [on the power, resp. integral, basis if T
     778             :  * is a polynomial, resp. an nf] is  <= B T_2(x)
     779             :  * den = multiplicative bound for denom(x) */
     780             : static GEN
     781        1368 : L2_bound(GEN nf, GEN den)
     782             : {
     783        1368 :   GEN M, L, prep, T = nf_get_pol(nf), tozk = nf_get_invzk(nf);
     784        1368 :   long bit = bit_accuracy(ZX_max_lg(T)) + bit_accuracy(ZM_max_lg(tozk));
     785        1368 :   long prec = nbits2prec(bit + degpol(T));
     786        1368 :   (void)initgaloisborne(nf, den, prec, &L, &prep, NULL);
     787        1368 :   M = vandermondeinverse(L, RgX_gtofp(T,prec), den, prep);
     788        1368 :   return RgM_fpnorml2(RgM_mul(tozk,M), DEFAULTPREC);
     789             : }
     790             : 
     791             : /* || L ||_p^p in dimension n (L may be a scalar) */
     792             : static GEN
     793        1902 : normlp(GEN L, long p, long n)
     794             : {
     795        1902 :   long i,l, t = typ(L);
     796             :   GEN z;
     797             : 
     798        1902 :   if (!is_vec_t(t)) return gmulsg(n, gpowgs(L, p));
     799             : 
     800        1026 :   l = lg(L); z = gen_0;
     801             :   /* assert(n == l-1); */
     802        3024 :   for (i=1; i<l; i++)
     803        1998 :     z = gadd(z, gpowgs(gel(L,i), p));
     804        1026 :   return z;
     805             : }
     806             : 
     807             : /* S = S0 + tS1, P = P0 + tP1 (Euclidean div. by t integer). For a true
     808             :  * factor (vS, vP), we have:
     809             :  *    | S vS + P vP |^2 < Btra
     810             :  * This implies | S1 vS + P1 vP |^2 < Bhigh, assuming t > sqrt(Btra).
     811             :  * d = dimension of low part (= [nf:Q])
     812             :  * n0 = bound for |vS|^2
     813             :  * */
     814             : static double
     815         366 : get_Bhigh(long n0, long d)
     816             : {
     817         366 :   double sqrtd = sqrt((double)d);
     818         366 :   double z = n0*sqrtd + sqrtd/2 * (d * (n0+1));
     819         366 :   z = 1. + 0.5 * z; return z * z;
     820             : }
     821             : 
     822             : typedef struct {
     823             :   GEN d;
     824             :   GEN dPinvS;   /* d P^(-1) S   [ integral ] */
     825             :   double **PinvSdbl; /* P^(-1) S as double */
     826             :   GEN S1, P1;   /* S = S0 + S1 q, idem P */
     827             : } trace_data;
     828             : 
     829             : /* S1 * u - P1 * round(P^-1 S u). K non-zero coords in u given by ind */
     830             : static GEN
     831       80874 : get_trace(GEN ind, trace_data *T)
     832             : {
     833       80874 :   long i, j, l, K = lg(ind)-1;
     834             :   GEN z, s, v;
     835             : 
     836       80874 :   s = gel(T->S1, ind[1]);
     837       80874 :   if (K == 1) return s;
     838             : 
     839             :   /* compute s = S1 u */
     840       79086 :   for (j=2; j<=K; j++) s = ZC_add(s, gel(T->S1, ind[j]));
     841             : 
     842             :   /* compute v := - round(P^1 S u) */
     843       79086 :   l = lg(s);
     844       79086 :   v = cgetg(l, t_VECSMALL);
     845      958920 :   for (i=1; i<l; i++)
     846             :   {
     847      879834 :     double r, t = 0.;
     848             :     /* quick approximate computation */
     849      879834 :     for (j=1; j<=K; j++) t += T->PinvSdbl[ ind[j] ][i];
     850      879834 :     r = floor(t + 0.5);
     851      879834 :     if (fabs(t + 0.5 - r) < 0.0001)
     852             :     { /* dubious, compute exactly */
     853          84 :       z = gen_0;
     854          84 :       for (j=1; j<=K; j++) z = addii(z, ((GEN**)T->dPinvS)[ ind[j] ][i]);
     855          84 :       v[i] = - itos( diviiround(z, T->d) );
     856             :     }
     857             :     else
     858      879750 :       v[i] = - (long)r;
     859             :   }
     860       79086 :   return ZC_add(s, ZM_zc_mul(T->P1, v));
     861             : }
     862             : 
     863             : static trace_data *
     864         648 : init_trace(trace_data *T, GEN S, nflift_t *L, GEN q)
     865             : {
     866         648 :   long e = gexpo(S), i,j, l,h;
     867             :   GEN qgood, S1, invd;
     868             : 
     869         648 :   if (e < 0) return NULL; /* S = 0 */
     870             : 
     871         612 :   qgood = int2n(e - 32); /* single precision check */
     872         612 :   if (cmpii(qgood, q) > 0) q = qgood;
     873             : 
     874         612 :   S1 = gdivround(S, q);
     875         612 :   if (gequal0(S1)) return NULL;
     876             : 
     877         192 :   invd = invr(itor(L->den, DEFAULTPREC));
     878             : 
     879         192 :   T->dPinvS = ZM_mul(L->iprk, S);
     880         192 :   l = lg(S);
     881         192 :   h = lgcols(T->dPinvS);
     882         192 :   T->PinvSdbl = (double**)cgetg(l, t_MAT);
     883        2856 :   for (j = 1; j < l; j++)
     884             :   {
     885        2664 :     double *t = (double *) stack_malloc_align(h * sizeof(double), sizeof(double));
     886        2664 :     GEN c = gel(T->dPinvS,j);
     887        2664 :     pari_sp av = avma;
     888        2664 :     T->PinvSdbl[j] = t;
     889        2664 :     for (i=1; i < h; i++) t[i] = rtodbl(mulri(invd, gel(c,i)));
     890        2664 :     avma = av;
     891             :   }
     892             : 
     893         192 :   T->d  = L->den;
     894         192 :   T->P1 = gdivround(L->prk, q);
     895         192 :   T->S1 = S1; return T;
     896             : }
     897             : 
     898             : static void
     899       13392 : update_trace(trace_data *T, long k, long i)
     900             : {
     901       26784 :   if (!T) return;
     902        7800 :   gel(T->S1,k)     = gel(T->S1,i);
     903        7800 :   gel(T->dPinvS,k) = gel(T->dPinvS,i);
     904        7800 :   T->PinvSdbl[k]   = T->PinvSdbl[i];
     905             : }
     906             : 
     907             : /* reduce coeffs mod (T,pk), then center mod pk */
     908             : static GEN
     909       13314 : FqX_centermod(GEN z, GEN T, GEN pk, GEN pks2)
     910             : {
     911             :   long i, l;
     912             :   GEN y;
     913       13314 :   if (!T) return centermod_i(z, pk, pks2);
     914       12240 :   y = FpXQX_red(z, T, pk); l = lg(y);
     915      115362 :   for (i = 2; i < l; i++)
     916             :   {
     917      103122 :     GEN c = gel(y,i);
     918      103122 :     if (typ(c) == t_INT)
     919       68772 :       c = centermodii(c, pk, pks2);
     920             :     else
     921       34350 :       c = FpX_center(c, pk, pks2);
     922      103122 :     gel(y,i) = c;
     923             :   }
     924       12240 :   return y;
     925             : }
     926             : 
     927             : typedef struct {
     928             :   GEN lt, C, Clt, C2lt, C2ltpol;
     929             : } div_data;
     930             : 
     931             : static void
     932        1392 : init_div_data(div_data *D, GEN pol, nflift_t *L)
     933             : {
     934        1392 :   GEN C = mul_content(L->topowden, L->dn);
     935        1392 :   GEN C2lt, Clt, lc = leading_coeff(pol), lt = is_pm1(lc)? NULL: absi(lc);
     936        1392 :   if (C)
     937             :   {
     938        1392 :     GEN C2 = sqri(C);
     939        1392 :     if (lt) {
     940         306 :       C2lt = mulii(C2, lt);
     941         306 :       Clt = mulii(C,lt);
     942             :     } else {
     943        1086 :       C2lt = C2;
     944        1086 :       Clt = C;
     945             :     }
     946             :   }
     947             :   else
     948           0 :     C2lt = Clt = lt;
     949        1392 :   D->lt = lt;
     950        1392 :   D->C = C;
     951        1392 :   D->Clt = Clt;
     952        1392 :   D->C2lt = C2lt;
     953        1392 :   D->C2ltpol = C2lt? RgX_Rg_mul(pol, C2lt): pol;
     954        1392 : }
     955             : static void
     956         978 : update_target(div_data *D, GEN pol)
     957         978 : { D->C2ltpol = D->Clt? RgX_Rg_mul(pol, D->Clt): pol; }
     958             : 
     959             : /* nb = number of modular factors; return a "good" K such that naive
     960             :  * recombination of up to maxK modular factors is not too costly */
     961             : long
     962        8310 : cmbf_maxK(long nb)
     963             : {
     964        8310 :   if (nb >  10) return 3;
     965        7830 :   return nb-1;
     966             : }
     967             : /* Naive recombination of modular factors: combine up to maxK modular
     968             :  * factors, degree <= klim
     969             :  *
     970             :  * target = polynomial we want to factor
     971             :  * famod = array of modular factors.  Product should be congruent to
     972             :  * target/lc(target) modulo p^a
     973             :  * For true factors: S1,S2 <= p^b, with b <= a and p^(b-a) < 2^31 */
     974             : /* set *done = 1 if factorisation is known to be complete */
     975             : static GEN
     976         324 : nfcmbf(nfcmbf_t *T, long klim, long *pmaxK, int *done)
     977             : {
     978         324 :   GEN nf = T->nf, famod = T->fact, bound = T->bound;
     979         324 :   GEN ltdn, nfpol = nf_get_pol(nf);
     980         324 :   long K = 1, cnt = 1, i,j,k, curdeg, lfamod = lg(famod)-1, dnf = degpol(nfpol);
     981         324 :   pari_sp av0 = avma;
     982         324 :   GEN Tpk = T->L->Tpk, pk = T->L->pk, pks2 = shifti(pk,-1);
     983         324 :   GEN ind      = cgetg(lfamod+1, t_VECSMALL);
     984         324 :   GEN deg      = cgetg(lfamod+1, t_VECSMALL);
     985         324 :   GEN degsofar = cgetg(lfamod+1, t_VECSMALL);
     986         324 :   GEN fa       = cgetg(lfamod+1, t_VEC);
     987         324 :   const double Bhigh = get_Bhigh(lfamod, dnf);
     988             :   trace_data _T1, _T2, *T1, *T2;
     989             :   div_data D;
     990             :   pari_timer ti;
     991             : 
     992         324 :   timer_start(&ti);
     993             : 
     994         324 :   *pmaxK = cmbf_maxK(lfamod);
     995         324 :   init_div_data(&D, T->pol, T->L);
     996         324 :   ltdn = mul_content(D.lt, T->L->dn);
     997             :   {
     998         324 :     GEN q = ceil_safe(sqrtr(T->BS_2));
     999         324 :     GEN t1,t2, lt2dn = mul_content(ltdn, D.lt);
    1000         324 :     GEN trace1   = cgetg(lfamod+1, t_MAT);
    1001         324 :     GEN trace2   = cgetg(lfamod+1, t_MAT);
    1002        2778 :     for (i=1; i <= lfamod; i++)
    1003             :     {
    1004        2454 :       pari_sp av = avma;
    1005        2454 :       GEN P = gel(famod,i);
    1006        2454 :       long d = degpol(P);
    1007             : 
    1008        2454 :       deg[i] = d; P += 2;
    1009        2454 :       t1 = gel(P,d-1);/* = - S_1 */
    1010        2454 :       t2 = Fq_sqr(t1, Tpk, pk);
    1011        2454 :       if (d > 1) t2 = Fq_sub(t2, gmul2n(gel(P,d-2), 1), Tpk, pk);
    1012             :       /* t2 = S_2 Newton sum */
    1013        2454 :       if (ltdn)
    1014             :       {
    1015         108 :         t1 = Fq_Fp_mul(t1, ltdn, Tpk, pk);
    1016         108 :         t2 = Fq_Fp_mul(t2, lt2dn, Tpk, pk);
    1017             :       }
    1018        2454 :       gel(trace1,i) = gclone( nf_bestlift(t1, NULL, T->L) );
    1019        2454 :       gel(trace2,i) = gclone( nf_bestlift(t2, NULL, T->L) ); avma = av;
    1020             :     }
    1021         324 :     T1 = init_trace(&_T1, trace1, T->L, q);
    1022         324 :     T2 = init_trace(&_T2, trace2, T->L, q);
    1023        2778 :     for (i=1; i <= lfamod; i++) {
    1024        2454 :       gunclone(gel(trace1,i));
    1025        2454 :       gunclone(gel(trace2,i));
    1026             :     }
    1027             :   }
    1028         324 :   degsofar[0] = 0; /* sentinel */
    1029             : 
    1030             :   /* ind runs through strictly increasing sequences of length K,
    1031             :    * 1 <= ind[i] <= lfamod */
    1032             : nextK:
    1033         600 :   if (K > *pmaxK || 2*K > lfamod) goto END;
    1034         498 :   if (DEBUGLEVEL > 3)
    1035           0 :     err_printf("\n### K = %d, %Ps combinations\n", K,binomial(utoipos(lfamod), K));
    1036         498 :   setlg(ind, K+1); ind[1] = 1;
    1037         498 :   i = 1; curdeg = deg[ind[1]];
    1038             :   for(;;)
    1039             :   { /* try all combinations of K factors */
    1040       95478 :     for (j = i; j < K; j++)
    1041             :     {
    1042       10764 :       degsofar[j] = curdeg;
    1043       10764 :       ind[j+1] = ind[j]+1; curdeg += deg[ind[j+1]];
    1044             :     }
    1045       84714 :     if (curdeg <= klim) /* trial divide */
    1046             :     {
    1047             :       GEN t, y, q;
    1048             :       pari_sp av;
    1049             : 
    1050       84714 :       av = avma;
    1051       84714 :       if (T1)
    1052             :       { /* d-1 test */
    1053       41700 :         t = get_trace(ind, T1);
    1054       41700 :         if (rtodbl(_norml2(t)) > Bhigh)
    1055             :         {
    1056       40722 :           if (DEBUGLEVEL>6) err_printf(".");
    1057       40722 :           avma = av; goto NEXT;
    1058             :         }
    1059             :       }
    1060       43992 :       if (T2)
    1061             :       { /* d-2 test */
    1062       39174 :         t = get_trace(ind, T2);
    1063       39174 :         if (rtodbl(_norml2(t)) > Bhigh)
    1064             :         {
    1065       38646 :           if (DEBUGLEVEL>3) err_printf("|");
    1066       38646 :           avma = av; goto NEXT;
    1067             :         }
    1068             :       }
    1069        5346 :       avma = av;
    1070        5346 :       y = ltdn; /* full computation */
    1071       18660 :       for (i=1; i<=K; i++)
    1072             :       {
    1073       13314 :         GEN q = gel(famod, ind[i]);
    1074       13314 :         if (y) q = gmul(y, q);
    1075       13314 :         y = FqX_centermod(q, Tpk, pk, pks2);
    1076             :       }
    1077        5346 :       y = nf_pol_lift(y, bound, T->L);
    1078        5346 :       if (!y)
    1079             :       {
    1080        4386 :         if (DEBUGLEVEL>3) err_printf("@");
    1081        4386 :         avma = av; goto NEXT;
    1082             :       }
    1083             :       /* y = topowden*dn*lt*\prod_{i in ind} famod[i] is apparently in O_K[X],
    1084             :        * in fact in (Z[Y]/nf.pol)[X] due to multiplication by C = topowden*dn.
    1085             :        * Try out this candidate factor */
    1086         960 :       q = RgXQX_divrem(D.C2ltpol, y, nfpol, ONLY_DIVIDES);
    1087         960 :       if (!q)
    1088             :       {
    1089          24 :         if (DEBUGLEVEL>3) err_printf("*");
    1090          24 :         avma = av; goto NEXT;
    1091             :       }
    1092             :       /* Original T->pol in O_K[X] with leading coeff lt in Z,
    1093             :        * y = C*lt \prod famod[i] is in O_K[X] with leading coeff in Z
    1094             :        * q = C^2*lt*pol / y = C * (lt*pol) / (lt*\prod famod[i]) is a
    1095             :        * K-rational factor, in fact in Z[Y]/nf.pol)[X] as above, with
    1096             :        * leading term C*lt. */
    1097         936 :       update_target(&D, q);
    1098         936 :       gel(fa,cnt++) = D.C2lt? RgX_int_normalize(y): y; /* make monic */
    1099        8844 :       for (i=j=k=1; i <= lfamod; i++)
    1100             :       { /* remove used factors */
    1101        7908 :         if (j <= K && i == ind[j]) j++;
    1102             :         else
    1103             :         {
    1104        6696 :           gel(famod,k) = gel(famod,i);
    1105        6696 :           update_trace(T1, k, i);
    1106        6696 :           update_trace(T2, k, i);
    1107        6696 :           deg[k] = deg[i]; k++;
    1108             :         }
    1109             :       }
    1110         936 :       lfamod -= K;
    1111         936 :       *pmaxK = cmbf_maxK(lfamod);
    1112         936 :       if (lfamod < 2*K) goto END;
    1113         714 :       i = 1; curdeg = deg[ind[1]];
    1114         714 :       if (DEBUGLEVEL > 2)
    1115             :       {
    1116           0 :         err_printf("\n"); timer_printf(&ti, "to find factor %Ps",y);
    1117           0 :         err_printf("remaining modular factor(s): %ld\n", lfamod);
    1118             :       }
    1119         714 :       continue;
    1120             :     }
    1121             : 
    1122             : NEXT:
    1123       83778 :     for (i = K+1;;)
    1124             :     {
    1125       94542 :       if (--i == 0) { K++; goto nextK; }
    1126       94266 :       if (++ind[i] <= lfamod - K + i)
    1127             :       {
    1128       83502 :         curdeg = degsofar[i-1] + deg[ind[i]];
    1129       83502 :         if (curdeg <= klim) break;
    1130             :       }
    1131       10764 :     }
    1132       84216 :   }
    1133             : END:
    1134         324 :   *done = 1;
    1135         324 :   if (degpol(D.C2ltpol) > 0)
    1136             :   { /* leftover factor */
    1137         324 :     GEN q = D.C2ltpol;
    1138         324 :     if (D.C2lt) q = RgX_int_normalize(q);
    1139         324 :     if (lfamod >= 2*K)
    1140             :     { /* restore leading coefficient [#930] */
    1141          42 :       if (D.lt) q = RgX_Rg_mul(q, D.lt);
    1142          42 :       *done = 0; /* ... may still be reducible */
    1143             :     }
    1144         324 :     setlg(famod, lfamod+1);
    1145         324 :     gel(fa,cnt++) = q;
    1146             :   }
    1147         324 :   if (DEBUGLEVEL>6) err_printf("\n");
    1148         324 :   if (cnt == 2) {
    1149          36 :     avma = av0;
    1150          36 :     return mkvec(T->pol);
    1151             :   }
    1152             :   else
    1153             :   {
    1154         288 :     setlg(fa, cnt);
    1155         288 :     return gerepilecopy(av0, fa);
    1156             :   }
    1157             : }
    1158             : 
    1159             : static GEN
    1160          24 : nf_chk_factors(nfcmbf_t *T, GEN P, GEN M_L, GEN famod, GEN pk)
    1161             : {
    1162          24 :   GEN nf = T->nf, bound = T->bound;
    1163          24 :   GEN nfT = nf_get_pol(nf);
    1164             :   long i, r;
    1165          24 :   GEN pol = P, list, piv, y;
    1166          24 :   GEN Tpk = T->L->Tpk;
    1167             :   div_data D;
    1168             : 
    1169          24 :   piv = ZM_hnf_knapsack(M_L);
    1170          24 :   if (!piv) return NULL;
    1171          12 :   if (DEBUGLEVEL>3) err_printf("ZM_hnf_knapsack output:\n%Ps\n",piv);
    1172             : 
    1173          12 :   r  = lg(piv)-1;
    1174          12 :   list = cgetg(r+1, t_VEC);
    1175          12 :   init_div_data(&D, pol, T->L);
    1176          12 :   for (i = 1;;)
    1177             :   {
    1178          54 :     pari_sp av = avma;
    1179          54 :     if (DEBUGLEVEL) err_printf("nf_LLL_cmbf: checking factor %ld\n", i);
    1180          54 :     y = chk_factors_get(D.lt, famod, gel(piv,i), Tpk, pk);
    1181             : 
    1182          54 :     if (! (y = nf_pol_lift(y, bound, T->L)) ) return NULL;
    1183          54 :     y = gerepilecopy(av, y);
    1184             :     /* y is the candidate factor */
    1185          54 :     pol = RgXQX_divrem(D.C2ltpol, y, nfT, ONLY_DIVIDES);
    1186          54 :     if (!pol) return NULL;
    1187             : 
    1188          54 :     if (D.C2lt) y = RgX_int_normalize(y);
    1189          54 :     gel(list,i) = y;
    1190          54 :     if (++i >= r) break;
    1191             : 
    1192          42 :     update_target(&D, pol);
    1193          42 :   }
    1194          12 :   gel(list,i) = RgX_int_normalize(pol); return list;
    1195             : }
    1196             : 
    1197             : static GEN
    1198       11688 : nf_to_Zq(GEN x, GEN T, GEN pk, GEN pks2, GEN proj)
    1199             : {
    1200             :   GEN y;
    1201       11688 :   if (typ(x) != t_COL) return centermodii(x, pk, pks2);
    1202        1812 :   if (!T)
    1203             :   {
    1204        1752 :     y = ZV_dotproduct(proj, x);
    1205        1752 :     return centermodii(y, pk, pks2);
    1206             :   }
    1207          60 :   y = ZM_ZC_mul(proj, x);
    1208          60 :   y = RgV_to_RgX(y, varn(T));
    1209          60 :   return FpX_center(FpX_rem(y, T, pk), pk, pks2);
    1210             : }
    1211             : 
    1212             : /* Assume P in nfX form, lc(P) != 0 mod p. Reduce P to Zp[X]/(T) mod p^a */
    1213             : static GEN
    1214        1014 : ZqX(GEN P, GEN pk, GEN T, GEN proj)
    1215             : {
    1216        1014 :   long i, l = lg(P);
    1217        1014 :   GEN z, pks2 = shifti(pk,-1);
    1218             : 
    1219        1014 :   z = cgetg(l,t_POL); z[1] = P[1];
    1220        1014 :   for (i=2; i<l; i++) gel(z,i) = nf_to_Zq(gel(P,i),T,pk,pks2,proj);
    1221        1014 :   return normalizepol_lg(z, l);
    1222             : }
    1223             : 
    1224             : static GEN
    1225        1014 : ZqX_normalize(GEN P, GEN lt, nflift_t *L)
    1226             : {
    1227        1014 :   GEN R = lt? RgX_Rg_mul(P, Fp_inv(lt, L->pk)): P;
    1228        1014 :   return ZqX(R, L->pk, L->Tpk, L->ZqProj);
    1229             : }
    1230             : 
    1231             : /* k allowing to reconstruct x, |x|^2 < C, from x mod pr^k */
    1232             : /* return log [  2sqrt(C/d) * ( (3/2)sqrt(gamma) )^(d-1) ] ^d / log N(pr)
    1233             :  * cf. Belabas relative van Hoeij algorithm, lemma 3.12 */
    1234             : static double
    1235        1014 : bestlift_bound(GEN C, long d, double alpha, GEN Npr)
    1236             : {
    1237        1014 :   const double y = 1 / (alpha - 0.25); /* = 2 if alpha = 3/4 */
    1238             :   double t;
    1239        1014 :   C = gtofp(C,DEFAULTPREC);
    1240             :   /* (1/2)log (4C/d) + (d-1)(log 3/2 sqrt(gamma)) */
    1241        1014 :   t = rtodbl(mplog(gmul2n(divru(C,d), 2))) * 0.5 + (d-1) * log(1.5 * sqrt(y));
    1242        1014 :   return ceil((t * d) / log(gtodouble(Npr)));
    1243             : }
    1244             : 
    1245             : static GEN
    1246        1386 : get_R(GEN M)
    1247             : {
    1248             :   GEN R;
    1249        1386 :   long i, l, prec = nbits2prec( gexpo(M) + 64 );
    1250             : 
    1251             :   for(;;)
    1252             :   {
    1253        1386 :     R = gaussred_from_QR(M, prec);
    1254        1386 :     if (R) break;
    1255           0 :     prec = precdbl(prec);
    1256           0 :   }
    1257        1386 :   l = lg(R);
    1258        1386 :   for (i=1; i<l; i++) gcoeff(R,i,i) = gen_1;
    1259        1386 :   return R;
    1260             : }
    1261             : 
    1262             : static void
    1263        1368 : init_proj(nflift_t *L, GEN nfT)
    1264             : {
    1265        1368 :   if (L->Tp)
    1266             :   {
    1267         102 :     GEN coTp = FpX_div(FpX_red(nfT, L->p), L->Tp,  L->p); /* Tp's cofactor */
    1268             :     GEN z, proj;
    1269         102 :     z = ZpX_liftfact(nfT, mkvec2(L->Tp, coTp), L->pk, L->p, L->k);
    1270         102 :     L->Tpk = gel(z,1);
    1271         102 :     proj = QXQV_to_FpM(L->topow, L->Tpk, L->pk);
    1272         102 :     if (L->topowden)
    1273         102 :       proj = FpM_red(ZM_Z_mul(proj, Fp_inv(L->topowden, L->pk)), L->pk);
    1274         102 :     L->ZqProj = proj;
    1275             :   }
    1276             :   else
    1277             :   {
    1278        1266 :     L->Tpk = NULL;
    1279        1266 :     L->ZqProj = dim1proj(L->prkHNF);
    1280             :   }
    1281        1368 : }
    1282             : 
    1283             : /* Square of the radius of largest ball inscript in PRK's fundamental domain,
    1284             :  *   whose orthogonalized vector's norms are the Bi
    1285             :  * Rmax ^2 =  min 1/4T_i where T_i = sum ( s_ij^2 / B_j) */
    1286             : static GEN
    1287        1386 : max_radius(GEN PRK, GEN B)
    1288             : {
    1289        1386 :   GEN S, smax = gen_0;
    1290        1386 :   pari_sp av = avma;
    1291        1386 :   long i, j, d = lg(PRK)-1;
    1292             : 
    1293        1386 :   S = RgM_inv( get_R(PRK) ); if (!S) pari_err_PREC("max_radius");
    1294        9420 :   for (i=1; i<=d; i++)
    1295             :   {
    1296        8034 :     GEN s = gen_0;
    1297      123264 :     for (j=1; j<=d; j++)
    1298      115230 :       s = mpadd(s, mpdiv( mpsqr(gcoeff(S,i,j)), gel(B,j)));
    1299        8034 :     if (mpcmp(s, smax) > 0) smax = s;
    1300             :   }
    1301        1386 :   return gerepileupto(av, ginv(gmul2n(smax, 2)));
    1302             : }
    1303             : 
    1304             : static void
    1305        1368 : bestlift_init(long a, GEN nf, GEN pr, GEN C, nflift_t *L)
    1306             : {
    1307        1368 :   const double alpha = 0.99; /* LLL parameter */
    1308        1368 :   const long d = nf_get_degree(nf);
    1309        1368 :   pari_sp av = avma, av2;
    1310             :   GEN prk, PRK, B, GSmin, pk;
    1311             :   pari_timer ti;
    1312             : 
    1313        1368 :   timer_start(&ti);
    1314        1368 :   if (!a) a = (long)bestlift_bound(C, d, alpha, pr_norm(pr));
    1315             : 
    1316          18 :   for (;; avma = av, a += (a==1)? 1: (a>>1)) /* roughly a *= 1.5 */
    1317             :   {
    1318        1386 :     if (DEBUGLEVEL>2) err_printf("exponent %ld\n",a);
    1319        1386 :     prk = idealpows(nf, pr, a);
    1320        1386 :     av2 = avma;
    1321        1386 :     pk = gcoeff(prk,1,1);
    1322        1386 :     PRK = ZM_lll_norms(prk, alpha, LLL_INPLACE, &B);
    1323        1386 :     GSmin = max_radius(PRK, B);
    1324        1386 :     if (gcmp(GSmin, C) >= 0) break;
    1325          18 :   }
    1326        1368 :   gerepileall(av2, 2, &PRK, &GSmin);
    1327        1368 :   if (DEBUGLEVEL>2)
    1328           0 :     err_printf("for this exponent, GSmin = %Ps\nTime reduction: %ld\n",
    1329             :       GSmin, timer_delay(&ti));
    1330        1368 :   L->k = a;
    1331        1368 :   L->den = L->pk = pk;
    1332        1368 :   L->prk = PRK;
    1333        1368 :   L->iprk = ZM_inv(PRK, pk);
    1334        1368 :   L->GSmin= GSmin;
    1335        1368 :   L->prkHNF = prk;
    1336        1368 :   init_proj(L, nf_get_pol(nf));
    1337        1368 : }
    1338             : 
    1339             : /* Let X = Tra * M_L, Y = bestlift(X) return V s.t Y = X - PRK V
    1340             :  * and set *eT2 = gexpo(Y)  [cf nf_bestlift, but memory efficient] */
    1341             : static GEN
    1342         210 : get_V(GEN Tra, GEN M_L, GEN PRK, GEN PRKinv, GEN pk, long *eT2)
    1343             : {
    1344         210 :   long i, e = 0, l = lg(M_L);
    1345         210 :   GEN V = cgetg(l, t_MAT);
    1346         210 :   *eT2 = 0;
    1347        2526 :   for (i = 1; i < l; i++)
    1348             :   { /* cf nf_bestlift(Tra * c) */
    1349        2316 :     pari_sp av = avma, av2;
    1350        2316 :     GEN v, T2 = ZM_ZC_mul(Tra, gel(M_L,i));
    1351             : 
    1352        2316 :     v = gdivround(ZM_ZC_mul(PRKinv, T2), pk); /* small */
    1353        2316 :     av2 = avma;
    1354        2316 :     T2 = ZC_sub(T2, ZM_ZC_mul(PRK, v));
    1355        2316 :     e = gexpo(T2); if (e > *eT2) *eT2 = e;
    1356        2316 :     avma = av2;
    1357        2316 :     gel(V,i) = gerepileupto(av, v); /* small */
    1358             :   }
    1359         210 :   return V;
    1360             : }
    1361             : 
    1362             : static GEN
    1363          42 : nf_LLL_cmbf(nfcmbf_t *T, long rec)
    1364             : {
    1365          42 :   const double BitPerFactor = 0.4; /* nb bits / modular factor */
    1366          42 :   nflift_t *L = T->L;
    1367          42 :   GEN famod = T->fact, ZC = T->ZC, Br = T->Br, P = T->pol, dn = T->L->dn;
    1368          42 :   long dnf = nf_get_degree(T->nf), dP = degpol(P);
    1369             :   long i, C, tmax, n0;
    1370             :   GEN lP, Bnorm, Tra, T2, TT, CM_L, m, list, ZERO, Btra;
    1371             :   double Bhigh;
    1372             :   pari_sp av, av2;
    1373          42 :   long ti_LLL = 0, ti_CF = 0;
    1374             :   pari_timer ti2, TI;
    1375             : 
    1376          42 :   lP = absi(leading_coeff(P));
    1377          42 :   if (is_pm1(lP)) lP = NULL;
    1378             : 
    1379          42 :   n0 = lg(famod) - 1;
    1380             :  /* Lattice: (S PRK), small vector (vS vP). To find k bound for the image,
    1381             :   * write S = S1 q + S0, P = P1 q + P0
    1382             :   * |S1 vS + P1 vP|^2 <= Bhigh for all (vS,vP) assoc. to true factors */
    1383          42 :   Btra = mulrr(ZC, mulur(dP*dP, normlp(Br, 2, dnf)));
    1384          42 :   Bhigh = get_Bhigh(n0, dnf);
    1385          42 :   C = (long)ceil(sqrt(Bhigh/n0)) + 1; /* C^2 n0 ~ Bhigh */
    1386          42 :   Bnorm = dbltor( n0 * C * C + Bhigh );
    1387          42 :   ZERO = zeromat(n0, dnf);
    1388             : 
    1389          42 :   av = avma;
    1390          42 :   TT = cgetg(n0+1, t_VEC);
    1391          42 :   Tra  = cgetg(n0+1, t_MAT);
    1392          42 :   for (i=1; i<=n0; i++) TT[i] = 0;
    1393          42 :   CM_L = scalarmat_s(C, n0);
    1394             :   /* tmax = current number of traces used (and computed so far) */
    1395         156 :   for(tmax = 0;; tmax++)
    1396             :   {
    1397         156 :     long a, b, bmin, bgood, delta, tnew = tmax + 1, r = lg(CM_L)-1;
    1398             :     GEN M_L, q, CM_Lp, oldCM_L, S1, P1, VV;
    1399         156 :     int first = 1;
    1400             : 
    1401             :     /* bound for f . S_k(genuine factor) = ZC * bound for T_2(S_tnew) */
    1402         156 :     Btra = mulrr(ZC, mulur(dP*dP, normlp(Br, 2*tnew, dnf)));
    1403         156 :     bmin = logint(ceil_safe(sqrtr(Btra)), gen_2) + 1;
    1404         156 :     if (DEBUGLEVEL>2)
    1405           0 :       err_printf("\nLLL_cmbf: %ld potential factors (tmax = %ld, bmin = %ld)\n",
    1406             :                  r, tmax, bmin);
    1407             : 
    1408             :     /* compute Newton sums (possibly relifting first) */
    1409         156 :     if (gcmp(L->GSmin, Btra) < 0)
    1410             :     {
    1411             :       GEN polred;
    1412             : 
    1413           0 :       bestlift_init((L->k)<<1, T->nf, T->pr, Btra, L);
    1414           0 :       polred = ZqX_normalize(T->polbase, lP, L);
    1415           0 :       famod = ZqX_liftfact(polred, famod, L->Tpk, L->pk, L->p, L->k);
    1416           0 :       for (i=1; i<=n0; i++) TT[i] = 0;
    1417             :     }
    1418        2940 :     for (i=1; i<=n0; i++)
    1419             :     {
    1420        2784 :       GEN h, lPpow = lP? powiu(lP, tnew): NULL;
    1421        2784 :       GEN z = polsym_gen(gel(famod,i), gel(TT,i), tnew, L->Tpk, L->pk);
    1422        2784 :       gel(TT,i) = z;
    1423        2784 :       h = gel(z,tnew+1);
    1424             :       /* make Newton sums integral */
    1425        2784 :       lPpow = mul_content(lPpow, dn);
    1426        2784 :       if (lPpow)
    1427           0 :         h = (typ(h) == t_INT)? Fp_mul(h, lPpow, L->pk): FpX_Fp_mul(h, lPpow, L->pk);
    1428        2784 :       gel(Tra,i) = nf_bestlift(h, NULL, L); /* S_tnew(famod) */
    1429             :     }
    1430             : 
    1431             :     /* compute truncation parameter */
    1432         156 :     if (DEBUGLEVEL>2) { timer_start(&ti2); timer_start(&TI); }
    1433         156 :     oldCM_L = CM_L;
    1434         156 :     av2 = avma;
    1435         156 :     b = delta = 0; /* -Wall */
    1436             : AGAIN:
    1437         210 :     M_L = Q_div_to_int(CM_L, utoipos(C));
    1438         210 :     VV = get_V(Tra, M_L, L->prk, L->iprk, L->pk, &a);
    1439         210 :     if (first)
    1440             :     { /* initialize lattice, using few p-adic digits for traces */
    1441         156 :       bgood = (long)(a - maxss(32, (long)(BitPerFactor * r)));
    1442         156 :       b = maxss(bmin, bgood);
    1443         156 :       delta = a - b;
    1444             :     }
    1445             :     else
    1446             :     { /* add more p-adic digits and continue reduction */
    1447          54 :       if (a < b) b = a;
    1448          54 :       b = maxss(b-delta, bmin);
    1449          54 :       if (b - delta/2 < bmin) b = bmin; /* near there. Go all the way */
    1450             :     }
    1451             : 
    1452             :     /* restart with truncated entries */
    1453         210 :     q = int2n(b);
    1454         210 :     P1 = gdivround(L->prk, q);
    1455         210 :     S1 = gdivround(Tra, q);
    1456         210 :     T2 = ZM_sub(ZM_mul(S1, M_L), ZM_mul(P1, VV));
    1457         210 :     m = vconcat( CM_L, T2 );
    1458         210 :     if (first)
    1459             :     {
    1460         156 :       first = 0;
    1461         156 :       m = shallowconcat( m, vconcat(ZERO, P1) );
    1462             :       /*     [ C M_L   0  ]
    1463             :        * m = [            ]   square matrix
    1464             :        *     [  T2'   PRK ]   T2' = Tra * M_L  truncated
    1465             :        */
    1466             :     }
    1467         210 :     CM_L = LLL_check_progress(Bnorm, n0, m, b == bmin, /*dbg:*/ &ti_LLL);
    1468         210 :     if (DEBUGLEVEL>2)
    1469           0 :       err_printf("LLL_cmbf: (a,b) =%4ld,%4ld; r =%3ld -->%3ld, time = %ld\n",
    1470           0 :                  a,b, lg(m)-1, CM_L? lg(CM_L)-1: 1, timer_delay(&TI));
    1471         252 :     if (!CM_L) { list = mkcol(RgX_int_normalize(P)); break; }
    1472         180 :     if (b > bmin)
    1473             :     {
    1474          54 :       CM_L = gerepilecopy(av2, CM_L);
    1475          54 :       goto AGAIN;
    1476             :     }
    1477         126 :     if (DEBUGLEVEL>2) timer_printf(&ti2, "for this trace");
    1478             : 
    1479         126 :     i = lg(CM_L) - 1;
    1480         126 :     if (i == r && ZM_equal(CM_L, oldCM_L))
    1481             :     {
    1482          60 :       CM_L = oldCM_L;
    1483          60 :       avma = av2; continue;
    1484             :     }
    1485             : 
    1486          66 :     CM_Lp = FpM_image(CM_L, utoipos(27449)); /* inexpensive test */
    1487          66 :     if (lg(CM_Lp) != lg(CM_L))
    1488             :     {
    1489           0 :       if (DEBUGLEVEL>2) err_printf("LLL_cmbf: rank decrease\n");
    1490           0 :       CM_L = ZM_hnf(CM_L);
    1491             :     }
    1492             : 
    1493          66 :     if (i <= r && i*rec < n0)
    1494             :     {
    1495             :       pari_timer ti;
    1496          24 :       if (DEBUGLEVEL>2) timer_start(&ti);
    1497          24 :       list = nf_chk_factors(T, P, Q_div_to_int(CM_L,utoipos(C)), famod, L->pk);
    1498          24 :       if (DEBUGLEVEL>2) ti_CF += timer_delay(&ti);
    1499          24 :       if (list) break;
    1500             :     }
    1501          54 :     CM_L = gerepilecopy(av2, CM_L);
    1502          54 :     if (gc_needed(av,1))
    1503             :     {
    1504           0 :       if(DEBUGMEM>1) pari_warn(warnmem,"nf_LLL_cmbf");
    1505           0 :       gerepileall(av, L->Tpk? 9: 8,
    1506             :                       &CM_L,&TT,&Tra,&famod,&L->pk,&L->GSmin,&L->prk,&L->iprk,&L->Tpk);
    1507             :     }
    1508         114 :   }
    1509          42 :   if (DEBUGLEVEL>2)
    1510           0 :     err_printf("* Time LLL: %ld\n* Time Check Factor: %ld\n",ti_LLL,ti_CF);
    1511          42 :   return list;
    1512             : }
    1513             : 
    1514             : static GEN
    1515         324 : nf_combine_factors(nfcmbf_t *T, GEN polred, long klim)
    1516             : {
    1517         324 :   nflift_t *L = T->L;
    1518             :   GEN res;
    1519             :   long maxK;
    1520             :   int done;
    1521             :   pari_timer ti;
    1522             : 
    1523         324 :   if (DEBUGLEVEL>2) timer_start(&ti);
    1524         324 :   T->fact = ZqX_liftfact(polred, T->fact, L->Tpk, L->pk, L->p, L->k);
    1525         324 :   if (DEBUGLEVEL>2) timer_printf(&ti, "Hensel lift");
    1526         324 :   res = nfcmbf(T, klim, &maxK, &done);
    1527         324 :   if (DEBUGLEVEL>2) timer_printf(&ti, "Naive recombination");
    1528         324 :   if (!done)
    1529             :   {
    1530          42 :     long l = lg(res)-1;
    1531             :     GEN v;
    1532          42 :     if (l > 1)
    1533             :     {
    1534           6 :       T->pol = gel(res,l);
    1535           6 :       T->polbase = RgX_to_nfX(T->nf, T->pol);
    1536             :     }
    1537          42 :     v = nf_LLL_cmbf(T, maxK);
    1538             :     /* remove last elt, possibly unfactored. Add all new ones. */
    1539          42 :     setlg(res, l); res = shallowconcat(res, v);
    1540             :   }
    1541         324 :   return res;
    1542             : }
    1543             : 
    1544             : static GEN
    1545         690 : nf_DDF_roots(GEN pol, GEN polred, GEN nfpol, long fl, nflift_t *L)
    1546             : {
    1547             :   GEN z, Cltx_r, ltdn;
    1548             :   long i, m;
    1549             :   div_data D;
    1550             : 
    1551         690 :   init_div_data(&D, pol, L);
    1552         690 :   ltdn = mul_content(D.lt, L->dn);
    1553         690 :   z = ZqX_roots(polred, L->Tpk, L->p, L->k);
    1554         690 :   Cltx_r = deg1pol_shallow(D.Clt? D.Clt: gen_1, NULL, varn(pol));
    1555        2682 :   for (m=1,i=1; i<lg(z); i++)
    1556             :   {
    1557        1992 :     GEN q, r = gel(z,i);
    1558             :     pari_sp av;
    1559             :     /* lt*dn*topowden * r = Clt * r */
    1560        1992 :     r = nf_bestlift_to_pol(ltdn? gmul(ltdn,r): r, NULL, L);
    1561        1992 :     av = avma;
    1562        1992 :     gel(Cltx_r,2) = gneg(r); /* check P(r) == 0 */
    1563        1992 :     q = RgXQX_divrem(D.C2ltpol, Cltx_r, nfpol, ONLY_DIVIDES); /* integral */
    1564        1992 :     avma = av;
    1565             :     /* don't go on with q, usually much larger that C2ltpol */
    1566        1992 :     if (q) {
    1567        1860 :       if (D.Clt) r = gdiv(r, D.Clt);
    1568        1860 :       gel(z,m++) = r;
    1569             :     }
    1570         132 :     else if (fl == ROOTS_SPLIT) return cgetg(1, t_VEC);
    1571             :   }
    1572         690 :   z[0] = evaltyp(t_VEC) | evallg(m);
    1573         690 :   return z;
    1574             : }
    1575             : 
    1576             : /* returns a factor of T in Fp of degree <= maxf, NULL if none exist */
    1577             : static GEN
    1578       23340 : get_good_factor(GEN T, ulong p, long maxf)
    1579             : {
    1580       23340 :   pari_sp av = avma;
    1581       23340 :   GEN r, list = gel(Flx_factor(ZX_to_Flx(T,p),p), 1);
    1582       23340 :   if (maxf == 1)
    1583             :   { /* deg.1 factors are best */
    1584       22404 :     r = gel(list,1);
    1585       22404 :     if (degpol(r) == 1) return r;
    1586             :   }
    1587             :   else
    1588             :   { /* otherwise, pick factor of largish degree */
    1589             :     long i, dr;
    1590        1740 :     for (i = lg(list)-1; i > 0; i--)
    1591             :     {
    1592        1266 :       r = gel(list,i); dr = degpol(r);
    1593        1266 :       if (dr <= maxf) return r;
    1594             :     }
    1595             :   }
    1596       16986 :   avma = av; return NULL; /* failure */
    1597             : }
    1598             : 
    1599             : /* Optimization problem: factorization of polynomials over large Fq is slow,
    1600             :  * BUT bestlift correspondingly faster.
    1601             :  * Return maximal residue degree to be considered when picking a prime ideal */
    1602             : static long
    1603        1530 : get_maxf(long nfdeg)
    1604             : {
    1605        1530 :   long maxf = 1;
    1606        1530 :   if      (nfdeg >= 45) maxf =32;
    1607        1518 :   else if (nfdeg >= 30) maxf =16;
    1608        1506 :   else if (nfdeg >= 15) maxf = 8;
    1609        1530 :   return maxf;
    1610             : }
    1611             : 
    1612             : /* Select a prime ideal pr over which to factor polbase.
    1613             :  * Return the number of factors (or roots, according to flag fl) mod pr,
    1614             :  * Input:
    1615             :  *   ct: number of attempts to find best
    1616             :  * Set:
    1617             :  *   lt: leading term of polbase (t_INT or NULL [ for 1 ])
    1618             :  *   pr: a suitable maximal ideal
    1619             :  *   Fa: factors found mod pr
    1620             :  *   Tp: polynomial defining Fq/Fp */
    1621             : static long
    1622        1176 : nf_pick_prime(long ct, GEN nf, GEN polbase, long fl,
    1623             :               GEN *lt, GEN *pr, GEN *Tp)
    1624             : {
    1625        1176 :   GEN nfpol = nf_get_pol(nf), bad = mulii(nf_get_disc(nf), nf_get_index(nf));
    1626        1176 :   long maxf, nfdeg = degpol(nfpol), dpol = degpol(polbase), nbf = 0;
    1627             :   ulong pp;
    1628             :   forprime_t S;
    1629             :   pari_timer ti_pr;
    1630             : 
    1631        1176 :   if (DEBUGLEVEL>3) timer_start(&ti_pr);
    1632        1176 :   *lt  = leading_coeff(polbase); /* t_INT */
    1633        1176 :   if (gequal1(*lt)) *lt = NULL;
    1634        1176 :   *pr = NULL;
    1635        1176 :   *Tp = NULL;
    1636             : 
    1637        1176 :   maxf = get_maxf(nfdeg);
    1638        1176 :   (void)u_forprime_init(&S, 2, ULONG_MAX);
    1639             :   /* select pr such that pol has the smallest number of factors, ct attempts */
    1640        1176 :   while ((pp = u_forprime_next(&S)))
    1641             :   {
    1642             :     GEN aT, apr, ap, amodpr, red, r;
    1643             :     long anbf;
    1644       23160 :     ulong ltp = 0;
    1645       23160 :     pari_sp av2 = avma;
    1646             : 
    1647             :     /* first step : select prime of high inertia degree */
    1648       40950 :     if (! umodiu(bad,pp)) continue;
    1649       20652 :     if (*lt) { ltp = umodiu(*lt, pp); if (!ltp) continue; }
    1650       20478 :     r = get_good_factor(nfpol, pp, maxf);
    1651       20478 :     if (!r) continue;
    1652             : 
    1653        6000 :     ap = utoipos(pp);
    1654        6000 :     apr = idealprimedec_kummer(nf, Flx_to_ZX(r), 1, ap);
    1655        6000 :     amodpr = zk_to_Fq_init(nf,&apr,&aT,&ap);
    1656             : 
    1657             :     /* second step : evaluate factorisation mod apr */
    1658        6000 :     red = nfX_to_FqX(polbase, nf, amodpr);
    1659        6000 :     if (!aT)
    1660             :     { /* degree 1 */
    1661        5586 :       red = ZX_to_Flx(red, pp);
    1662        5586 :       if (ltp) red = Flx_normalize(red, pp);
    1663        5586 :       if (!Flx_is_squarefree(red, pp)) { avma = av2; continue; }
    1664        5028 :       anbf = fl == FACTORS? Flx_nbfact(red, pp): Flx_nbroots(red, pp);
    1665             :     }
    1666             :     else
    1667             :     {
    1668         414 :       if (ltp) red = FqX_normalize(red, aT,ap);
    1669         414 :       if (!FqX_is_squarefree(red,aT,ap)) { avma = av2; continue; }
    1670         564 :       anbf = fl == FACTORS? FqX_nbfact(red, aT, ap)
    1671         564 :                           : FqX_nbroots(red, aT, ap);
    1672             :     }
    1673        5532 :     if (fl == ROOTS_SPLIT && anbf < dpol) return anbf;
    1674        5334 :     if (anbf <= 1)
    1675             :     {
    1676         390 :       if (fl == FACTORS) return anbf; /* irreducible */
    1677         354 :       if (!anbf) return 0; /* no root */
    1678             :     }
    1679        5208 :     if (DEBUGLEVEL>3)
    1680           0 :       err_printf("%3ld %s at prime\n  %Ps\nTime: %ld\n",
    1681             :                  anbf, fl == FACTORS?"factors": "roots", apr, timer_delay(&ti_pr));
    1682             : 
    1683        5208 :     if (!nbf || anbf < nbf
    1684        4014 :              || (anbf == nbf && pr_get_f(apr) > pr_get_f(*pr)))
    1685             :     {
    1686        1206 :       nbf = anbf;
    1687        1206 :       *pr = apr;
    1688        1206 :       *Tp = aT;
    1689             :     }
    1690        4002 :     else avma = av2;
    1691        5208 :     if (--ct <= 0) break;
    1692             :   }
    1693        1014 :   if (!nbf) pari_err_OVERFLOW("nf_pick_prime [ran out of primes]");
    1694        1014 :   return nbf;
    1695             : }
    1696             : 
    1697             : /* assume lt(T) is a t_INT and T square free */
    1698             : static GEN
    1699         192 : nfsqff_trager(GEN u, GEN T, GEN dent)
    1700             : {
    1701         192 :   long k = 0, i, lx;
    1702         192 :   GEN U, P, x0, mx0, fa, n = ZX_ZXY_rnfequation(T, u, &k);
    1703             :   int tmonic;
    1704         192 :   if (DEBUGLEVEL>4) err_printf("nfsqff_trager: choosing k = %ld\n",k);
    1705             : 
    1706             :   /* n guaranteed to be squarefree */
    1707         192 :   fa = ZX_DDF(Q_primpart(n)); lx = lg(fa);
    1708         192 :   if (lx == 2) return mkcol(u);
    1709             : 
    1710         126 :   tmonic = is_pm1(leading_coeff(T));
    1711         126 :   P = cgetg(lx,t_COL);
    1712         126 :   x0 = deg1pol_shallow(stoi(-k), gen_0, varn(T));
    1713         126 :   mx0 = deg1pol_shallow(stoi(k), gen_0, varn(T));
    1714         126 :   U = RgXQX_translate(u, mx0, T);
    1715         126 :   if (!tmonic) U = Q_primpart(U);
    1716         540 :   for (i=lx-1; i>0; i--)
    1717             :   {
    1718         414 :     GEN f = gel(fa,i), F = nfgcd(U, f, T, dent);
    1719         414 :     F = RgXQX_translate(F, x0, T);
    1720             :     /* F = gcd(f, u(t - x0)) [t + x0] = gcd(f(t + x0), u), more efficient */
    1721         414 :     if (typ(F) != t_POL || degpol(F) == 0)
    1722           0 :       pari_err_IRREDPOL("factornf [modulus]",T);
    1723         414 :     gel(P,i) = QXQX_normalize(F, T);
    1724             :   }
    1725         126 :   gen_sort_inplace(P, (void*)&cmp_RgX, &gen_cmp_RgX, NULL);
    1726         126 :   return P;
    1727             : }
    1728             : 
    1729             : /* Factor polynomial a on the number field defined by polynomial T, using
    1730             :  * Trager's trick */
    1731             : GEN
    1732          96 : polfnf(GEN a, GEN T)
    1733             : {
    1734          96 :   GEN rep = cgetg(3, t_MAT), A, B, y, dent, bad;
    1735             :   long dA;
    1736             :   int tmonic;
    1737             : 
    1738          96 :   if (typ(a)!=t_POL) pari_err_TYPE("polfnf",a);
    1739          96 :   if (typ(T)!=t_POL) pari_err_TYPE("polfnf",T);
    1740          96 :   T = Q_primpart(T); tmonic = is_pm1(leading_coeff(T));
    1741          96 :   RgX_check_ZX(T,"polfnf");
    1742          96 :   A = Q_primpart( QXQX_normalize(RgX_nffix("polfnf",T,a,1), T) );
    1743          96 :   dA = degpol(A);
    1744          96 :   if (dA <= 0)
    1745             :   {
    1746           0 :     avma = (pari_sp)(rep + 3);
    1747           0 :     return (dA == 0)? trivial_fact(): zerofact(varn(A));
    1748             :   }
    1749          96 :   bad = dent = ZX_disc(T);
    1750          96 :   if (tmonic) dent = indexpartial(T, dent);
    1751          96 :   (void)nfgcd_all(A,RgX_deriv(A), T, dent, &B);
    1752          96 :   if (degpol(B) != dA) B = Q_primpart( QXQX_normalize(B, T) );
    1753          96 :   ensure_lt_INT(B);
    1754          96 :   y = nfsqff_trager(B, T, dent);
    1755          96 :   fact_from_sqff(rep, A, B, y, T, bad);
    1756          96 :   return sort_factor_pol(rep, cmp_RgX);
    1757             : }
    1758             : 
    1759             : static int
    1760        2478 : nfsqff_use_Trager(long n, long dpol)
    1761             : {
    1762        2478 :   return dpol*3<n;
    1763             : }
    1764             : 
    1765             : /* return the factorization of the square-free polynomial pol. Not memory-clean
    1766             :    The coeffs of pol are in Z_nf and its leading term is a rational integer.
    1767             :    deg(pol) > 0, deg(nfpol) > 1
    1768             :    fl is either FACTORS (return factors), or ROOTS / ROOTS_SPLIT (return roots):
    1769             :      - ROOTS, return only the roots of x in nf
    1770             :      - ROOTS_SPLIT, as ROOTS if pol splits, [] otherwise
    1771             :    den is usually 1, otherwise nf.zk is doubtful, and den bounds the
    1772             :    denominator of an arbitrary element of Z_nf on nf.zk */
    1773             : static GEN
    1774        1608 : nfsqff(GEN nf, GEN pol, long fl, GEN den)
    1775             : {
    1776        1608 :   long n, nbf, dpol = degpol(pol);
    1777             :   GEN pr, C0, polbase;
    1778        1608 :   GEN N2, res, polred, lt, nfpol = typ(nf)==t_POL?nf:nf_get_pol(nf);
    1779             :   nfcmbf_t T;
    1780             :   nflift_t L;
    1781             :   pari_timer ti, ti_tot;
    1782             : 
    1783        1608 :   if (DEBUGLEVEL>2) { timer_start(&ti); timer_start(&ti_tot); }
    1784        1608 :   n = degpol(nfpol);
    1785             :   /* deg = 1 => irreducible */
    1786        1608 :   if (dpol == 1) {
    1787         336 :     if (fl == FACTORS) return mkvec(QXQX_normalize(pol, nfpol));
    1788         330 :     return mkvec(gneg(gdiv(gel(pol,2),gel(pol,3))));
    1789             :   }
    1790        1272 :   if (typ(nf)==t_POL || nfsqff_use_Trager(n,dpol))
    1791             :   {
    1792             :     GEN z;
    1793          96 :     if (DEBUGLEVEL>2) err_printf("Using Trager's method\n");
    1794          96 :     if (typ(nf) != t_POL) den =  mulii(den, nf_get_index(nf));
    1795          96 :     z = nfsqff_trager(Q_primpart(pol), nfpol, den);
    1796          96 :     if (fl != FACTORS) {
    1797          78 :       long i, l = lg(z);
    1798         246 :       for (i = 1; i < l; i++)
    1799             :       {
    1800         186 :         GEN LT, t = gel(z,i); if (degpol(t) > 1) break;
    1801         168 :         LT = gel(t,3);
    1802         168 :         if (typ(LT) == t_POL) LT = gel(LT,2); /* constant */
    1803         168 :         gel(z,i) = gdiv(gel(t,2), negi(LT));
    1804             :       }
    1805          78 :       setlg(z, i);
    1806          78 :       if (fl == ROOTS_SPLIT && i != l) return cgetg(1,t_VEC);
    1807             :     }
    1808          96 :     return z;
    1809             :   }
    1810             : 
    1811        1176 :   polbase = RgX_to_nfX(nf, pol);
    1812        1176 :   nbf = nf_pick_prime(5, nf, polbase, fl, &lt, &pr, &L.Tp);
    1813        1176 :   if (fl == ROOTS_SPLIT && nbf < dpol) return cgetg(1,t_VEC);
    1814        1140 :   if (nbf <= 1)
    1815             :   {
    1816         198 :     if (fl == FACTORS) return mkvec(QXQX_normalize(pol, nfpol)); /* irred. */
    1817         162 :     if (!nbf) return cgetg(1,t_VEC); /* no root */
    1818             :   }
    1819             : 
    1820        1014 :   if (DEBUGLEVEL>2) {
    1821           0 :     timer_printf(&ti, "choice of a prime ideal");
    1822           0 :     err_printf("Prime ideal chosen: %Ps\n", pr);
    1823             :   }
    1824        1014 :   L.tozk = nf_get_invzk(nf);
    1825        1014 :   L.topow= nf_get_zkprimpart(nf);
    1826        1014 :   L.topowden = nf_get_zkden(nf);
    1827        1014 :   if (is_pm1(den)) den = NULL;
    1828        1014 :   L.dn = den;
    1829        1014 :   T.ZC = L2_bound(nf, den);
    1830        1014 :   T.Br = nf_root_bounds(pol, nf); if (lt) T.Br = gmul(T.Br, lt);
    1831             : 
    1832             :   /* C0 = bound for T_2(Q_i), Q | P */
    1833        1014 :   if (fl != FACTORS) C0 = normlp(T.Br, 2, n);
    1834         324 :   else               C0 = nf_factor_bound(nf, polbase);
    1835        1014 :   T.bound = mulrr(T.ZC, C0); /* bound for |Q_i|^2 in Z^n on chosen Z-basis */
    1836             : 
    1837        1014 :   N2 = mulur(dpol*dpol, normlp(T.Br, 4, n)); /* bound for T_2(lt * S_2) */
    1838        1014 :   T.BS_2 = mulrr(T.ZC, N2); /* bound for |S_2|^2 on chosen Z-basis */
    1839             : 
    1840        1014 :   if (DEBUGLEVEL>2) {
    1841           0 :     timer_printf(&ti, "bound computation");
    1842           0 :     err_printf("  1) T_2 bound for %s: %Ps\n",
    1843             :                fl == FACTORS?"factor": "root", C0);
    1844           0 :     err_printf("  2) Conversion from T_2 --> | |^2 bound : %Ps\n", T.ZC);
    1845           0 :     err_printf("  3) Final bound: %Ps\n", T.bound);
    1846             :   }
    1847             : 
    1848        1014 :   L.p = pr_get_p(pr);
    1849        1014 :   if (L.Tp && degpol(L.Tp) == 1) L.Tp = NULL;
    1850        1014 :   bestlift_init(0, nf, pr, T.bound, &L);
    1851        1014 :   if (DEBUGLEVEL>2) timer_start(&ti);
    1852        1014 :   polred = ZqX_normalize(polbase, lt, &L); /* monic */
    1853             : 
    1854        1014 :   if (fl != FACTORS) {
    1855         690 :     GEN z = nf_DDF_roots(pol, polred, nfpol, fl, &L);
    1856         690 :     if (lg(z) == 1) return cgetg(1, t_VEC);
    1857         678 :     return z;
    1858             :   }
    1859             : 
    1860         324 :   T.fact = gel(FqX_factor(polred, L.Tp, L.p), 1);
    1861         324 :   if (DEBUGLEVEL>2) timer_printf(&ti, "splitting mod %Ps", pr);
    1862         324 :   T.pr = pr;
    1863         324 :   T.L  = &L;
    1864         324 :   T.polbase = polbase;
    1865         324 :   T.pol   = pol;
    1866         324 :   T.nf    = nf;
    1867         324 :   res = nf_combine_factors(&T, polred, dpol-1);
    1868         324 :   if (DEBUGLEVEL>2)
    1869           0 :     err_printf("Total Time: %ld\n===========\n", timer_delay(&ti_tot));
    1870         324 :   return res;
    1871             : }
    1872             : 
    1873             : /* assume pol monic in nf.zk[X] */
    1874             : GEN
    1875          72 : nfroots_if_split(GEN *pnf, GEN pol)
    1876             : {
    1877          72 :   GEN T = get_nfpol(*pnf,pnf), den = fix_nf(pnf, &T, &pol);
    1878          72 :   pari_sp av = avma;
    1879          72 :   GEN z = nfsqff(*pnf, pol, ROOTS_SPLIT, den);
    1880          72 :   if (lg(z) == 1) { avma = av; return NULL; }
    1881          36 :   return gerepilecopy(av, z);
    1882             : }
    1883             : 
    1884             : /*******************************************************************/
    1885             : /*                                                                 */
    1886             : /*              Roots of unity in a number field                   */
    1887             : /*     (alternative to nfrootsof1 using factorization in K[X])     */
    1888             : /*                                                                 */
    1889             : /*******************************************************************/
    1890             : /* Code adapted from nffactor. Structure of the algorithm; only step 1 is
    1891             :  * specific to roots of unity.
    1892             :  *
    1893             :  * [Step 1]: guess roots via ramification. If trivial output this.
    1894             :  * [Step 2]: select prime [p] unramified and ideal [pr] above
    1895             :  * [Step 3]: evaluate the maximal exponent [k] such that the fondamental domain
    1896             :  *           of a LLL-reduction of [prk] = pr^k contains a ball of radius larger
    1897             :  *           than the norm of any root of unity.
    1898             :  * [Step 3]: select a heuristic exponent,
    1899             :  *           LLL reduce prk=pr^k and verify the exponent is sufficient,
    1900             :  *           otherwise try a larger one.
    1901             :  * [Step 4]: factor the cyclotomic polynomial mod [pr],
    1902             :  *           Hensel lift to pr^k and find the representative in the ball
    1903             :  *           If there is it is a primitive root */
    1904             : 
    1905             : typedef struct {
    1906             :   GEN q;
    1907             :   GEN modpr;
    1908             :   GEN pr;
    1909             :   nflift_t *L;
    1910             : } prklift_t;
    1911             : 
    1912             : /* Choose prime ideal unramified with "large" inertia degree */
    1913             : static void
    1914         354 : nf_pick_prime_for_units(GEN nf, prklift_t *P)
    1915             : {
    1916         354 :   GEN nfpol = nf_get_pol(nf), bad = mulii(nf_get_disc(nf), nf_get_index(nf));
    1917         354 :   GEN aT, amodpr, apr, ap = NULL, r = NULL;
    1918         354 :   long nfdeg = degpol(nfpol), maxf = get_maxf(nfdeg);
    1919             :   ulong pp;
    1920             :   forprime_t S;
    1921             : 
    1922         354 :   (void)u_forprime_init(&S, 2, ULONG_MAX);
    1923         354 :   while ( (pp = u_forprime_next(&S)) )
    1924             :   {
    1925        3786 :     if (! umodiu(bad,pp)) continue;
    1926        2862 :     r = get_good_factor(nfpol, pp, maxf);
    1927        2862 :     if (r) break;
    1928             :   }
    1929         354 :   if (!r) pari_err_OVERFLOW("nf_pick_prime [ran out of primes]");
    1930         354 :   ap = utoipos(pp);
    1931         354 :   apr = idealprimedec_kummer(nf, Flx_to_ZX(r), 1, ap);
    1932         354 :   amodpr = zk_to_Fq_init(nf,&apr,&aT,&ap);
    1933         354 :   P->pr = apr;
    1934         354 :   P->q = pr_norm(apr);
    1935         354 :   P->modpr = amodpr;
    1936         354 :   P->L->p = ap;
    1937         354 :   P->L->Tp = aT;
    1938         354 :   P->L->tozk = nf_get_invzk(nf);
    1939         354 :   P->L->topow = nf_get_zkprimpart(nf);
    1940         354 :   P->L->topowden = nf_get_zkden(nf);
    1941         354 : }
    1942             : 
    1943             : /* *Heuristic* exponent k such that the fundamental domain of pr^k
    1944             :  * should contain the ball of radius C */
    1945             : static double
    1946         354 : mybestlift_bound(GEN C)
    1947             : {
    1948         354 :   C = gtofp(C,DEFAULTPREC);
    1949             : #if 0 /* d = nf degree, Npr = Norm(pr) */
    1950             :   const double alpha = 0.99; /* LLL parameter */
    1951             :   const double y = 1 / (alpha - 0.25); /* = 2 if alpha = 3/4 */
    1952             :   double t;
    1953             :   t = rtodbl(mplog(gmul2n(divru(C,d), 4))) * 0.5 + (d-1) * log(1.5 * sqrt(y));
    1954             :   return ceil((t * d) / log(gtodouble(Npr))); /* proved upper bound */
    1955             : #endif
    1956         354 :   return ceil(log(gtodouble(C)) / 0.2) + 3;
    1957             : }
    1958             : 
    1959             : /* simplified nf_DDF_roots: polcyclo(n) monic in ZX either splits or has no
    1960             :  * root in nf.
    1961             :  * Return a root or NULL (no root) */
    1962             : static GEN
    1963         366 : nfcyclo_root(long n, GEN nfpol, nflift_t *L)
    1964             : {
    1965         366 :   GEN q, r, Cltx_r, pol = polcyclo(n,0), gn = utoipos(n);
    1966             :   div_data D;
    1967             : 
    1968         366 :   init_div_data(&D, pol, L);
    1969         366 :   (void)Fq_sqrtn(gen_1, gn, L->Tp, L->p, &r);
    1970             :   /* r primitive n-th root of 1 in Fq */
    1971         366 :   r = Zq_sqrtnlift(gen_1, gn, r, L->Tpk, L->p, L->k);
    1972             :   /* lt*dn*topowden * r = Clt * r */
    1973         366 :   r = nf_bestlift_to_pol(r, NULL, L);
    1974         366 :   Cltx_r = deg1pol_shallow(D.Clt? D.Clt: gen_1, gneg(r), varn(pol));
    1975             :   /* check P(r) == 0 */
    1976         366 :   q = RgXQX_divrem(D.C2ltpol, Cltx_r, nfpol, ONLY_DIVIDES); /* integral */
    1977         366 :   if (!q) return NULL;
    1978         342 :   if (D.Clt) r = gdiv(r, D.Clt);
    1979         342 :   return r;
    1980             : }
    1981             : 
    1982             : /* Guesses the number of roots of unity in number field [nf].
    1983             :  * Computes gcd of N(P)-1 for some primes. The value returned is a proven
    1984             :  * multiple of the correct value. */
    1985             : static long
    1986         822 : guess_roots(GEN nf)
    1987             : {
    1988         822 :   long c = 0, nfdegree = nf_get_degree(nf), B = nfdegree + 20, l;
    1989         822 :   ulong p = 2;
    1990         822 :   GEN T = nf_get_pol(nf), D = nf_get_disc(nf), index = nf_get_index(nf);
    1991         822 :   GEN nbroots = NULL;
    1992             :   forprime_t S;
    1993             :   pari_sp av;
    1994             : 
    1995         822 :   (void)u_forprime_init(&S, 3, ULONG_MAX);
    1996         822 :   av = avma;
    1997             :   /* result must be stationary (counter c) for at least B loops */
    1998       28320 :   for (l=1; (p = u_forprime_next(&S)); l++)
    1999             :   {
    2000             :     GEN old, F, pf_1, Tp;
    2001       28320 :     long i, nb, gcdf = 0;
    2002             : 
    2003       28320 :     if (!umodiu(D,p) || !umodiu(index,p)) continue;
    2004       26850 :     Tp = ZX_to_Flx(T,p); /* squarefree */
    2005       26850 :     F = Flx_nbfact_by_degree(Tp, &nb, p);
    2006             :     /* the gcd of the p^f - 1 is p^(gcd of the f's) - 1 */
    2007      233688 :     for (i = 1; i <= nfdegree; i++)
    2008      214020 :       if (F[i]) {
    2009       27174 :         gcdf = gcdf? cgcd(gcdf, i): i;
    2010       27174 :         if (gcdf == 1) break;
    2011             :       }
    2012       26850 :     pf_1 = subiu(powuu(p, gcdf), 1);
    2013       26850 :     old = nbroots;
    2014       26850 :     nbroots = nbroots? gcdii(pf_1, nbroots): pf_1;
    2015       26850 :     if (DEBUGLEVEL>5)
    2016           0 :       err_printf("p=%lu; gcf(f(P/p))=%ld; nbroots | %Ps",p, gcdf, nbroots);
    2017             :     /* if same result go on else reset the stop counter [c] */
    2018       26850 :     if (old && equalii(nbroots,old))
    2019       24996 :     { if (!is_bigint(nbroots) && ++c > B) break; }
    2020             :     else
    2021        1854 :       c = 0;
    2022             :   }
    2023         822 :   if (!nbroots) pari_err_OVERFLOW("guess_roots [ran out of primes]");
    2024         822 :   if (DEBUGLEVEL>5) err_printf("%ld loops\n",l);
    2025         822 :   avma = av; return itos(nbroots);
    2026             : }
    2027             : 
    2028             : /* T(x) an irreducible ZX. Is it of the form Phi_n(c \pm x) ?
    2029             :  * Return NULL if not, and a root of 1 of maximal order in Z[x]/(T) otherwise
    2030             :  *
    2031             :  * N.B. Set n_squarefree = 1 if n is squarefree, and 0 otherwise.
    2032             :  * This last parameter is inconvenient, but it allows a cheap
    2033             :  * stringent test. (n guessed from guess_roots())*/
    2034             : static GEN
    2035         270 : ZXirred_is_cyclo_translate(GEN T, long n_squarefree)
    2036             : {
    2037         270 :   long r, m, d = degpol(T);
    2038         270 :   GEN T1, c = divis_rem(gel(T, d+1), d, &r); /* d-1 th coeff divisible by d ? */
    2039             :   /* The trace coefficient of polcyclo(n) is \pm1 if n is square free, and 0
    2040             :    * otherwise. */
    2041         270 :   if (!n_squarefree)
    2042         126 :   { if (r) return NULL; }
    2043             :   else
    2044             :   {
    2045         144 :     if (r < -1)
    2046             :     {
    2047           0 :       r += d;
    2048           0 :       c = subiu(c, 1);
    2049             :     }
    2050         144 :     else if (r == d-1)
    2051             :     {
    2052          30 :       r = -1;
    2053          30 :       c = addiu(c, 1);
    2054             :     }
    2055         144 :     if (r != 1 && r != -1) return NULL;
    2056             :   }
    2057         228 :   if (signe(c)) /* presumably Phi_guess(c \pm x) */
    2058          30 :     T = RgX_translate(T, negi(c));
    2059         228 :   if (!n_squarefree) T = RgX_deflate_max(T, &m);
    2060             :   /* presumably Phi_core(guess)(\pm x), cyclotomic iff original T was */
    2061         228 :   T1 = ZX_graeffe(T);
    2062         228 :   if (ZX_equal(T, T1)) /* T = Phi_n, n odd */
    2063          30 :     return deg1pol_shallow(gen_m1, negi(c), varn(T));
    2064         198 :   else if (ZX_equal(T1, ZX_unscale(T, gen_m1))) /* T = Phi_{2n}, nodd */
    2065         180 :     return deg1pol_shallow(gen_1, c, varn(T));
    2066          18 :   return NULL;
    2067             : }
    2068             : 
    2069             : static GEN
    2070        1326 : trivroots(void) { return mkvec2(gen_2, gen_m1); }
    2071             : /* Number of roots of unity in number field [nf]. */
    2072             : GEN
    2073        1890 : rootsof1(GEN nf)
    2074             : {
    2075             :   prklift_t P;
    2076             :   nflift_t L;
    2077             :   GEN fa, LP, LE, C0, z, prim_root, disc, nfpol;
    2078             :   pari_timer ti;
    2079             :   long i, l, nbguessed, nbroots, nfdegree;
    2080             :   pari_sp av;
    2081             : 
    2082        1890 :   nf = checknf(nf);
    2083        1890 :   if (nf_get_r1(nf)) return trivroots();
    2084             : 
    2085             :   /* Step 1 : guess number of roots and discard trivial case 2 */
    2086         822 :   if (DEBUGLEVEL>2) timer_start(&ti);
    2087         822 :   nbguessed = guess_roots(nf);
    2088         822 :   if (DEBUGLEVEL>2)
    2089           0 :     timer_printf(&ti, "guessing roots of 1 [guess = %ld]", nbguessed);
    2090         822 :   if (nbguessed == 2) return trivroots();
    2091             : 
    2092         564 :   nfdegree = nf_get_degree(nf);
    2093         564 :   fa = factoru(nbguessed);
    2094         564 :   LP = gel(fa,1); l = lg(LP);
    2095         564 :   LE = gel(fa,2);
    2096         564 :   disc = nf_get_disc(nf);
    2097        1542 :   for (i = 1; i < l; i++)
    2098             :   {
    2099         978 :     long p = LP[i];
    2100             :     /* Degree and ramification test: find largest k such that Q(zeta_{p^k})
    2101             :      * may be a subfield of K. Q(zeta_p^k) has degree (p-1)p^(k-1)
    2102             :      * and v_p(discriminant) = ((p-1)k-1)p^(k-1); so we must have
    2103             :      * v_p(disc_K) >= ((p-1)k-1) * n / (p-1) = kn - q, where q = n/(p-1) */
    2104         978 :     if (p == 2)
    2105             :     { /* the test simplifies a little in that case */
    2106             :       long v, vnf, k;
    2107         564 :       if (LE[i] == 1) continue;
    2108         198 :       vnf = vals(nfdegree);
    2109         198 :       v = vali(disc);
    2110         222 :       for (k = minss(LE[i], vnf+1); k >= 1; k--)
    2111         222 :         if (v >= nfdegree*(k-1)) { nbguessed >>= LE[i]-k; LE[i] = k; break; }
    2112             :       /* N.B the test above always works for k = 1: LE[i] >= 1 */
    2113             :     }
    2114             :     else
    2115             :     {
    2116             :       long v, vnf, k;
    2117         414 :       ulong r, q = udivuu_rem(nfdegree, p-1, &r);
    2118         414 :       if (r) { nbguessed /= upowuu(p, LE[i]); LE[i] = 0; continue; }
    2119             :       /* q = n/(p-1) */
    2120         414 :       vnf = u_lval(q, p);
    2121         414 :       v = Z_lval(disc, p);
    2122         414 :       for (k = minss(LE[i], vnf+1); k >= 0; k--)
    2123         414 :         if (v >= nfdegree*k-(long)q)
    2124         414 :         { nbguessed /= upowuu(p, LE[i]-k); LE[i] = k; break; }
    2125             :       /* N.B the test above always works for k = 0: LE[i] >= 0 */
    2126             :     }
    2127             :   }
    2128         564 :   if (DEBUGLEVEL>2)
    2129           0 :     timer_printf(&ti, "after ramification conditions [guess = %ld]", nbguessed);
    2130         564 :   if (nbguessed == 2) return trivroots();
    2131         564 :   av = avma;
    2132             : 
    2133             :   /* Step 1.5 : test if nf.pol == subst(polcyclo(nbguessed), x, \pm x+c) */
    2134         564 :   if (eulerphiu_fact(fa) == (ulong)nfdegree)
    2135             :   {
    2136         270 :     GEN elt = ZXirred_is_cyclo_translate(nf_get_pol(nf),
    2137             :                                          uissquarefree_fact(fa));
    2138         270 :     if (elt)
    2139             :     {
    2140         210 :       if (DEBUGLEVEL>2)
    2141           0 :         timer_printf(&ti, "checking for cyclotomic polynomial [yes]");
    2142         210 :       return gerepilecopy(av, mkvec2(utoipos(nbguessed), elt));
    2143             :     }
    2144          60 :     avma = av;
    2145             :   }
    2146         354 :   if (DEBUGLEVEL>2)
    2147           0 :     timer_printf(&ti, "checking for cyclotomic polynomial [no]");
    2148             : 
    2149             :   /* Step 2 : choose a prime ideal for local lifting */
    2150         354 :   P.L = &L; nf_pick_prime_for_units(nf, &P);
    2151         354 :   if (DEBUGLEVEL>2)
    2152           0 :     timer_printf(&ti, "choosing prime %Ps, degree %ld",
    2153           0 :              P.L->p, P.L->Tp? degpol(P.L->Tp): 1);
    2154             : 
    2155             :   /* Step 3 : compute a reduced pr^k allowing lifting of local solutions */
    2156             :   /* evaluate maximum L2 norm of a root of unity in nf */
    2157         354 :   C0 = gmulsg(nfdegree, L2_bound(nf, gen_1));
    2158             :   /* lift and reduce pr^k */
    2159         354 :   if (DEBUGLEVEL>2) err_printf("Lift pr^k; GSmin wanted: %Ps\n",C0);
    2160         354 :   bestlift_init((long)mybestlift_bound(C0), nf, P.pr, C0, P.L);
    2161         354 :   P.L->dn = NULL;
    2162         354 :   if (DEBUGLEVEL>2) timer_start(&ti);
    2163             : 
    2164             :   /* Step 4 : actual computation of roots */
    2165         354 :   nbroots = 2; prim_root = gen_m1;
    2166         354 :   nfpol = nf_get_pol(nf);
    2167        1002 :   for (i = 1; i < l; i++)
    2168             :   { /* for all prime power factors of nbguessed, find a p^k-th root of unity */
    2169         648 :     long k, p = LP[i];
    2170         954 :     for (k = minss(LE[i], Z_lval(subiu(P.q,1UL),p)); k > 0; k--)
    2171             :     { /* find p^k-th roots */
    2172         648 :       pari_sp av = avma;
    2173         648 :       long pk = upowuu(p,k);
    2174         648 :       if (pk==2) continue; /* no need to test second roots ! */
    2175         366 :       z = nfcyclo_root(pk, nfpol, P.L);
    2176         366 :       if (DEBUGLEVEL>2) timer_printf(&ti, "for factoring Phi_%ld^%ld", p,k);
    2177         366 :       if (z) {
    2178         342 :         if (DEBUGLEVEL>2) err_printf("  %ld-th root of unity found.\n", pk);
    2179         342 :         if (p==2) { nbroots = pk; prim_root = z; }
    2180         276 :         else     { nbroots *= pk; prim_root = nfmul(nf, prim_root,z); }
    2181         342 :         break;
    2182             :       }
    2183          24 :       avma = av;
    2184          24 :       if (DEBUGLEVEL) pari_warn(warner,"rootsof1: wrong guess");
    2185             :     }
    2186             :   }
    2187         354 :   return gerepilecopy(av, mkvec2(utoi(nbroots), prim_root));
    2188             : }
    2189             : 
    2190             : static long
    2191           0 : zk_equal1(GEN y)
    2192             : {
    2193           0 :   if (typ(y) == t_INT) return equali1(y);
    2194           0 :   return equali1(gel(y,1)) && ZV_isscalar(y);
    2195             : }
    2196             : /* x^w = 1 */
    2197             : static GEN
    2198           0 : is_primitive_root(GEN nf, GEN fa, GEN x, long w)
    2199             : {
    2200           0 :   GEN P = gel(fa,1);
    2201           0 :   long i, l = lg(P);
    2202             : 
    2203           0 :   for (i = 1; i < l; i++)
    2204             :   {
    2205           0 :     long p = itos(gel(P,i));
    2206           0 :     GEN y = nfpow_u(nf,x, w/p);
    2207           0 :     if (zk_equal1(y) > 0) /* y = 1 */
    2208             :     {
    2209           0 :       if (p != 2 || !equali1(gcoeff(fa,i,2))) return NULL;
    2210           0 :       x = gneg_i(x);
    2211             :     }
    2212             :   }
    2213           0 :   return x;
    2214             : }
    2215             : GEN
    2216           0 : rootsof1_kannan(GEN nf)
    2217             : {
    2218           0 :   pari_sp av = avma;
    2219             :   long N, k, i, ws, prec;
    2220             :   GEN z, y, d, list, w;
    2221             : 
    2222           0 :   nf = checknf(nf);
    2223           0 :   if ( nf_get_r1(nf) ) return trivroots();
    2224             : 
    2225           0 :   N = nf_get_degree(nf); prec = nf_get_prec(nf);
    2226             :   for (;;)
    2227             :   {
    2228           0 :     GEN R = R_from_QR(nf_get_G(nf), prec);
    2229           0 :     if (R)
    2230             :     {
    2231           0 :       y = fincke_pohst(mkvec(R), utoipos(N), N * N, 0, NULL);
    2232           0 :       if (y) break;
    2233             :     }
    2234           0 :     prec = precdbl(prec);
    2235           0 :     if (DEBUGLEVEL) pari_warn(warnprec,"rootsof1",prec);
    2236           0 :     nf = nfnewprec_shallow(nf,prec);
    2237           0 :   }
    2238           0 :   if (itos(ground(gel(y,2))) != N) pari_err_BUG("rootsof1 (bug1)");
    2239           0 :   w = gel(y,1); ws = itos(w);
    2240           0 :   if (ws == 2) { avma = av; return trivroots(); }
    2241             : 
    2242           0 :   d = Z_factor(w); list = gel(y,3); k = lg(list);
    2243           0 :   for (i=1; i<k; i++)
    2244             :   {
    2245           0 :     z = is_primitive_root(nf, d, gel(list,i), ws);
    2246           0 :     if (z) return gerepilecopy(av, mkvec2(w, z));
    2247             :   }
    2248           0 :   pari_err_BUG("rootsof1");
    2249             :   return NULL; /* LCOV_EXCL_LINE */
    2250             : }

Generated by: LCOV version 1.11