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 - kummer.c (source / functions) Hit Total Coverage
Test: PARI/GP v2.10.0 lcov report (development 21059-cbe0d6a) Lines: 719 836 86.0 %
Date: 2017-09-22 06:24:58 Functions: 51 56 91.1 %
Legend: Lines: hit not hit

          Line data    Source code
       1             : /* Copyright (C) 2000  The PARI group.
       2             : 
       3             : This file is part of the PARI/GP package.
       4             : 
       5             : PARI/GP is free software; you can redistribute it and/or modify it under the
       6             : terms of the GNU General Public License as published by the Free Software
       7             : Foundation. It is distributed in the hope that it will be useful, but WITHOUT
       8             : ANY WARRANTY WHATSOEVER.
       9             : 
      10             : Check the License for details. You should have received a copy of it, along
      11             : with the package; see the file 'COPYING'. If not, write to the Free Software
      12             : Foundation, Inc., 51 Franklin Street, Fifth Floor, Boston, MA 02110-1301 USA. */
      13             : 
      14             : /*******************************************************************/
      15             : /*                                                                 */
      16             : /*                      KUMMER EXTENSIONS                          */
      17             : /*                                                                 */
      18             : /*******************************************************************/
      19             : #include "pari.h"
      20             : #include "paripriv.h"
      21             : 
      22             : typedef struct {
      23             :   GEN R; /* nf.pol */
      24             :   GEN x; /* tau ( Mod(x, R) ) */
      25             :   GEN zk;/* action of tau on nf.zk (as t_MAT) */
      26             : } tau_s;
      27             : 
      28             : typedef struct {
      29             :   GEN polnf, invexpoteta1, powg;
      30             :   tau_s *tau;
      31             :   long m;
      32             : } toK_s;
      33             : 
      34             : typedef struct {
      35             :   GEN R; /* ZX, compositum(P,Q) */
      36             :   GEN p; /* QX, Mod(p,R) root of P */
      37             :   GEN q; /* QX, Mod(q,R) root of Q */
      38             :   long k; /* Q[X]/R generated by q + k p */
      39             :   GEN rev;
      40             : } compo_s;
      41             : 
      42             : static long
      43        1085 : prank(GEN cyc, long ell)
      44             : {
      45             :   long i;
      46        3101 :   for (i=1; i<lg(cyc); i++)
      47        2338 :     if (smodis(gel(cyc,i),ell)) break;
      48        1085 :   return i-1;
      49             : }
      50             : 
      51             : /* increment y, which runs through [0,d-1]^(k-1). Return 0 when done. */
      52             : static int
      53         140 : increment(GEN y, long k, long d)
      54             : {
      55         140 :   long i = k, j;
      56             :   do
      57             :   {
      58         161 :     if (--i == 0) return 0;
      59         133 :     y[i]++;
      60         133 :   } while (y[i] >= d);
      61         112 :   for (j = i+1; j < k; j++) y[j] = 0;
      62         112 :   return 1;
      63             : }
      64             : 
      65             : static int
      66         539 : ok_congruence(GEN X, ulong ell, long lW, GEN vecMsup)
      67             : {
      68             :   long i, l;
      69         539 :   if (zv_equal0(X)) return 0;
      70         539 :   l = lg(X);
      71         924 :   for (i=lW; i<l; i++)
      72         413 :     if (X[i] == 0) return 0;
      73         511 :   l = lg(vecMsup);
      74         798 :   for (i=1; i<l; i++)
      75         287 :     if (zv_equal0(Flm_Flc_mul(gel(vecMsup,i),X, ell))) return 0;
      76         511 :   return 1;
      77             : }
      78             : 
      79             : static int
      80         259 : ok_sign(GEN X, GEN msign, GEN arch)
      81             : {
      82         259 :   return zv_equal(Flm_Flc_mul(msign, X, 2), arch);
      83             : }
      84             : 
      85             : /* REDUCTION MOD ell-TH POWERS */
      86             : 
      87             : #if 0
      88             : static GEN
      89             : logarch2arch(GEN x, long r1, long prec)
      90             : {
      91             :   long i, lx;
      92             :   GEN y = cgetg_copy(x, &lx);
      93             :   if (typ(x) == t_MAT)
      94             :   {
      95             :     for (i=1; i<lx; i++) gel(y,i) = logarch2arch(gel(x,i), r1, prec);
      96             :   }
      97             :   else
      98             :   {
      99             :     for (i=1; i<=r1;i++) gel(y,i) = gexp(gel(x,i),prec);
     100             :     for (   ; i<lx; i++) gel(y,i) = gexp(gmul2n(gel(x,i),-1),prec);
     101             :   }
     102             :   return y;
     103             : }
     104             : #endif
     105             : 
     106             : /* make be integral by multiplying by t in (Q^*)^ell */
     107             : static GEN
     108         434 : reduce_mod_Qell(GEN bnfz, GEN be, GEN gell)
     109             : {
     110             :   GEN c;
     111         434 :   be = nf_to_scalar_or_basis(bnfz, be);
     112         434 :   be = Q_primitive_part(be, &c);
     113         434 :   if (c)
     114             :   {
     115         238 :     GEN d, fa = factor(c);
     116         238 :     gel(fa,2) = FpC_red(gel(fa,2), gell);
     117         238 :     d = factorback(fa);
     118         238 :     be = typ(be) == t_INT? mulii(be,d): ZC_Z_mul(be, d);
     119             :   }
     120         434 :   return be;
     121             : }
     122             : 
     123             : /* return q, q^n r = x, v_pr(r) < n for all pr. Insist q is a genuine n-th
     124             :  * root (i.e r = 1) if strict != 0. */
     125             : static GEN
     126        1260 : idealsqrtn(GEN nf, GEN x, GEN gn, int strict)
     127             : {
     128        1260 :   long i, l, n = itos(gn);
     129             :   GEN fa, q, Ex, Pr;
     130             : 
     131        1260 :   fa = idealfactor(nf, x);
     132        1260 :   Pr = gel(fa,1); l = lg(Pr);
     133        1260 :   Ex = gel(fa,2); q = NULL;
     134        3395 :   for (i=1; i<l; i++)
     135             :   {
     136        2135 :     long ex = itos(gel(Ex,i));
     137        2135 :     GEN e = stoi(ex / n);
     138        2135 :     if (strict && ex % n) pari_err_SQRTN("idealsqrtn", fa);
     139        2135 :     if (q) q = idealmulpowprime(nf, q, gel(Pr,i), e);
     140         469 :     else   q = idealpow(nf, gel(Pr,i), e);
     141             :   }
     142        1260 :   return q? q: gen_1;
     143             : }
     144             : 
     145             : static GEN
     146         434 : reducebeta(GEN bnfz, GEN b, GEN ell)
     147             : {
     148         434 :   long prec = nf_get_prec(bnfz);
     149         434 :   GEN y, elllogfu, nf = bnf_get_nf(bnfz), fu = bnf_get_fu_nocheck(bnfz);
     150             : 
     151         434 :   if (DEBUGLEVEL>1) err_printf("reducing beta = %Ps\n",b);
     152         434 :   b = reduce_mod_Qell(nf, b, ell);
     153             :   /* reduce l-th root */
     154         434 :   y = idealsqrtn(nf, b, ell, 0); /* (b) = y^ell I, I integral */
     155         434 :   if (typ(y) == t_MAT && !is_pm1(gcoeff(y,1,1)))
     156             :   {
     157         168 :     GEN T = idealred(nf, mkvec2(y, gen_1)), t = gel(T,2);
     158             :     /* (t)*T[1] = y, T[1] integral and small */
     159         168 :     if (gcmp(idealnorm(nf,t), gen_1) > 0)
     160         154 :       b = nfmul(nf, b, nfpow(nf, t, negi(ell)));
     161             :   }
     162         434 :   if (DEBUGLEVEL>1) err_printf("beta reduced via ell-th root = %Ps\n",b);
     163             :   /* log. embeddings of fu^ell */
     164         434 :   elllogfu = RgM_Rg_mul(real_i(bnf_get_logfu(bnfz)), ell);
     165             :   for (;;)
     166             :   {
     167         455 :     GEN emb, z = get_arch_real(nf, b, &emb, prec);
     168         455 :     if (z)
     169             :     {
     170         434 :       GEN ex = RgM_Babai(elllogfu, z);
     171         434 :       if (ex)
     172             :       {
     173         434 :         b = nfdiv(nf, b, nffactorback(nf, fu, RgC_Rg_mul(ex,ell)));
     174         434 :         break;
     175             :       }
     176             :     }
     177          21 :     prec = precdbl(prec);
     178          21 :     if (DEBUGLEVEL) pari_warn(warnprec,"reducebeta",prec);
     179          21 :     nf = nfnewprec_shallow(nf,prec);
     180          21 :   }
     181         434 :   if (DEBUGLEVEL>1) err_printf("beta LLL-reduced mod U^l = %Ps\n",b);
     182         434 :   return b;
     183             : }
     184             : 
     185             : /* FIXME: remove */
     186             : static GEN
     187         455 : tauofalg(GEN x, tau_s *tau) {
     188         455 :   long tx = typ(x);
     189         455 :   if (tx == t_POLMOD) { x = gel(x,2); tx = typ(x); }
     190         455 :   if (tx == t_POL) x = RgX_RgXQ_eval(x, tau->x, tau->R);
     191         455 :   return mkpolmod(x, tau->R);
     192             : }
     193             : 
     194             : /* compute Gal(K(\zeta_l)/K) */
     195             : static tau_s *
     196         231 : get_tau(tau_s *tau, GEN nf, compo_s *C, long g)
     197             : {
     198             :   GEN U;
     199             : 
     200             :   /* compute action of tau: q^g + kp */
     201         231 :   U = RgX_add(RgXQ_powu(C->q, g, C->R), RgX_muls(C->p, C->k));
     202         231 :   U = RgX_RgXQ_eval(C->rev, U, C->R);
     203             : 
     204         231 :   tau->x  = U;
     205         231 :   tau->R  = C->R;
     206         231 :   tau->zk = nfgaloismatrix(nf, U);
     207         231 :   return tau;
     208             : }
     209             : 
     210             : static GEN tauoffamat(GEN x, tau_s *tau);
     211             : 
     212             : static GEN
     213        9940 : tauofelt(GEN x, tau_s *tau)
     214             : {
     215        9940 :   switch(typ(x))
     216             :   {
     217        8323 :     case t_COL: return RgM_RgC_mul(tau->zk, x);
     218        1162 :     case t_MAT: return tauoffamat(x, tau);
     219         455 :     default: return tauofalg(x, tau);
     220             :   }
     221             : }
     222             : static GEN
     223        1330 : tauofvec(GEN x, tau_s *tau)
     224             : {
     225             :   long i, l;
     226        1330 :   GEN y = cgetg_copy(x, &l);
     227        1330 :   for (i=1; i<l; i++) gel(y,i) = tauofelt(gel(x,i), tau);
     228        1330 :   return y;
     229             : }
     230             : /* [x, tau(x), ..., tau^(m-1)(x)] */
     231             : static GEN
     232         630 : powtau(GEN x, long m, tau_s *tau)
     233             : {
     234         630 :   GEN y = cgetg(m+1, t_VEC);
     235             :   long i;
     236         630 :   gel(y,1) = x;
     237         630 :   for (i=2; i<=m; i++) gel(y,i) = tauofelt(gel(y,i-1), tau);
     238         630 :   return y;
     239             : }
     240             : /* x^lambda */
     241             : static GEN
     242         539 : lambdaofelt(GEN x, toK_s *T)
     243             : {
     244         539 :   tau_s *tau = T->tau;
     245         539 :   long i, m = T->m;
     246         539 :   GEN y = cgetg(1, t_MAT), powg = T->powg; /* powg[i] = g^i */
     247        1372 :   for (i=1; i<m; i++)
     248             :   {
     249         833 :     y = famat_mulpow_shallow(y, x, gel(powg,m-i));
     250         833 :     x = tauofelt(x, tau);
     251             :   }
     252         539 :   return famat_mul_shallow(y, x);
     253             : }
     254             : static GEN
     255         448 : lambdaofvec(GEN x, toK_s *T)
     256             : {
     257             :   long i, l;
     258         448 :   GEN y = cgetg_copy(x, &l);
     259         448 :   for (i=1; i<l; i++) gel(y,i) = lambdaofelt(gel(x,i), T);
     260         448 :   return y;
     261             : }
     262             : 
     263             : static GEN
     264        1162 : tauoffamat(GEN x, tau_s *tau)
     265             : {
     266        1162 :   return mkmat2(tauofvec(gel(x,1), tau), gel(x,2));
     267             : }
     268             : 
     269             : static GEN
     270         140 : tauofideal(GEN id, tau_s *tau)
     271             : {
     272         140 :   return ZM_hnfmodid(RgM_mul(tau->zk, id), gcoeff(id, 1,1));
     273             : }
     274             : 
     275             : static int
     276         644 : isprimeidealconj(GEN P, GEN Q, tau_s *tau)
     277             : {
     278         644 :   GEN p = pr_get_p(P);
     279         644 :   GEN x = pr_get_gen(P);
     280         644 :   if (!equalii(p, pr_get_p(Q))
     281         490 :    || pr_get_e(P) != pr_get_e(Q)
     282         490 :    || pr_get_f(P) != pr_get_f(Q)) return 0;
     283         483 :   if (ZV_equal(x, pr_get_gen(Q))) return 1;
     284             :   for(;;)
     285             :   {
     286        1295 :     if (ZC_prdvd(x,Q)) return 1;
     287         980 :     x = FpC_red(tauofelt(x, tau), p);
     288         980 :     if (ZC_prdvd(x,P)) return 0;
     289         812 :   }
     290             : }
     291             : 
     292             : static int
     293        1092 : isconjinprimelist(GEN S, GEN pr, tau_s *tau)
     294             : {
     295             :   long i, l;
     296             : 
     297        1092 :   if (!tau) return 0;
     298         756 :   l = lg(S);
     299        1085 :   for (i=1; i<l; i++)
     300         644 :     if (isprimeidealconj(gel(S,i),pr,tau)) return 1;
     301         441 :   return 0;
     302             : }
     303             : 
     304             : /* assume x in basistoalg form */
     305             : static GEN
     306        1036 : downtoK(toK_s *T, GEN x)
     307             : {
     308        1036 :   long degKz = lg(T->invexpoteta1) - 1;
     309        1036 :   GEN y = gmul(T->invexpoteta1, Rg_to_RgC(lift_shallow(x), degKz));
     310        1036 :   return gmodulo(gtopolyrev(y,varn(T->polnf)), T->polnf);
     311             : }
     312             : 
     313             : static GEN
     314           0 : no_sol(long all, long i)
     315             : {
     316           0 :   if (!all) pari_err_BUG(stack_sprintf("kummer [bug%ld]", i));
     317           0 :   return cgetg(1,t_VEC);
     318             : }
     319             : 
     320             : static GEN
     321         420 : get_gell(GEN bnr, GEN subgp, long all)
     322             : {
     323             :   GEN gell;
     324         420 :   if (all && all != -1) return utoipos(labs(all));
     325         399 :   if (!subgp) return ZV_prod(bnr_get_cyc(bnr));
     326         399 :   gell = det(subgp);
     327         399 :   if (typ(gell) != t_INT) pari_err_TYPE("rnfkummer",gell);
     328         399 :   return gell;
     329             : }
     330             : 
     331             : typedef struct {
     332             :   GEN Sm, Sml1, Sml2, Sl, ESml2;
     333             : } primlist;
     334             : 
     335             : static int
     336         413 : build_list_Hecke(primlist *L, GEN nfz, GEN fa, GEN gothf, GEN gell, tau_s *tau)
     337             : {
     338             :   GEN listpr, listex, pr, factell;
     339         413 :   long vp, i, l, ell = itos(gell), degKz = nf_get_degree(nfz);
     340             : 
     341         413 :   if (!fa) fa = idealfactor(nfz, gothf);
     342         413 :   listpr = gel(fa,1);
     343         413 :   listex = gel(fa,2); l = lg(listpr);
     344         413 :   L->Sm  = vectrunc_init(l);
     345         413 :   L->Sml1= vectrunc_init(l);
     346         413 :   L->Sml2= vectrunc_init(l);
     347         413 :   L->Sl  = vectrunc_init(l+degKz);
     348         413 :   L->ESml2=vecsmalltrunc_init(l);
     349        1274 :   for (i=1; i<l; i++)
     350             :   {
     351         861 :     pr = gel(listpr,i);
     352         861 :     vp = itos(gel(listex,i));
     353         861 :     if (!equalii(pr_get_p(pr), gell))
     354             :     {
     355         574 :       if (vp != 1) return 1;
     356         574 :       if (!isconjinprimelist(L->Sm,pr,tau)) vectrunc_append(L->Sm,pr);
     357             :     }
     358             :     else
     359             :     {
     360         287 :       long e = pr_get_e(pr), vd = (vp-1)*(ell-1)-ell*e;
     361         287 :       if (vd > 0) return 4;
     362         287 :       if (vd==0)
     363             :       {
     364          63 :         if (!isconjinprimelist(L->Sml1,pr,tau)) vectrunc_append(L->Sml1, pr);
     365             :       }
     366             :       else
     367             :       {
     368         224 :         if (vp==1) return 2;
     369         224 :         if (!isconjinprimelist(L->Sml2,pr,tau))
     370             :         {
     371         224 :           vectrunc_append(L->Sml2, pr);
     372         224 :           vecsmalltrunc_append(L->ESml2, vp);
     373             :         }
     374             :       }
     375             :     }
     376             :   }
     377         413 :   factell = idealprimedec(nfz,gell); l = lg(factell);
     378         931 :   for (i=1; i<l; i++)
     379             :   {
     380         518 :     pr = gel(factell,i);
     381         518 :     if (!idealval(nfz,gothf,pr) && !isconjinprimelist(L->Sl,pr,tau))
     382         224 :       vectrunc_append(L->Sl, pr);
     383             :   }
     384         413 :   return 0; /* OK */
     385             : }
     386             : 
     387             : /* Return a Flm */
     388             : static GEN
     389         672 : logall(GEN nf, GEN vec, long lW, long mginv, long ell, GEN pr, long ex)
     390             : {
     391         672 :   GEN m, M, sprk = zlog_pr_init(nf, pr, ex);
     392         672 :   long ellrank, i, l = lg(vec);
     393             : 
     394         672 :   ellrank = prank(gel(sprk,1), ell);
     395         672 :   M = cgetg(l,t_MAT);
     396        2835 :   for (i=1; i<l; i++)
     397             :   {
     398        2163 :     m = zlog_pr(nf, gel(vec,i), sprk);
     399        2163 :     setlg(m, ellrank+1);
     400        2163 :     if (i < lW) m = gmulsg(mginv, m);
     401        2163 :     gel(M,i) = ZV_to_Flv(m, ell);
     402             :   }
     403         672 :   return M;
     404             : }
     405             : 
     406             : /* compute the u_j (see remark 5.2.15.) */
     407             : static GEN
     408         413 : get_u(GEN cyc, long rc, GEN gell)
     409             : {
     410         413 :   long i, l = lg(cyc);
     411         413 :   GEN u = cgetg(l,t_VEC);
     412         413 :   for (i=1; i<=rc; i++) gel(u,i) = gen_0;
     413         413 :   for (   ; i<  l; i++) gel(u,i) = Fp_inv(gel(cyc,i), gell);
     414         413 :   return u;
     415             : }
     416             : 
     417             : /* alg. 5.2.15. with remark */
     418             : static GEN
     419         469 : isprincipalell(GEN bnfz, GEN id, GEN cycgen, GEN u, GEN gell, long rc)
     420             : {
     421         469 :   long i, l = lg(cycgen);
     422         469 :   GEN v, b, db, y = bnfisprincipal0(bnfz, id, nf_FORCE|nf_GENMAT);
     423             : 
     424         469 :   v = FpC_red(gel(y,1), gell);
     425         469 :   b = gel(y,2);
     426         469 :   if (typ(b) == t_COL)
     427             :   {
     428         406 :     b = Q_remove_denom(gel(y,2), &db);
     429         406 :     if (db) b = famat_mulpow_shallow(b, db, gen_m1);
     430             :   }
     431         672 :   for (i=rc+1; i<l; i++)
     432             :   {
     433         203 :     GEN e = Fp_mul(gel(v,i), gel(u,i), gell);
     434         203 :     b = famat_mulpow_shallow(b, gel(cycgen,i), e);
     435             :   }
     436         469 :   setlg(v,rc+1); return mkvec2(v, b);
     437             : }
     438             : 
     439             : static GEN
     440        1330 : famat_factorback(GEN v, GEN e)
     441             : {
     442        1330 :   long i, l = lg(e);
     443        1330 :   GEN V = cgetg(1, t_MAT);
     444        1330 :   for (i=1; i<l; i++) V = famat_mulpow_shallow(V, gel(v,i), gel(e,i));
     445        1330 :   return V;
     446             : }
     447             : 
     448             : static GEN
     449         434 : compute_beta(GEN X, GEN vecWB, GEN ell, GEN bnfz)
     450             : {
     451             :   GEN BE, be;
     452         434 :   BE = famat_reduce(famat_factorback(vecWB, zv_to_ZV(X)));
     453         434 :   gel(BE,2) = centermod(gel(BE,2), ell);
     454         434 :   be = nffactorback(bnfz, BE, NULL);
     455         434 :   be = reducebeta(bnfz, be, ell);
     456         434 :   if (DEBUGLEVEL>1) err_printf("beta reduced = %Ps\n",be);
     457         434 :   return be;
     458             : }
     459             : 
     460             : static GEN
     461         413 : get_Selmer(GEN bnf, GEN cycgen, long rc)
     462             : {
     463         413 :   GEN U = bnf_build_units(bnf), tu = gel(U,1), fu = vecslice(U, 2, lg(U)-1);
     464         413 :   return shallowconcat(shallowconcat(fu,mkvec(tu)), vecslice(cycgen,1,rc));
     465             : }
     466             : 
     467             : GEN
     468       34433 : lift_if_rational(GEN x)
     469             : {
     470             :   long lx, i;
     471             :   GEN y;
     472             : 
     473       34433 :   switch(typ(x))
     474             :   {
     475        4648 :     default: break;
     476             : 
     477             :     case t_POLMOD:
     478       20335 :       y = gel(x,2);
     479       20335 :       if (typ(y) == t_POL)
     480             :       {
     481        6538 :         long d = degpol(y);
     482        6538 :         if (d > 0) return x;
     483        1197 :         return (d < 0)? gen_0: gel(y,2);
     484             :       }
     485       13797 :       return y;
     486             : 
     487        3948 :     case t_POL: lx = lg(x);
     488        3948 :       for (i=2; i<lx; i++) gel(x,i) = lift_if_rational(gel(x,i));
     489        3948 :       break;
     490        5502 :     case t_VEC: case t_COL: case t_MAT: lx = lg(x);
     491        5502 :       for (i=1; i<lx; i++) gel(x,i) = lift_if_rational(gel(x,i));
     492             :   }
     493       14098 :   return x;
     494             : }
     495             : 
     496             : /* A column vector representing a subgroup of prime index */
     497             : static GEN
     498           0 : grptocol(GEN H)
     499             : {
     500           0 :   long i, j, l = lg(H);
     501           0 :   GEN col = cgetg(l, t_VECSMALL);
     502           0 :   for (i = 1; i < l; i++)
     503             :   {
     504           0 :     ulong ell = itou( gcoeff(H,i,i) );
     505           0 :     if (ell == 1) col[i] = 0; else { col[i] = ell-1; break; }
     506             :   }
     507           0 :   for (j=i; ++j < l; ) col[j] = itou( gcoeff(H,i,j) );
     508           0 :   return col;
     509             : }
     510             : 
     511             : /* Reorganize kernel basis so that the tests of ok_congruence can be ok
     512             :  * for y[ncyc]=1 and y[1..ncyc]=1 */
     513             : static GEN
     514           0 : fix_kernel(GEN K, GEN M, GEN vecMsup, long lW, long ell)
     515             : {
     516           0 :   pari_sp av = avma;
     517           0 :   long i, j, idx, ffree, dK = lg(K)-1;
     518           0 :   GEN Ki, Kidx = cgetg(dK+1, t_VECSMALL);
     519             : 
     520             :   /* First step: Gauss elimination on vectors lW...lg(M) */
     521           0 :   for (idx = lg(K), i=lg(M); --i >= lW; )
     522             :   {
     523           0 :     for (j=dK; j > 0; j--) if (coeff(K, i, j)) break;
     524           0 :     if (!j)
     525             :     { /* Do our best to ensure that K[dK,i] != 0 */
     526           0 :       if (coeff(K, i, dK)) continue;
     527           0 :       for (j = idx; j < dK; j++)
     528           0 :         if (coeff(K, i, j) && coeff(K, Kidx[j], dK) != ell - 1)
     529           0 :           Flv_add_inplace(gel(K,dK), gel(K,j), ell);
     530             :     }
     531           0 :     if (j != --idx) swap(gel(K, j), gel(K, idx));
     532           0 :     Kidx[idx] = i;
     533           0 :     if (coeff(K,i,idx) != 1)
     534           0 :       Flv_Fl_div_inplace(gel(K,idx), coeff(K,i,idx), ell);
     535           0 :     Ki = gel(K,idx);
     536           0 :     if (coeff(K,i,dK) != 1)
     537             :     {
     538           0 :       ulong t = Fl_sub(coeff(K,i,dK), 1, ell);
     539           0 :       Flv_sub_inplace(gel(K,dK), Flv_Fl_mul(Ki, t, ell), ell);
     540             :     }
     541           0 :     for (j = dK; --j > 0; )
     542             :     {
     543           0 :       if (j == idx) continue;
     544           0 :       if (coeff(K,i,j))
     545           0 :         Flv_sub_inplace(gel(K,j), Flv_Fl_mul(Ki, coeff(K,i,j), ell), ell);
     546             :     }
     547             :   }
     548             :   /* ffree = first vector that is not "free" for the scalar products */
     549           0 :   ffree = idx;
     550             :   /* Second step: for each hyperplane equation in vecMsup, do the same
     551             :    * thing as before. */
     552           0 :   for (i=1; i < lg(vecMsup); i++)
     553             :   {
     554           0 :     GEN Msup = gel(vecMsup,i);
     555             :     ulong dotprod;
     556           0 :     if (lgcols(Msup) != 2) continue;
     557           0 :     Msup = zm_row(Msup, 1);
     558           0 :     for (j=ffree; --j > 0; )
     559             :     {
     560           0 :       dotprod = Flv_dotproduct(Msup, gel(K,j), ell);
     561           0 :       if (dotprod)
     562             :       {
     563           0 :         if (j != --ffree) swap(gel(K, j), gel(K, ffree));
     564           0 :         if (dotprod != 1) Flv_Fl_div_inplace(gel(K, ffree), dotprod, ell);
     565           0 :         break;
     566             :       }
     567             :     }
     568           0 :     if (!j)
     569             :     { /* Do our best to ensure that vecMsup.K[dK] != 0 */
     570           0 :       if (Flv_dotproduct(Msup, gel(K,dK), ell) == 0)
     571             :       {
     572           0 :         for (j = ffree-1; j <= dK; j++)
     573           0 :           if (Flv_dotproduct(Msup, gel(K,j), ell)
     574           0 :               && coeff(K,Kidx[j],dK) != ell-1)
     575           0 :             Flv_add_inplace(gel(K,dK), gel(K,j), ell);
     576             :       }
     577           0 :       continue;
     578             :     }
     579           0 :     Ki = gel(K,ffree);
     580           0 :     dotprod = Flv_dotproduct(Msup, gel(K,dK), ell);
     581           0 :     if (dotprod != 1)
     582             :     {
     583           0 :       ulong t = Fl_sub(dotprod,1,ell);
     584           0 :       Flv_sub_inplace(gel(K,dK), Flv_Fl_mul(Ki,t,ell), ell);
     585             :     }
     586           0 :     for (j = dK; --j > 0; )
     587             :     {
     588           0 :       if (j == ffree) continue;
     589           0 :       dotprod = Flv_dotproduct(Msup, gel(K,j), ell);
     590           0 :       if (dotprod) Flv_sub_inplace(gel(K,j), Flv_Fl_mul(Ki,dotprod,ell), ell);
     591             :     }
     592             :   }
     593           0 :   if (ell == 2)
     594             :   {
     595           0 :     for (i = ffree, j = ffree-1; i <= dK && j; i++, j--)
     596           0 :     { swap(gel(K,i), gel(K,j)); }
     597             :   }
     598             :   /* Try to ensure that y = vec_ei(n, i) gives a good candidate */
     599           0 :   for (i = 1; i < dK; i++) Flv_add_inplace(gel(K,i), gel(K,dK), ell);
     600           0 :   return gerepilecopy(av, K);
     601             : }
     602             : 
     603             : static GEN
     604           0 : Flm_init(long m, long n)
     605             : {
     606           0 :   GEN M = cgetg(n+1, t_MAT);
     607           0 :   long i; for (i = 1; i <= n; i++) gel(M,i) = cgetg(m+1, t_VECSMALL);
     608           0 :   return M;
     609             : }
     610             : static void
     611           0 : Flv_fill(GEN v, GEN y)
     612             : {
     613           0 :   long i, l = lg(y);
     614           0 :   for (i = 1; i < l; i++) v[i] = y[i];
     615           0 : }
     616             : 
     617             : static GEN
     618         616 : get_badbnf(GEN bnf)
     619             : {
     620             :   long i, l;
     621         616 :   GEN bad = gen_1, gen = bnf_get_gen(bnf);
     622         616 :   l = lg(gen);
     623        1064 :   for (i = 1; i < l; i++)
     624             :   {
     625         448 :     GEN g = gel(gen,i);
     626         448 :     bad = lcmii(bad, gcoeff(g,1,1));
     627             :   }
     628         616 :   return bad;
     629             : }
     630             : /* Let K base field, L/K described by bnr (conductor f) + H. Return a list of
     631             :  * primes coprime to f*ell of degree 1 in K whose images in Cl_f(K) generate H:
     632             :  * thus they all split in Lz/Kz; t in Kz is such that
     633             :  * t^(1/p) generates Lz => t is an ell-th power in k(pr) for all such primes.
     634             :  * Restrict to primes not dividing
     635             :  * - the index fz of the polynomial defining Kz, or
     636             :  * - the modulus, or
     637             :  * - ell, or
     638             :  * - a generator in bnf.gen or bnfz.gen */
     639             : static GEN
     640         399 : get_prlist(GEN bnr, GEN H, ulong ell, GEN bnfz)
     641             : {
     642         399 :   pari_sp av0 = avma;
     643             :   forprime_t T;
     644             :   ulong p;
     645             :   GEN L, nf, cyc, bad, cond, condZ, Hsofar;
     646         399 :   L = cgetg(1, t_VEC);
     647         399 :   cyc = bnr_get_cyc(bnr);
     648         399 :   nf = bnr_get_nf(bnr);
     649             : 
     650         399 :   cond = gel(bnr_get_mod(bnr), 1);
     651         399 :   condZ = gcoeff(cond,1,1);
     652         399 :   bad = get_badbnf(bnr_get_bnf(bnr));
     653         399 :   if (bnfz)
     654             :   {
     655         217 :     GEN badz = lcmii(get_badbnf(bnfz), nf_get_index(bnf_get_nf(bnfz)));
     656         217 :     bad = mulii(bad,badz);
     657             :   }
     658         399 :   bad = lcmii(muliu(condZ, ell), bad);
     659             :   /* restrict to primes not dividing bad */
     660             : 
     661         399 :   u_forprime_init(&T, 2, ULONG_MAX);
     662         399 :   Hsofar = cgetg(1, t_MAT);
     663        7266 :   while ((p = u_forprime_next(&T)))
     664             :   {
     665             :     GEN LP;
     666             :     long i, l;
     667        6867 :     if (p == ell || !umodiu(bad, p)) continue;
     668        5775 :     LP = idealprimedec_limit_f(nf, utoipos(p), 1);
     669        5775 :     l = lg(LP);
     670        9387 :     for (i = 1; i < l; i++)
     671             :     {
     672        4011 :       pari_sp av = avma;
     673        4011 :       GEN M, P = gel(LP,i), v = bnrisprincipal(bnr, P, 0);
     674        4011 :       if (!hnf_invimage(H, v)) { avma = av; continue; }
     675        1218 :       M = shallowconcat(Hsofar, v);
     676        1218 :       M = ZM_hnfmodid(M, cyc);
     677        1218 :       if (ZM_equal(M, Hsofar)) continue;
     678         833 :       L = shallowconcat(L, mkvec(P));
     679         833 :       Hsofar = M;
     680             :       /* the primes in L generate H */
     681         833 :       if (ZM_equal(M, H)) return gerepilecopy(av0, L);
     682             :     }
     683             :   }
     684           0 :   pari_err_BUG("rnfkummer [get_prlist]");
     685           0 :   return NULL;
     686             : }
     687             : /*Lprz list of prime ideals in Kz that must split completely in Lz/Kz, vecWA
     688             :  * generators for the S-units used to build the Kummer generators. Return
     689             :  * matsmall M such that \prod WA[j]^x[j] ell-th power mod pr[i] iff
     690             :  * \sum M[i,j] x[j] = 0 (mod ell) */
     691             : static GEN
     692         399 : subgroup_info(GEN bnfz, GEN Lprz, long ell, GEN vecWA)
     693             : {
     694         399 :   GEN nfz = bnf_get_nf(bnfz), M, gell = utoipos(ell), Lell = mkvec(gell);
     695         399 :   long i, j, l = lg(vecWA), lz = lg(Lprz);
     696         399 :   M = cgetg(l, t_MAT);
     697         399 :   for (j=1; j<l; j++) gel(M,j) = cgetg(lz, t_VECSMALL);
     698        1232 :   for (i=1; i < lz; i++)
     699             :   {
     700         833 :     GEN pr = gel(Lprz,i), EX = subiu(pr_norm(pr), 1);
     701         833 :     GEN N, g,T,p, prM = idealhnf(nfz, pr);
     702         833 :     GEN modpr = zk_to_Fq_init(nfz, &pr,&T,&p);
     703         833 :     long v = Z_lvalrem(divis(EX,ell), ell, &N) + 1; /* Norm(pr)-1 = N * ell^v */
     704         833 :     GEN ellv = powuu(ell, v);
     705         833 :     g = gener_Fq_local(T,p, Lell);
     706         833 :     g = Fq_pow(g,N, T,p); /* order ell^v */
     707        4186 :     for (j=1; j < l; j++)
     708             :     {
     709        3353 :       GEN logc, c = gel(vecWA,j);
     710        3353 :       if (typ(c) == t_MAT) /* famat */
     711        1071 :         c = famat_makecoprime(nfz, gel(c,1), gel(c,2), pr, prM, EX);
     712        3353 :       c = nf_to_Fq(nfz, c, modpr);
     713        3353 :       c = Fq_pow(c, N, T,p);
     714        3353 :       logc = Fq_log(c, g, ellv, T,p);
     715        3353 :       ucoeff(M, i,j) = umodiu(logc, ell);
     716             :     }
     717             :   }
     718         399 :   return M;
     719             : }
     720             : 
     721             : /* if all!=0, give all equations of degree 'all'. Assume bnr modulus is the
     722             :  * conductor */
     723             : static GEN
     724         182 : rnfkummersimple(GEN bnr, GEN subgroup, GEN gell, long all)
     725             : {
     726             :   long ell, i, j, degK, dK;
     727             :   long lSml2, lSl2, lSp, rc, lW;
     728             :   long prec;
     729         182 :   long rk=0, ncyc=0;
     730         182 :   GEN mat=NULL, matgrp=NULL, xell, be1 = NULL;
     731         182 :   long firstpass = all<0;
     732             : 
     733             :   GEN bnf,nf,bid,ideal,arch,cycgen;
     734             :   GEN cyc;
     735             :   GEN Sp,listprSp,matP;
     736         182 :   GEN res=NULL,u,M,K,y,vecMsup,vecW,vecWB,vecBp,msign;
     737             :   primlist L;
     738             : 
     739         182 :   bnf = bnr_get_bnf(bnr); (void)bnf_build_units(bnf);
     740         182 :   nf  = bnf_get_nf(bnf);
     741         182 :   degK = nf_get_degree(nf);
     742             : 
     743         182 :   bid = bnr_get_bid(bnr);
     744         182 :   ideal= bid_get_ideal(bid);
     745         182 :   arch = bid_get_arch(bid); /* this is the conductor */
     746         182 :   ell = itos(gell);
     747         182 :   i = build_list_Hecke(&L, nf, bid_get_fact2(bid), ideal, gell, NULL);
     748         182 :   if (i) return no_sol(all,i);
     749             : 
     750         182 :   lSml2 = lg(L.Sml2);
     751         182 :   Sp = shallowconcat(L.Sm, L.Sml1); lSp = lg(Sp);
     752         182 :   listprSp = shallowconcat(L.Sml2, L.Sl); lSl2 = lg(listprSp);
     753             : 
     754         182 :   cycgen = bnf_build_cycgen(bnf);
     755         182 :   cyc = bnf_get_cyc(bnf); rc = prank(cyc, ell);
     756             : 
     757         182 :   vecW = get_Selmer(bnf, cycgen, rc);
     758         182 :   u = get_u(cyc, rc, gell);
     759             : 
     760         182 :   vecBp = cgetg(lSp, t_VEC);
     761         182 :   matP  = cgetg(lSp, t_MAT);
     762         357 :   for (j = 1; j < lSp; j++)
     763             :   {
     764         175 :     GEN L = isprincipalell(bnf,gel(Sp,j), cycgen,u,gell,rc);
     765         175 :     gel( matP,j) = gel(L,1);
     766         175 :     gel(vecBp,j) = gel(L,2);
     767             :   }
     768         182 :   vecWB = shallowconcat(vecW, vecBp);
     769             : 
     770         182 :   prec = DEFAULTPREC +
     771         182 :       nbits2extraprec(((degK-1) * (gexpo(vecWB) + gexpo(nf_get_M(nf)))));
     772         182 :   if (nf_get_prec(nf) < prec) nf = nfnewprec_shallow(nf, prec);
     773         182 :   msign = nfsign(nf, vecWB);
     774         182 :   arch = ZV_to_zv(arch);
     775             : 
     776         182 :   vecMsup = cgetg(lSml2,t_VEC);
     777         182 :   M = NULL;
     778         343 :   for (i = 1; i < lSl2; i++)
     779             :   {
     780         161 :     GEN pr = gel(listprSp,i);
     781         161 :     long e = pr_get_e(pr), z = ell * (e / (ell-1));
     782             : 
     783         161 :     if (i < lSml2)
     784             :     {
     785          91 :       z += 1 - L.ESml2[i];
     786          91 :       gel(vecMsup,i) = logall(nf, vecWB, 0,0, ell, pr,z+1);
     787             :     }
     788         161 :     M = vconcat(M, logall(nf, vecWB, 0,0, ell, pr,z));
     789             :   }
     790         182 :   lW = lg(vecW);
     791         182 :   M = vconcat(M, shallowconcat(zero_Flm(rc,lW-1), ZM_to_Flm(matP, ell)));
     792         182 :   if (!all)
     793             :   { /* primes landing in subgroup must be totally split */
     794         182 :     GEN Lpr = get_prlist(bnr, subgroup, ell, NULL);
     795         182 :     GEN M2 = subgroup_info(bnf, Lpr, ell, vecWB);
     796         182 :     M = vconcat(M, M2);
     797             :   }
     798         182 :   K = Flm_ker(M, ell);
     799         182 :   dK = lg(K)-1;
     800             : 
     801         182 :   if (all < 0)
     802           0 :     K = fix_kernel(K, M, vecMsup, lW, ell);
     803             : 
     804         182 :   y = cgetg(dK+1,t_VECSMALL);
     805         182 :   if (all) res = cgetg(1,t_VEC); /* in case all = 1 */
     806         182 :   if (all < 0)
     807             :   {
     808           0 :     ncyc = dK;
     809           0 :     mat = Flm_init(dK, ncyc);
     810           0 :     if (all == -1) matgrp = Flm_init(lg(bnr_get_cyc(bnr)), ncyc+1);
     811           0 :     rk = 0;
     812             :   }
     813         182 :   xell = pol_xn(ell, 0);
     814             :   do {
     815         182 :     dK = lg(K)-1;
     816         371 :     while (dK)
     817             :     {
     818         189 :       for (i=1; i<dK; i++) y[i] = 0;
     819         189 :       y[i] = 1; /* y = [0,...,0,1,0,...,0], 1 at dK'th position */
     820             :       do
     821             :       {
     822         287 :         pari_sp av = avma;
     823         287 :         GEN be, P=NULL, X;
     824         287 :         if (all < 0)
     825             :         {
     826           0 :           Flv_fill(gel(mat, rk+1), y);
     827           0 :           setlg(mat, rk+2);
     828           0 :           if (Flm_rank(mat, ell) <= rk) continue;
     829             :         }
     830         287 : FOUND:  X = Flm_Flc_mul(K, y, ell);
     831         287 :         if (ok_congruence(X, ell, lW, vecMsup) && ok_sign(X, msign, arch))
     832             :         {/* be satisfies all congruences, x^ell - be is irreducible, signature
     833             :           * and relative discriminant are correct */
     834         182 :           if (all < 0) rk++;
     835         182 :           be = compute_beta(X, vecWB, gell, bnf);
     836         182 :           be = nf_to_scalar_or_alg(nf, be);
     837         182 :           if (typ(be) == t_POL) be = mkpolmod(be, nf_get_pol(nf));
     838         182 :           if (all == -1)
     839             :           {
     840           0 :             pari_sp av2 = avma;
     841           0 :             GEN Kgrp, colgrp = grptocol(rnfnormgroup(bnr, gsub(xell, be)));
     842           0 :             if (ell != 2)
     843             :             {
     844           0 :               if (rk == 1) be1 = be;
     845             :               else
     846             :               { /* Compute the pesky scalar */
     847           0 :                 GEN K2, C = cgetg(4, t_MAT);
     848           0 :                 gel(C,1) = gel(matgrp,1);
     849           0 :                 gel(C,2) = colgrp;
     850           0 :                 gel(C,3) = grptocol(rnfnormgroup(bnr, gsub(xell, gmul(be1,be))));
     851           0 :                 K2 = Flm_ker(C, ell);
     852           0 :                 if (lg(K2) != 2) pari_err_BUG("linear algebra");
     853           0 :                 K2 = gel(K2,1);
     854           0 :                 if (K2[1] != K2[2])
     855           0 :                   Flv_Fl_mul_inplace(colgrp, Fl_div(K2[2],K2[1],ell), ell);
     856             :               }
     857             :             }
     858           0 :             Flv_fill(gel(matgrp,rk), colgrp);
     859           0 :             setlg(matgrp, rk+1);
     860           0 :             Kgrp = Flm_ker(matgrp, ell);
     861           0 :             if (lg(Kgrp) == 2)
     862             :             {
     863           0 :               setlg(gel(Kgrp,1), rk+1);
     864           0 :               y = Flm_Flc_mul(mat, gel(Kgrp,1), ell);
     865           0 :               all = 0; goto FOUND;
     866             :             }
     867           0 :             avma = av2;
     868             :           }
     869             :           else
     870             :           {
     871         182 :             P = gsub(xell, be);
     872         182 :             if (all)
     873           0 :               res = shallowconcat(res, gerepileupto(av, P));
     874             :             else
     875             :             {
     876         182 :               if (ZM_equal(rnfnormgroup(bnr,P),subgroup)) return P; /*DONE*/
     877           0 :               avma = av; continue;
     878             :             }
     879             :           }
     880           0 :           if (all < 0 && rk == ncyc) return res;
     881           0 :           if (firstpass) break;
     882             :         }
     883         105 :         else avma = av;
     884         105 :       } while (increment(y, dK, ell));
     885           7 :       y[dK--] = 0;
     886             :     }
     887           0 :   } while (firstpass--);
     888           0 :   return all? res: gen_0;
     889             : }
     890             : 
     891             : /* alg. 5.3.11 (return only discrete log mod ell) */
     892             : static GEN
     893         826 : isvirtualunit(GEN bnf, GEN v, GEN cycgen, GEN cyc, GEN gell, long rc)
     894             : {
     895         826 :   GEN L, b, eps, y, q, nf = bnf_get_nf(bnf);
     896         826 :   GEN w = idealsqrtn(nf, v, gell, 1);
     897         826 :   long i, l = lg(cycgen);
     898             : 
     899         826 :   L = bnfisprincipal0(bnf, w, nf_GENMAT|nf_FORCE);
     900         826 :   q = gel(L,1);
     901         826 :   if (ZV_equal0(q)) { eps = v; y = q; }
     902             :   else
     903             :   {
     904         140 :     y = cgetg(l,t_COL);
     905         140 :     for (i=1; i<l; i++) gel(y,i) = diviiexact(mulii(gell,gel(q,i)), gel(cyc,i));
     906         140 :     eps = famat_mulpow_shallow(famat_factorback(cycgen,y), gel(L,2), gell);
     907         140 :     eps = famat_mul_shallow(famat_inv(eps), v);
     908             :   }
     909         826 :   setlg(y, rc+1);
     910         826 :   b = bnfisunit(bnf,eps);
     911         826 :   if (lg(b) == 1) pari_err_BUG("isvirtualunit");
     912         826 :   return shallowconcat(lift_shallow(b), y);
     913             : }
     914             : 
     915             : /* J a vector of elements in nfz = relative extension of nf by polrel,
     916             :  * return the Steinitz element attached to the module generated by J */
     917             : static GEN
     918         651 : Stelt(GEN nf, GEN J, GEN polrel)
     919             : {
     920         651 :   long i, l = lg(J), vx = varn(polrel);
     921         651 :   GEN A = cgetg(l, t_VEC), I = cgetg(l, t_VEC);
     922        4487 :   for (i = 1; i < l; i++)
     923             :   {
     924        3836 :     GEN v = gel(J,i);
     925        3836 :     if (typ(v) == t_POL) { v = RgX_rem(v, polrel); setvarn(v,vx); }
     926        3836 :     gel(A,i) = v;
     927        3836 :     gel(I,i) = gen_1;
     928             :   }
     929         651 :   A = RgV_to_RgM(A, degpol(polrel));
     930         651 :   return idealprod(nf, gel(nfhnf(nf, mkvec2(A,I)),2));
     931             : }
     932             : 
     933             : static GEN
     934         126 : polrelKzK(toK_s *T, GEN x)
     935             : {
     936         126 :   GEN P = roots_to_pol(powtau(x, T->m, T->tau), 0);
     937         126 :   long i, l = lg(P);
     938         126 :   for (i=2; i<l; i++) gel(P,i) = downtoK(T, gel(P,i));
     939         126 :   return P;
     940             : }
     941             : 
     942             : /* N: Cl_m(Kz) --> Cl_m(K), lift subgroup from bnr to bnrz using Algo 4.1.11 */
     943             : static GEN
     944         126 : invimsubgroup(GEN bnrz, GEN bnr, GEN subgroup, toK_s *T)
     945             : {
     946             :   long l, j;
     947             :   GEN P, cyc, gen, U, polrel, StZk;
     948         126 :   GEN nf = bnr_get_nf(bnr), nfz = bnr_get_nf(bnrz);
     949         126 :   GEN polz = nf_get_pol(nfz), zkzD = nf_get_zkprimpart(nfz);
     950             : 
     951         126 :   polrel = polrelKzK(T, pol_x(varn(polz)));
     952         126 :   StZk = Stelt(nf, zkzD, polrel);
     953         126 :   cyc = bnr_get_cyc(bnrz); l = lg(cyc);
     954         126 :   gen = bnr_get_gen(bnrz);
     955         126 :   P = cgetg(l,t_MAT);
     956         651 :   for (j=1; j<l; j++)
     957             :   {
     958         525 :     GEN g, id = idealhnf_shallow(nfz, gel(gen,j));
     959         525 :     g = Stelt(nf, RgV_RgM_mul(zkzD, id), polrel);
     960         525 :     g = idealdiv(nf, g, StZk); /* N_{Kz/K}(gen[j]) */
     961         525 :     gel(P,j) = isprincipalray(bnr, g);
     962             :   }
     963         126 :   (void)ZM_hnfall_i(shallowconcat(P, subgroup), &U, 1);
     964         126 :   setlg(U, l); for (j=1; j<l; j++) setlg(U[j], l);
     965         126 :   return ZM_hnfmodid(U, cyc);
     966             : }
     967             : 
     968             : static GEN
     969         252 : pol_from_Newton(GEN S)
     970             : {
     971         252 :   long i, k, l = lg(S);
     972         252 :   GEN C = cgetg(l+1, t_VEC), c = C + 1;
     973         252 :   gel(c,0) = gen_1;
     974         252 :   gel(c,1) = gel(S,1); /* gen_0 in our case */
     975         882 :   for (k = 2; k < l; k++)
     976             :   {
     977         630 :     GEN s = gel(S,k);
     978         630 :     for (i = 2; i < k-1; i++) s = gadd(s, gmul(gel(S,i), gel(c,k-i)));
     979         630 :     gel(c,k) = gdivgs(s, -k);
     980             :   }
     981         252 :   return gtopoly(C, 0);
     982             : }
     983             : 
     984             : /* - mu_b = sum_{0 <= i < m} floor(r_b r_{d-1-i} / ell) tau^i */
     985             : static GEN
     986         602 : get_mmu(long b, GEN r, long ell)
     987             : {
     988         602 :   long i, m = lg(r)-1;
     989         602 :   GEN M = cgetg(m+1, t_VEC);
     990         602 :   for (i = 0; i < m; i++) gel(M,i+1) = stoi((r[b + 1] * r[m - i]) / ell);
     991         602 :   return M;
     992             : }
     993             : 
     994             : /* coeffs(x, a..b) in variable v >= varn(x) */
     995             : static GEN
     996        5964 : split_pol(GEN x, long v, long a, long b)
     997             : {
     998        5964 :   long i, l = degpol(x);
     999        5964 :   GEN y = x + a, z;
    1000             : 
    1001        5964 :   if (l < b) b = l;
    1002        5964 :   if (a > b || varn(x) != v) return pol_0(v);
    1003        5320 :   l = b-a + 3;
    1004        5320 :   z = cgetg(l, t_POL); z[1] = x[1];
    1005        5320 :   for (i = 2; i < l; i++) gel(z,i) = gel(y,i);
    1006        5320 :   return normalizepol_lg(z, l);
    1007             : }
    1008             : 
    1009             : /* return (den_a * z) mod (v^ell - num_a/den_a), assuming deg(z) < 2*ell
    1010             :  * allow either num/den to be NULL (= 1) */
    1011             : static GEN
    1012        2982 : mod_Xell_a(GEN z, long v, long ell, GEN num_a, GEN den_a)
    1013             : {
    1014        2982 :   GEN z1 = split_pol(z, v, ell, degpol(z));
    1015        2982 :   GEN z0 = split_pol(z, v, 0,   ell-1); /* z = v^ell z1 + z0*/
    1016        2982 :   if (den_a) z0 = gmul(den_a, z0);
    1017        2982 :   if (num_a) z1 = gmul(num_a, z1);
    1018        2982 :   return gadd(z0, z1);
    1019             : }
    1020             : static GEN
    1021         854 : to_alg(GEN nfz, GEN c, long v)
    1022             : {
    1023             :   GEN z, D;
    1024         854 :   if (typ(c) != t_COL) return c;
    1025         602 :   z = gmul(nf_get_zkprimpart(nfz), c);
    1026         602 :   if (typ(z) == t_POL) setvarn(z, v);
    1027         602 :   D = nf_get_zkden(nfz);
    1028         602 :   if (!equali1(D)) z = RgX_Rg_div(z, D);
    1029         602 :   return z;
    1030             : }
    1031             : 
    1032             : /* th. 5.3.5. and prop. 5.3.9. */
    1033             : static GEN
    1034         252 : compute_polrel(GEN nfz, toK_s *T, GEN be, long g, long ell)
    1035             : {
    1036         252 :   long i, k, m = T->m, vT = fetch_var(), vz = fetch_var();
    1037             :   GEN r, powtaubet, S, p1, root, num_t, den_t, nfzpol, powtau_prim_invbe;
    1038             :   GEN prim_Rk, C_Rk, prim_root, C_root, prim_invbe, C_invbe;
    1039             :   pari_timer ti;
    1040             : 
    1041         252 :   r = cgetg(m+1,t_VECSMALL); /* r[i+1] = g^i mod ell */
    1042         252 :   r[1] = 1;
    1043         252 :   for (i=2; i<=m; i++) r[i] = (r[i-1] * g) % ell;
    1044         252 :   powtaubet = powtau(be, m, T->tau);
    1045         252 :   if (DEBUGLEVEL>1) { err_printf("Computing Newton sums: "); timer_start(&ti); }
    1046         252 :   prim_invbe = Q_primitive_part(nfinv(nfz, be), &C_invbe);
    1047         252 :   powtau_prim_invbe = powtau(prim_invbe, m, T->tau);
    1048             : 
    1049         252 :   root = cgetg(ell + 2, t_POL);
    1050         252 :   root[1] = evalsigne(1) | evalvarn(0);
    1051         252 :   for (i = 0; i < ell; i++) gel(root,2+i) = gen_0;
    1052         854 :   for (i = 0; i < m; i++)
    1053             :   { /* compute (1/be) ^ (-mu) instead of be^mu [mu << 0].
    1054             :      * 1/be = C_invbe * prim_invbe */
    1055         602 :     GEN mmu = get_mmu(i, r, ell);
    1056             :     /* p1 = prim_invbe ^ -mu */
    1057         602 :     p1 = to_alg(nfz, nffactorback(nfz, powtau_prim_invbe, mmu), vz);
    1058         602 :     if (C_invbe) p1 = gmul(p1, powgi(C_invbe, RgV_sumpart(mmu, m)));
    1059             :     /* root += zeta_ell^{r_i} T^{r_i} be^mu_i */
    1060         602 :     gel(root, 2 + r[i+1]) = monomial(p1, r[i+1], vT);
    1061             :   }
    1062             :   /* Other roots are as above with z_ell --> z_ell^j.
    1063             :    * Treat all contents (C_*) and principal parts (prim_*) separately */
    1064         252 :   prim_Rk = prim_root = Q_primitive_part(root, &C_root);
    1065         252 :   C_Rk = C_root;
    1066             : 
    1067         252 :   r = vecsmall_reverse(r); /* theta^ell = be^( sum tau^a r_{d-1-a} ) */
    1068             :   /* Compute modulo X^ell - 1, T^ell - t, nfzpol(vz) */
    1069         252 :   p1 = to_alg(nfz, nffactorback(nfz, powtaubet, r), vz);
    1070         252 :   num_t = Q_remove_denom(p1, &den_t);
    1071             : 
    1072         252 :   nfzpol = leafcopy(nf_get_pol(nfz));
    1073         252 :   setvarn(nfzpol, vz);
    1074         252 :   S = cgetg(ell+1, t_VEC); /* Newton sums */
    1075         252 :   gel(S,1) = gen_0;
    1076         882 :   for (k = 2; k <= ell; k++)
    1077             :   { /* compute the k-th Newton sum */
    1078         630 :     pari_sp av = avma;
    1079         630 :     GEN z, D, Rk = gmul(prim_Rk, prim_root);
    1080         630 :     C_Rk = mul_content(C_Rk, C_root);
    1081         630 :     Rk = mod_Xell_a(Rk, 0, ell, NULL, NULL); /* mod X^ell - 1 */
    1082        3010 :     for (i = 2; i < lg(Rk); i++)
    1083             :     {
    1084        2380 :       if (typ(gel(Rk,i)) != t_POL) continue;
    1085        2352 :       z = mod_Xell_a(gel(Rk,i), vT, ell, num_t,den_t); /* mod T^ell - t */
    1086        2352 :       gel(Rk,i) = RgXQX_red(z, nfzpol); /* mod nfz.pol */
    1087             :     }
    1088         630 :     if (den_t) C_Rk = mul_content(C_Rk, ginv(den_t));
    1089         630 :     prim_Rk = Q_primitive_part(Rk, &D);
    1090         630 :     C_Rk = mul_content(C_Rk, D); /* root^k = prim_Rk * C_Rk */
    1091             : 
    1092             :     /* Newton sum is ell * constant coeff (in X), which has degree 0 in T */
    1093         630 :     z = polcoeff_i(prim_Rk, 0, 0);
    1094         630 :     z = polcoeff_i(z      , 0,vT);
    1095         630 :     z = downtoK(T, gmulgs(z, ell));
    1096         630 :     if (C_Rk) z = gmul(z, C_Rk);
    1097         630 :     gerepileall(av, C_Rk? 3: 2, &z, &prim_Rk, &C_Rk);
    1098         630 :     if (DEBUGLEVEL>1) { err_printf("%ld(%ld) ", k, timer_delay(&ti)); err_flush(); }
    1099         630 :     gel(S,k) = z;
    1100             :   }
    1101         252 :   if (DEBUGLEVEL>1) err_printf("\n");
    1102         252 :   (void)delete_var();
    1103         252 :   (void)delete_var(); return pol_from_Newton(S);
    1104             : }
    1105             : 
    1106             : /* lift elt t in nf to nfz, algebraic form */
    1107             : static GEN
    1108         343 : lifttoKz(GEN nf, GEN t, compo_s *C)
    1109             : {
    1110         343 :   GEN x = nf_to_scalar_or_alg(nf, t);
    1111         343 :   if (typ(x) != t_POL) return x;
    1112         343 :   return RgX_RgXQ_eval(x, C->p, C->R);
    1113             : }
    1114             : /* lift ideal id in nf to nfz */
    1115             : static GEN
    1116         231 : ideallifttoKz(GEN nfz, GEN nf, GEN id, compo_s *C)
    1117             : {
    1118         231 :   GEN I = idealtwoelt(nf,id);
    1119         231 :   GEN x = nf_to_scalar_or_alg(nf, gel(I,2));
    1120         231 :   if (typ(x) != t_POL) return gel(I,1);
    1121         147 :   gel(I,2) = algtobasis(nfz, RgX_RgXQ_eval(x, C->p, C->R));
    1122         147 :   return idealhnf_two(nfz,I);
    1123             : }
    1124             : /* lift ideal pr in nf to ONE prime in nfz (the others are conjugate under tau
    1125             :  * and bring no further information on e_1 W). Assume pr coprime to
    1126             :  * index of both nf and nfz, and unramified in Kz/K (minor simplification) */
    1127             : static GEN
    1128         378 : prlifttoKz(GEN nfz, GEN nf, GEN pr, compo_s *C)
    1129             : {
    1130         378 :   GEN F, p = pr_get_p(pr), t = pr_get_gen(pr), T = nf_get_pol(nfz);
    1131         378 :   if (nf_get_degree(nf) != 1)
    1132             :   { /* restrict to primes above pr */
    1133         343 :     t = Q_primpart( lifttoKz(nf,t,C) );
    1134         343 :     T = FpX_gcd(FpX_red(T,p), FpX_red(t,p), p);
    1135         343 :     T = FpX_normalize(T, p);
    1136             :   }
    1137         378 :   F = FpX_factor(T, p);
    1138         378 :   return idealprimedec_kummer(nfz,gcoeff(F,1,1), pr_get_e(pr), p);
    1139             : }
    1140             : static GEN
    1141         217 : get_przlist(GEN L, GEN nfz, GEN nf, compo_s *C)
    1142             : {
    1143             :   long i, l;
    1144         217 :   GEN M = cgetg_copy(L, &l);
    1145         217 :   for (i = 1; i < l; i++) gel(M,i) = prlifttoKz(nfz, nf, gel(L,i), C);
    1146         217 :   return M;
    1147             : }
    1148             : 
    1149             : static void
    1150         231 : compositum_red(compo_s *C, GEN P, GEN Q)
    1151             : {
    1152         231 :   GEN p, q, a, z = gel(compositum2(P, Q),1);
    1153         231 :   a = gel(z,1);
    1154         231 :   p = gel(gel(z,2), 2);
    1155         231 :   q = gel(gel(z,3), 2);
    1156         231 :   C->k = itos( gel(z,4) );
    1157             :   /* reduce R. FIXME: should be polredbest(a, 1), but breaks rnfkummer bench */
    1158         231 :   z = polredabs0(a, nf_ORIG|nf_PARTIALFACT);
    1159         231 :   C->R = gel(z,1);
    1160         231 :   a = gel(gel(z,2), 2);
    1161         231 :   C->p = RgX_RgXQ_eval(p, a, C->R);
    1162         231 :   C->q = RgX_RgXQ_eval(q, a, C->R);
    1163         231 :   C->rev = QXQ_reverse(a, C->R);
    1164         231 :   if (DEBUGLEVEL>1) err_printf("polred(compositum) = %Ps\n",C->R);
    1165         231 : }
    1166             : 
    1167             : /* replace P->C^(-deg P) P(xC) for the largest integer C such that coefficients
    1168             :  * remain algebraic integers. Lift *rational* coefficients */
    1169             : static void
    1170         252 : nfX_Z_normalize(GEN nf, GEN P)
    1171             : {
    1172             :   long i, l;
    1173         252 :   GEN C, Cj, PZ = cgetg_copy(P, &l);
    1174         252 :   PZ[1] = P[1];
    1175        1386 :   for (i = 2; i < l; i++) /* minor variation on RgX_to_nfX (create PZ) */
    1176             :   {
    1177        1134 :     GEN z = nf_to_scalar_or_basis(nf, gel(P,i));
    1178        1134 :     if (typ(z) == t_INT)
    1179         735 :       gel(PZ,i) = gel(P,i) = z;
    1180             :     else
    1181         399 :       gel(PZ,i) = ZV_content(z);
    1182             :   }
    1183         252 :   (void)ZX_Z_normalize(PZ, &C);
    1184             : 
    1185         504 :   if (C == gen_1) return;
    1186          77 :   Cj = C;
    1187         322 :   for (i = l-2; i > 1; i--)
    1188             :   {
    1189         245 :     if (i != l-2) Cj = mulii(Cj, C);
    1190         245 :     gel(P,i) = gdiv(gel(P,i), Cj);
    1191             :   }
    1192             : }
    1193             : 
    1194             : /* v[i] = g^i mod ell, 1 <= i < m */
    1195             : static GEN
    1196         231 : Fl_powers_FpV(ulong g, long m, ulong ell)
    1197             : {
    1198         231 :   GEN v = cgetg(m, t_VEC);
    1199         231 :   ulong gi = g % ell;
    1200             :   long i;
    1201         231 :   for (i=1; i<m; i++) { gel(v,i) = utoi(gi); gi = (gi*g) % ell; }
    1202         231 :   return v;
    1203             : }
    1204             : 
    1205             : static GEN
    1206         420 : _rnfkummer(GEN bnr, GEN subgroup, long all, long prec)
    1207             : {
    1208             :   long ell, i, j, m, d, dK, dc, rc, ru, rv, g, mginv, degK, degKz, vnf;
    1209             :   long lSp, lSml2, lSl2, lW;
    1210             :   GEN polnf,bnf,nf,bnfz,nfz,bid,ideal,cycgen,gell,p1,p2,vselmer;
    1211             :   GEN cyc,gen;
    1212             :   GEN Q,idealz,gothf;
    1213         420 :   GEN res=NULL,u,M,K,y,vecMsup,vecW,vecWA,vecWB,vecB,vecC,vecAp,vecBp;
    1214             :   GEN matP, Sp, listprSp, Tc, Tv, P;
    1215             :   primlist L;
    1216             :   toK_s T;
    1217             :   tau_s _tau, *tau;
    1218             :   compo_s COMPO;
    1219             :   pari_timer t;
    1220         420 :   long rk=0, ncyc=0;
    1221         420 :   GEN mat = NULL;
    1222         420 :   long firstpass = all<0;
    1223             : 
    1224         420 :   if (DEBUGLEVEL) timer_start(&t);
    1225         420 :   checkbnr(bnr);
    1226         420 :   bnf = bnr_get_bnf(bnr);
    1227         420 :   nf  = bnf_get_nf(bnf);
    1228         420 :   polnf = nf_get_pol(nf); vnf = varn(polnf);
    1229         420 :   if (!vnf) pari_err_PRIORITY("rnfkummer", polnf, "=", 0);
    1230             :   /* step 7 */
    1231         420 :   p1 = bnrconductor_i(bnr, subgroup, 2);
    1232         420 :   if (DEBUGLEVEL) timer_printf(&t, "[rnfkummer] conductor");
    1233         420 :   bnr      = gel(p1,2);
    1234         420 :   subgroup = gel(p1,3);
    1235         420 :   gell = get_gell(bnr,subgroup,all);
    1236         420 :   ell = itos(gell);
    1237         420 :   if (ell == 1) return pol_x(0);
    1238         420 :   if (!uisprime(ell)) pari_err_IMPL("kummer for composite relative degree");
    1239         420 :   if (all && all != -1 && umodiu(bnr_get_no(bnr), ell))
    1240           7 :     return cgetg(1, t_VEC);
    1241         413 :   if (bnf_get_tuN(bnf) % ell == 0)
    1242         182 :     return rnfkummersimple(bnr, subgroup, gell, all);
    1243             : 
    1244         231 :   if (all == -1) all = 0;
    1245         231 :   bid = bnr_get_bid(bnr);
    1246         231 :   ideal = bid_get_ideal(bid);
    1247             :   /* step 1 of alg 5.3.5. */
    1248         231 :   if (DEBUGLEVEL>2) err_printf("Step 1\n");
    1249         231 :   compositum_red(&COMPO, polnf, polcyclo(ell,vnf));
    1250             :   /* step 2 */
    1251         231 :   if (DEBUGLEVEL>2) err_printf("Step 2\n");
    1252         231 :   if (DEBUGLEVEL) timer_printf(&t, "[rnfkummer] compositum");
    1253         231 :   degK  = degpol(polnf);
    1254         231 :   degKz = degpol(COMPO.R);
    1255         231 :   m = degKz / degK;
    1256         231 :   d = (ell-1) / m;
    1257         231 :   g = (long)Fl_powu(pgener_Fl(ell), d, ell);
    1258         231 :   if (Fl_powu((ulong)g, m, ell*ell) == 1) g += ell;
    1259             :   /* ord(g) = m in all (Z/ell^k)^* */
    1260             :   /* step 3 */
    1261         231 :   if (DEBUGLEVEL>2) err_printf("Step 3\n");
    1262             :   /* could factor disc(R) using th. 2.1.6. */
    1263         231 :   bnfz = Buchall(COMPO.R, nf_FORCE, maxss(prec,BIGDEFAULTPREC));
    1264         231 :   if (DEBUGLEVEL) timer_printf(&t, "[rnfkummer] bnfinit(Kz)");
    1265         231 :   cycgen = bnf_build_cycgen(bnfz);
    1266         231 :   nfz = bnf_get_nf(bnfz);
    1267         231 :   cyc = bnf_get_cyc(bnfz); rc = prank(cyc,ell);
    1268         231 :   gen = bnf_get_gen(bnfz);
    1269         231 :   u = get_u(cyc, rc, gell);
    1270             : 
    1271         231 :   vselmer = get_Selmer(bnfz, cycgen, rc);
    1272         231 :   if (DEBUGLEVEL) timer_printf(&t, "[rnfkummer] Selmer group");
    1273         231 :   ru = (degKz>>1)-1;
    1274         231 :   rv = rc+ru+1;
    1275         231 :   tau = get_tau(&_tau, nfz, &COMPO, g);
    1276             : 
    1277             :   /* step 4 */
    1278         231 :   if (DEBUGLEVEL>2) err_printf("Step 4\n");
    1279         231 :   vecB=cgetg(rc+1,t_VEC);
    1280         231 :   Tc=cgetg(rc+1,t_MAT);
    1281         371 :   for (j=1; j<=rc; j++)
    1282             :   {
    1283         140 :     p1 = tauofideal(gel(gen,j), tau);
    1284         140 :     p1 = isprincipalell(bnfz, p1, cycgen,u,gell,rc);
    1285         140 :     gel(Tc,j)  = gel(p1,1);
    1286         140 :     gel(vecB,j)= gel(p1,2);
    1287             :   }
    1288             : 
    1289         231 :   vecC = cgetg(rc+1,t_VEC);
    1290         231 :   if (rc)
    1291             :   {
    1292         126 :     for (j=1; j<=rc; j++) gel(vecC,j) = cgetg(1, t_MAT);
    1293         126 :     p1 = cgetg(m,t_VEC);
    1294         126 :     gel(p1,1) = matid(rc);
    1295         126 :     for (j=2; j<=m-1; j++) gel(p1,j) = gmul(gel(p1,j-1),Tc);
    1296         126 :     p2 = vecB;
    1297         294 :     for (j=1; j<=m-1; j++)
    1298             :     {
    1299         168 :       GEN z = FpM_red(gmulsg((j*d)%ell,gel(p1,m-j)), gell);
    1300         168 :       p2 = tauofvec(p2, tau);
    1301         364 :       for (i=1; i<=rc; i++)
    1302         392 :         gel(vecC,i) = famat_mul_shallow(gel(vecC,i),
    1303         196 :                                         famat_factorback(p2,gel(z,i)));
    1304             :     }
    1305         126 :     for (i=1; i<=rc; i++) gel(vecC,i) = famat_reduce(gel(vecC,i));
    1306             :   }
    1307             :   /* step 5 */
    1308         231 :   if (DEBUGLEVEL>2) err_printf("Step 5\n");
    1309         231 :   Tv = cgetg(rv+1,t_MAT);
    1310        1057 :   for (j=1; j<=rv; j++)
    1311             :   {
    1312         826 :     p1 = tauofelt(gel(vselmer,j), tau);
    1313         826 :     if (typ(p1) == t_MAT) /* famat */
    1314         140 :       p1 = nffactorback(nfz, gel(p1,1), FpC_red(gel(p1,2),gell));
    1315         826 :     gel(Tv,j) = isvirtualunit(bnfz, p1, cycgen,cyc,gell,rc);
    1316             :   }
    1317         231 :   P = FpM_ker(RgM_Rg_add_shallow(Tv, stoi(-g)), gell);
    1318         231 :   lW = lg(P); vecW = cgetg(lW,t_VEC);
    1319         231 :   for (j=1; j<lW; j++) gel(vecW,j) = famat_factorback(vselmer, gel(P,j));
    1320             :   /* step 6 */
    1321         231 :   if (DEBUGLEVEL>2) err_printf("Step 6\n");
    1322         231 :   Q = FpM_ker(RgM_Rg_add_shallow(shallowtrans(Tc), stoi(-g)), gell);
    1323             :   /* step 8 */
    1324         231 :   if (DEBUGLEVEL>2) err_printf("Step 8\n");
    1325         231 :   p1 = RgXQ_matrix_pow(COMPO.p, degKz, degK, COMPO.R);
    1326         231 :   T.invexpoteta1 = RgM_inv(p1); /* left inverse */
    1327         231 :   T.polnf = polnf;
    1328         231 :   T.tau = tau;
    1329         231 :   T.m = m;
    1330         231 :   T.powg = Fl_powers_FpV(g, m, ell);
    1331             : 
    1332         231 :   idealz = ideallifttoKz(nfz, nf, ideal, &COMPO);
    1333         231 :   if (umodiu(gcoeff(ideal,1,1), ell)) gothf = idealz;
    1334             :   else
    1335             :   { /* ell | N(ideal) */
    1336         126 :     GEN bnrz = Buchray(bnfz, idealz, nf_INIT|nf_GEN);
    1337         126 :     GEN subgroupz = invimsubgroup(bnrz, bnr, subgroup, &T);
    1338         126 :     gothf = bnrconductor_i(bnrz,subgroupz,0);
    1339             :   }
    1340             :   /* step 9, 10, 11 */
    1341         231 :   if (DEBUGLEVEL>2) err_printf("Step 9, 10 and 11\n");
    1342         231 :   i = build_list_Hecke(&L, nfz, NULL, gothf, gell, tau);
    1343         231 :   if (i) return no_sol(all,i);
    1344             : 
    1345         231 :   lSml2 = lg(L.Sml2);
    1346         231 :   Sp = shallowconcat(L.Sm, L.Sml1); lSp = lg(Sp);
    1347         231 :   listprSp = shallowconcat(L.Sml2, L.Sl); lSl2 = lg(listprSp);
    1348             : 
    1349             :   /* step 12 */
    1350         231 :   if (DEBUGLEVEL>2) err_printf("Step 12\n");
    1351         231 :   vecAp = cgetg(lSp, t_VEC);
    1352         231 :   vecBp = cgetg(lSp, t_VEC);
    1353         231 :   matP  = cgetg(lSp, t_MAT);
    1354             : 
    1355         385 :   for (j = 1; j < lSp; j++)
    1356             :   {
    1357             :     GEN e, a;
    1358         154 :     p1 = isprincipalell(bnfz, gel(Sp,j), cycgen,u,gell,rc);
    1359         154 :     e = gel(p1,1); gel(matP,j) = e;
    1360         154 :     a = gel(p1,2);
    1361         154 :     gel(vecBp,j) = famat_mul_shallow(famat_factorback(vecC, gneg(e)), a);
    1362             :   }
    1363         231 :   vecAp = lambdaofvec(vecBp, &T);
    1364             :   /* step 13 */
    1365         231 :   if (DEBUGLEVEL>2) err_printf("Step 13\n");
    1366         231 :   vecWA = shallowconcat(vecW, vecAp);
    1367         231 :   vecWB = shallowconcat(vecW, vecBp);
    1368             : 
    1369             :   /* step 14, 15, and 17 */
    1370         231 :   if (DEBUGLEVEL>2) err_printf("Step 14, 15 and 17\n");
    1371         231 :   mginv = (m * Fl_inv(g,ell)) % ell;
    1372         231 :   vecMsup = cgetg(lSml2,t_VEC);
    1373         231 :   M = NULL;
    1374         518 :   for (i = 1; i < lSl2; i++)
    1375             :   {
    1376         287 :     GEN pr = gel(listprSp,i);
    1377         287 :     long e = pr_get_e(pr), z = ell * (e / (ell-1));
    1378             : 
    1379         287 :     if (i < lSml2)
    1380             :     {
    1381         133 :       z += 1 - L.ESml2[i];
    1382         133 :       gel(vecMsup,i) = logall(nfz, vecWA,lW,mginv,ell, pr,z+1);
    1383             :     }
    1384         287 :     M = vconcat(M, logall(nfz, vecWA,lW,mginv,ell, pr,z));
    1385             :   }
    1386         231 :   dc = lg(Q)-1;
    1387         231 :   if (dc)
    1388             :   {
    1389         105 :     GEN QtP = gmul(shallowtrans(Q), matP);
    1390         105 :     M = vconcat(M, shallowconcat(zero_Flm(dc,lW-1), ZM_to_Flm(QtP,ell)));
    1391             :   }
    1392         231 :   if (!M) M = zero_Flm(1, lSp-1 + lW-1);
    1393             : 
    1394         231 :   if (!all)
    1395             :   { /* primes landing in subgroup must be totally split */
    1396         217 :     GEN lambdaWB = shallowconcat(lambdaofvec(vecW, &T), vecAp);/*vecWB^lambda*/
    1397         217 :     GEN Lpr = get_prlist(bnr, subgroup, ell, bnfz);
    1398         217 :     GEN Lprz= get_przlist(Lpr, nfz, nf, &COMPO);
    1399         217 :     GEN M2 = subgroup_info(bnfz, Lprz, ell, lambdaWB);
    1400         217 :     M = vconcat(M, M2);
    1401             :   }
    1402             :   /* step 16 */
    1403         231 :   if (DEBUGLEVEL>2) err_printf("Step 16\n");
    1404         231 :   K = Flm_ker(M, ell);
    1405         231 :   if (all < 0)
    1406           0 :     K = fix_kernel(K, M, vecMsup, lW, ell);
    1407             :   /* step 18 & ff */
    1408         231 :   if (DEBUGLEVEL>2) err_printf("Step 18\n");
    1409         231 :   dK = lg(K)-1;
    1410         231 :   y = cgetg(dK+1,t_VECSMALL);
    1411         231 :   if (all) res = cgetg(1, t_VEC);
    1412         231 :   if (DEBUGLEVEL) timer_printf(&t, "[rnfkummer] candidate list");
    1413         231 :   if (all < 0) { ncyc = dK; rk = 0; mat = zero_Flm(lg(M)-1, ncyc); }
    1414             : 
    1415             :   do {
    1416         231 :     dK = lg(K)-1;
    1417         483 :     while (dK)
    1418             :     {
    1419         238 :       for (i=1; i<dK; i++) y[i] = 0;
    1420         238 :       y[i] = 1; /* y = [0,...,0,1,0,...,0], 1 at dK'th position */
    1421             :       do
    1422             :       { /* cf. algo 5.3.18 */
    1423         252 :         GEN H, be, P, X = Flm_Flc_mul(K, y, ell);
    1424         252 :         if (ok_congruence(X, ell, lW, vecMsup))
    1425             :         {
    1426         252 :           pari_sp av = avma;
    1427         252 :           if (all < 0)
    1428             :           {
    1429           0 :             gel(mat, rk+1) = X;
    1430           0 :             if (Flm_rank(mat,ell) <= rk) continue;
    1431           0 :             rk++;
    1432             :           }
    1433         252 :           be = compute_beta(X, vecWB, gell, bnfz);
    1434         252 :           P = compute_polrel(nfz, &T, be, g, ell);
    1435         252 :           nfX_Z_normalize(nf, P);
    1436         252 :           if (DEBUGLEVEL>1) err_printf("polrel(beta) = %Ps\n", P);
    1437         252 :           if (!all) {
    1438         217 :             H = rnfnormgroup(bnr, P);
    1439         217 :             if (ZM_equal(subgroup, H)) return P; /* DONE */
    1440           0 :             avma = av; continue;
    1441             :           } else {
    1442          35 :             GEN P0 = Q_primpart(lift_shallow(P));
    1443          35 :             GEN g = nfgcd(P0, RgX_deriv(P0), polnf, nf_get_index(nf));
    1444          35 :             if (degpol(g)) continue;
    1445          35 :             H = rnfnormgroup(bnr, P);
    1446          35 :             if (!ZM_equal(subgroup,H) && !bnrisconductor(bnr,H)) continue;
    1447             :           }
    1448          35 :           P = gerepilecopy(av, P);
    1449          35 :           res = shallowconcat(res, P);
    1450          35 :           if (all < 0 && rk == ncyc) return res;
    1451          35 :           if (firstpass) break;
    1452             :         }
    1453          35 :       } while (increment(y, dK, ell));
    1454          21 :       y[dK--] = 0;
    1455             :     }
    1456          14 :   } while (firstpass--);
    1457          14 :   if (all) return res;
    1458           0 :   return gen_0; /* FAIL */
    1459             : }
    1460             : 
    1461             : GEN
    1462         420 : rnfkummer(GEN bnr, GEN subgroup, long all, long prec)
    1463             : {
    1464         420 :   pari_sp av = avma;
    1465         420 :   return gerepilecopy(av, _rnfkummer(bnr, subgroup, all, prec));
    1466             : }

Generated by: LCOV version 1.11