Code coverage tests

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

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

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

LCOV - code coverage report
Current view: top level - basemath - buch3.c (source / functions) Hit Total Coverage
Test: PARI/GP v2.12.1 lcov report (development 25819-e703fe1174) Lines: 1513 1606 94.2 %
Date: 2020-09-18 06:10:04 Functions: 118 124 95.2 %
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             : /*                       RAY CLASS FIELDS                          */
      17             : /*                                                                 */
      18             : /*******************************************************************/
      19             : #include "pari.h"
      20             : #include "paripriv.h"
      21             : 
      22             : static GEN
      23      387631 : bnr_get_El(GEN bnr) { return gel(bnr,3); }
      24             : static GEN
      25      386104 : bnr_get_U(GEN bnr) { return gel(bnr,4); }
      26             : static GEN
      27        9107 : bnr_get_Ui(GEN bnr) { return gmael(bnr,4,3); }
      28             : 
      29             : /* faster than Buchray */
      30             : GEN
      31          35 : bnfnarrow(GEN bnf)
      32             : {
      33             :   GEN nf, cyc, gen, Cyc, Gen, A, GD, v, w, H, invpi, L, R, u, U0, Uoo, archp, sarch;
      34             :   long r1, j, l, t, RU;
      35             :   pari_sp av;
      36             : 
      37          35 :   bnf = checkbnf(bnf);
      38          35 :   nf = bnf_get_nf(bnf);
      39          35 :   r1 = nf_get_r1(nf); if (!r1) return gcopy( bnf_get_clgp(bnf) );
      40             : 
      41             :   /* simplified version of nfsign_units; r1 > 0 so bnf.tu = -1 */
      42          35 :   av = avma; archp = identity_perm(r1);
      43          35 :   A = bnf_get_logfu(bnf); RU = lg(A)+1;
      44          35 :   invpi = invr( mppi(nf_get_prec(nf)) );
      45          35 :   v = cgetg(RU,t_MAT); gel(v, 1) = const_vecsmall(r1, 1); /* nfsign(-1) */
      46          98 :   for (j=2; j<RU; j++) gel(v,j) = nfsign_from_logarch(gel(A,j-1), invpi, archp);
      47             :   /* up to here */
      48             : 
      49          35 :   v = Flm_image(v, 2); t = lg(v)-1;
      50          35 :   if (t == r1) { set_avma(av); return gcopy( bnf_get_clgp(bnf) ); }
      51             : 
      52          28 :   v = Flm_suppl(v,2); /* v = (sgn(U)|H) in GL_r1(F_2) */
      53          28 :   H = zm_to_ZM( vecslice(v, t+1, r1) ); /* supplement H of sgn(U) */
      54          28 :   w = rowslice(Flm_inv(v,2), t+1, r1); /* H*w*z = proj of z on H // sgn(U) */
      55             : 
      56          28 :   sarch = nfarchstar(nf, NULL, archp);
      57          28 :   cyc = bnf_get_cyc(bnf);
      58          28 :   gen = bnf_get_gen(bnf); l = lg(gen);
      59          28 :   L = cgetg(l,t_MAT); GD = gmael(bnf,9,3);
      60          63 :   for (j=1; j<l; j++)
      61             :   {
      62          35 :     GEN z = nfsign_from_logarch(gel(GD,j), invpi, archp);
      63          35 :     gel(L,j) = zc_to_ZC( Flm_Flc_mul(w, z, 2) );
      64             :   }
      65             :   /* [cyc, 0; L, 2] = relation matrix for Cl_f */
      66          28 :   R = shallowconcat(
      67             :     vconcat(diagonal_shallow(cyc), L),
      68             :     vconcat(zeromat(l-1, r1-t), scalarmat_shallow(gen_2,r1-t)));
      69          28 :   Cyc = ZM_snf_group(R, NULL, &u);
      70          28 :   U0 = rowslice(u, 1, l-1);
      71          28 :   Uoo = ZM_mul(H, rowslice(u, l, nbrows(u)));
      72          28 :   l = lg(Cyc); Gen = cgetg(l,t_VEC);
      73          91 :   for (j = 1; j < l; j++)
      74             :   {
      75          63 :     GEN g = gel(U0,j), s = gel(Uoo,j);
      76          63 :     g = (lg(g) == 1)? gen_1: Q_primpart( idealfactorback(nf,gen,g,0) );
      77          63 :     if (!ZV_equal0(s))
      78             :     {
      79          28 :       GEN a = set_sign_mod_divisor(nf, ZV_to_Flv(s,2), gen_1, sarch);
      80          28 :       g = is_pm1(g)? a: idealmul(nf, a, g);
      81             :     }
      82          63 :     gel(Gen,j) = g;
      83             :   }
      84          28 :   return gerepilecopy(av, mkvec3(shifti(bnf_get_no(bnf),r1-t), Cyc, Gen));
      85             : }
      86             : 
      87             : /********************************************************************/
      88             : /**                                                                **/
      89             : /**                  REDUCTION MOD IDELE                           **/
      90             : /**                                                                **/
      91             : /********************************************************************/
      92             : 
      93             : static GEN
      94       25830 : compute_fact(GEN nf, GEN U, GEN gen)
      95             : {
      96       25830 :   long i, j, l = lg(U), h = lgcols(U); /* l > 1 */
      97       25830 :   GEN basecl = cgetg(l,t_VEC), G;
      98             : 
      99       25830 :   G = mkvec2(NULL, trivial_fact());
     100       55573 :   for (j = 1; j < l; j++)
     101             :   {
     102       29743 :     GEN z = NULL;
     103       99540 :     for (i = 1; i < h; i++)
     104             :     {
     105       69797 :       GEN g, e = gcoeff(U,i,j); if (!signe(e)) continue;
     106             : 
     107       32151 :       g = gel(gen,i);
     108       32151 :       if (typ(g) != t_MAT)
     109             :       {
     110       21266 :         if (z)
     111        2247 :           gel(z,2) = famat_mulpow_shallow(gel(z,2), g, e);
     112             :         else
     113       19019 :           z = mkvec2(NULL, to_famat_shallow(g, e));
     114       21266 :         continue;
     115             :       }
     116       10885 :       gel(G,1) = g;
     117       10885 :       g = idealpowred(nf,G,e);
     118       10885 :       z = z? idealmulred(nf,z,g): g;
     119             :     }
     120       29743 :     gel(z,2) = famat_reduce(gel(z,2));
     121       29743 :     gel(basecl,j) = z;
     122             :   }
     123       25830 :   return basecl;
     124             : }
     125             : 
     126             : static int
     127       15281 : too_big(GEN nf, GEN bet)
     128             : {
     129       15281 :   GEN x = nfnorm(nf,bet);
     130       15281 :   switch (typ(x))
     131             :   {
     132        8988 :     case t_INT: return abscmpii(x, gen_1);
     133        6293 :     case t_FRAC: return abscmpii(gel(x,1), gel(x,2));
     134             :   }
     135           0 :   pari_err_BUG("wrong type in too_big");
     136             :   return 0; /* LCOV_EXCL_LINE */
     137             : }
     138             : 
     139             : /* true nf; GTM 193: Algo 4.3.4. Reduce x mod divisor */
     140             : static GEN
     141       14938 : idealmoddivisor_aux(GEN nf, GEN x, GEN f, GEN sarch)
     142             : {
     143       14938 :   pari_sp av = avma;
     144             :   GEN a, A;
     145             : 
     146       14938 :   if ( is_pm1(gcoeff(f,1,1)) ) /* f = 1 */
     147             :   {
     148         455 :     A = idealred(nf, mkvec2(x, gen_1));
     149         455 :     A = nfinv(nf, gel(A,2));
     150             :   }
     151             :   else
     152             :   {/* given coprime integral ideals x and f (f HNF), compute "small"
     153             :     * G in x, such that G = 1 mod (f). GTM 193: Algo 4.3.3 */
     154       14483 :     GEN G = idealaddtoone_raw(nf, x, f);
     155       14483 :     GEN D = idealaddtoone_i(nf, idealdiv(nf,G,x), f);
     156       14483 :     A = nfdiv(nf,D,G);
     157             :   }
     158       14938 :   if (too_big(nf,A) > 0) return gc_const(av, x);
     159       13461 :   a = set_sign_mod_divisor(nf, NULL, A, sarch);
     160       13461 :   if (a != A && too_big(nf,A) > 0) return gc_const(av, x);
     161       13461 :   return idealmul(nf, a, x);
     162             : }
     163             : 
     164             : GEN
     165        4214 : idealmoddivisor(GEN bnr, GEN x)
     166             : {
     167        4214 :   GEN nf = bnr_get_nf(bnr), bid = bnr_get_bid(bnr);
     168        4214 :   return idealmoddivisor_aux(nf, x, bid_get_ideal(bid), bid_get_sarch(bid));
     169             : }
     170             : 
     171             : /* v_pr(L0 * cx) */
     172             : static long
     173       10668 : fast_val(GEN L0, GEN cx, GEN pr)
     174             : {
     175       10668 :   pari_sp av = avma;
     176       10668 :   long v = typ(L0) == t_INT? 0: ZC_nfval(L0,pr);
     177       10668 :   if (cx)
     178             :   {
     179        9800 :     long w = Q_pval(cx, pr_get_p(pr));
     180        9800 :     if (w) v += w * pr_get_e(pr);
     181             :   }
     182       10668 :   return gc_long(av,v);
     183             : }
     184             : 
     185             : /* x coprime to fZ, return y = x mod fZ, y integral */
     186             : static GEN
     187        3703 : make_integral_Z(GEN x, GEN fZ)
     188             : {
     189        3703 :   GEN d, y = Q_remove_denom(x, &d);
     190        3703 :   if (d) y = FpC_Fp_mul(y, Fp_inv(d, fZ), fZ);
     191        3703 :   return y;
     192             : }
     193             : 
     194             : /* p pi^(-1) mod f */
     195             : static GEN
     196        3885 : get_pinvpi(GEN nf, GEN fZ, GEN p, GEN pi, GEN *v)
     197             : {
     198        3885 :   if (!*v) {
     199        3703 :     GEN invpi = nfinv(nf, pi);
     200        3703 :     *v = make_integral_Z(RgC_Rg_mul(invpi, p), mulii(p, fZ));
     201             :   }
     202        3885 :   return *v;
     203             : }
     204             : /* uniformizer pi for pr, coprime to F/p */
     205             : static GEN
     206        7868 : get_pi(GEN F, GEN pr, GEN *v)
     207             : {
     208        7868 :   if (!*v) *v = pr_uniformizer(pr, F);
     209        7868 :   return *v;
     210             : }
     211             : 
     212             : /* true nf */
     213             : static GEN
     214       32592 : bnr_grp(GEN nf, GEN U, GEN gen, GEN cyc, GEN bid)
     215             : {
     216       32592 :   GEN h = ZV_prod(cyc);
     217             :   GEN f, fZ, basecl, fa, pr, t, EX, sarch, F, P, vecpi, vecpinvpi;
     218             :   long i,j,l,lp;
     219             : 
     220       32592 :   if (lg(U) == 1) return mkvec3(h, cyc, cgetg(1, t_VEC));
     221       25830 :   basecl = compute_fact(nf, U, gen); /* generators in factored form */
     222       25830 :   EX = gel(bid_get_cyc(bid),1); /* exponent of (O/f)^* */
     223       25830 :   f  = bid_get_ideal(bid); fZ = gcoeff(f,1,1);
     224       25830 :   fa = bid_get_fact(bid);
     225       25830 :   sarch = bid_get_sarch(bid);
     226       25830 :   P = gel(fa,1); F = prV_lcm_capZ(P);
     227             : 
     228       25830 :   lp = lg(P);
     229       25830 :   vecpinvpi = cgetg(lp, t_VEC);
     230       25830 :   vecpi  = cgetg(lp, t_VEC);
     231       64617 :   for (i=1; i<lp; i++)
     232             :   {
     233       38787 :     pr = gel(P,i);
     234       38787 :     gel(vecpi,i)    = NULL; /* to be computed if needed */
     235       38787 :     gel(vecpinvpi,i) = NULL; /* to be computed if needed */
     236             :   }
     237             : 
     238       25830 :   l = lg(basecl);
     239       55573 :   for (i=1; i<l; i++)
     240             :   {
     241             :     GEN p, pi, pinvpi, dmulI, mulI, G, I, A, e, L, newL;
     242             :     long la, v, k;
     243             :     pari_sp av;
     244             :     /* G = [I, A=famat(L,e)] is a generator, I integral */
     245       29743 :     G = gel(basecl,i);
     246       29743 :     I = gel(G,1);
     247       29743 :     A = gel(G,2); L = gel(A,1); e = gel(A,2);
     248             :     /* if no reduction took place in compute_fact, everybody is still coprime
     249             :      * to f + no denominators */
     250       29743 :     if (!I) { gel(basecl,i) = famat_to_nf_moddivisor(nf, L, e, bid); continue; }
     251       10724 :     if (lg(A) == 1) { gel(basecl,i) = I; continue; }
     252             : 
     253             :     /* compute mulI so that mulI * I coprime to f
     254             :      * FIXME: use idealcoprime ??? (Less efficient. Fix idealcoprime!) */
     255       10724 :     dmulI = mulI = NULL;
     256       25816 :     for (j=1; j<lp; j++)
     257             :     {
     258       15092 :       pr = gel(P,j);
     259       15092 :       v  = idealval(nf, I, pr);
     260       15092 :       if (!v) continue;
     261        3766 :       p  = pr_get_p(pr);
     262        3766 :       pi = get_pi(F, pr, &gel(vecpi,j));
     263        3766 :       pinvpi = get_pinvpi(nf, fZ, p, pi, &gel(vecpinvpi,j));
     264        3766 :       t = nfpow_u(nf, pinvpi, (ulong)v);
     265        3766 :       mulI = mulI? nfmuli(nf, mulI, t): t;
     266        3766 :       t = powiu(p, v);
     267        3766 :       dmulI = dmulI? mulii(dmulI, t): t;
     268             :     }
     269             : 
     270             :     /* make all components of L coprime to f.
     271             :      * Assuming (L^e * I, f) = 1, then newL^e * mulI = L^e */
     272       10724 :     la = lg(e); newL = cgetg(la, t_VEC);
     273       17339 :     for (k=1; k<la; k++)
     274             :     {
     275        6615 :       GEN cx, LL = nf_to_scalar_or_basis(nf, gel(L,k));
     276        6615 :       GEN L0 = Q_primitive_part(LL, &cx); /* LL = L0*cx (faster nfval) */
     277       17283 :       for (j=1; j<lp; j++)
     278             :       {
     279       10668 :         pr = gel(P,j);
     280       10668 :         v  = fast_val(L0,cx, pr); /* = val_pr(LL) */
     281       10668 :         if (!v) continue;
     282        4102 :         p  = pr_get_p(pr);
     283        4102 :         pi = get_pi(F, pr, &gel(vecpi,j));
     284        4102 :         if (v > 0)
     285             :         {
     286         119 :           pinvpi = get_pinvpi(nf, fZ, p, pi, &gel(vecpinvpi,j));
     287         119 :           t = nfpow_u(nf,pinvpi, (ulong)v);
     288         119 :           LL = nfmul(nf, LL, t);
     289         119 :           LL = gdiv(LL, powiu(p, v));
     290             :         }
     291             :         else
     292             :         {
     293        3983 :           t = nfpow_u(nf,pi,(ulong)(-v));
     294        3983 :           LL = nfmul(nf, LL, t);
     295             :         }
     296             :       }
     297        6615 :       LL = make_integral(nf,LL,f,P);
     298        6615 :       gel(newL,k) = typ(LL) == t_INT? LL: FpC_red(LL, fZ);
     299             :     }
     300             : 
     301       10724 :     av = avma;
     302             :     /* G in nf, = L^e mod f */
     303       10724 :     G = famat_to_nf_modideal_coprime(nf, newL, e, f, EX);
     304       10724 :     if (mulI)
     305             :     {
     306        3752 :       G = nfmuli(nf, G, mulI);
     307        3752 :       G = typ(G) == t_COL? ZC_hnfrem(G, ZM_Z_mul(f, dmulI))
     308        3752 :                          : modii(G, mulii(fZ,dmulI));
     309        3752 :       G = RgC_Rg_div(G, dmulI);
     310             :     }
     311       10724 :     G = set_sign_mod_divisor(nf,A,G,sarch);
     312       10724 :     I = idealmul(nf,I,G);
     313             :     /* more or less useless, but cheap at this point */
     314       10724 :     I = idealmoddivisor_aux(nf,I,f,sarch);
     315       10724 :     gel(basecl,i) = gerepilecopy(av, I);
     316             :   }
     317       25830 :   return mkvec3(h, cyc, basecl);
     318             : }
     319             : 
     320             : /********************************************************************/
     321             : /**                                                                **/
     322             : /**                   INIT RAY CLASS GROUP                         **/
     323             : /**                                                                **/
     324             : /********************************************************************/
     325             : GEN
     326       83636 : bnr_subgroup_check(GEN bnr, GEN H, GEN *pdeg)
     327             : {
     328       83636 :   GEN no = bnr_get_no(bnr);
     329       83636 :   if (H && isintzero(H)) H = NULL;
     330       83636 :   if (H)
     331             :   {
     332       79961 :     GEN h, cyc = bnr_get_cyc(bnr);
     333       79961 :     switch(typ(H))
     334             :     {
     335          98 :       case t_INT:
     336          98 :         H = scalarmat_shallow(H, lg(cyc)-1);
     337             :         /* fall through */
     338        2807 :       case t_MAT:
     339        2807 :         RgM_check_ZM(H, "bnr_subgroup_check");
     340        2807 :         H = ZM_hnfmodid(H, cyc);
     341        2807 :         break;
     342       77154 :       case t_VEC:
     343       77154 :         if (char_check(cyc, H)) { H = charker(cyc, H); break; }
     344           0 :       default: pari_err_TYPE("bnr_subgroup_check", H);
     345             :     }
     346       79961 :     h = ZM_det_triangular(H);
     347       79961 :     if (equalii(h, no)) H = NULL; else no = h;
     348             :   }
     349       83636 :   if (pdeg) *pdeg = no;
     350       83636 :   return H;
     351             : }
     352             : 
     353             : void
     354        1617 : bnr_subgroup_sanitize(GEN *pbnr, GEN *pH)
     355             : {
     356        1617 :   GEN D, cnd, mod, cyc, bnr = *pbnr, H = *pH;
     357             : 
     358        1617 :   if (nftyp(bnr)==typ_BNF) bnr = Buchray(bnr, gen_1, nf_INIT);
     359        1561 :   else checkbnr(bnr);
     360        1603 :   cyc = bnr_get_cyc(bnr);
     361        1603 :   if (!H) mod = cyc_get_expo(cyc);
     362        1239 :   else switch(typ(H))
     363             :   {
     364         315 :     case t_INT: mod = H; break;
     365           7 :     case t_VEC:
     366           7 :       if (!char_check(cyc, H))
     367           0 :         pari_err_TYPE("bnr_subgroup_sanitize [character]", H);
     368           7 :       H = charker(cyc, H); /* character -> subgroup */
     369         917 :     case t_MAT:
     370         917 :       H = hnfmodid(H, cyc); /* make sure H is a left divisor of Mat(cyc) */
     371         903 :       D = ZM_snf(H); /* structure of Cl_f / H */
     372         903 :       mod = lg(D) == 1? gen_1: gel(D,1);
     373         903 :       break;
     374           7 :     default: pari_err_TYPE("bnr_subroup_sanitize [subgroup]", H);
     375           0 :       mod = NULL;
     376             :   }
     377        1582 :   cnd = bnrconductormod(bnr, H, mod);
     378        1582 :   *pbnr = gel(cnd,2); *pH = gel(cnd,3);
     379        1582 : }
     380             : void
     381        1141 : bnr_char_sanitize(GEN *pbnr, GEN *pchi)
     382             : {
     383        1141 :   GEN cnd, cyc, bnr = *pbnr, chi = *pchi;
     384             : 
     385        1141 :   if (nftyp(bnr)==typ_BNF) bnr = Buchray(bnr, gen_1, nf_INIT);
     386        1141 :   else checkbnr(bnr);
     387        1141 :   cyc = bnr_get_cyc(bnr);
     388        1141 :   if (typ(chi) != t_VEC || !char_check(cyc, chi))
     389           0 :     pari_err_TYPE("bnr_char_sanitize [character]", chi);
     390        1141 :   cnd = bnrconductormod(bnr, chi, charorder(cyc, chi));
     391        1141 :   *pbnr = gel(cnd,2); *pchi = gel(cnd,3);
     392        1141 : }
     393             : 
     394             : 
     395             : /* c a rational content (NULL or t_INT or t_FRAC), return u*c as a ZM/d */
     396             : static GEN
     397       37198 : ZM_content_mul(GEN u, GEN c, GEN *pd)
     398             : {
     399       37198 :   *pd = gen_1;
     400       37198 :   if (c)
     401             :   {
     402       26229 :     if (typ(c) == t_FRAC) { *pd = gel(c,2); c = gel(c,1); }
     403       26229 :     if (!is_pm1(c)) u = ZM_Z_mul(u, c);
     404             :   }
     405       37198 :   return u;
     406             : }
     407             : 
     408             : /* bnr natural generators: bnf gens made coprime to modulus + bid gens */
     409             : static GEN
     410       45633 : get_Gen(GEN bnf, GEN bid, GEN El)
     411             : {
     412       45633 :   GEN Gen = shallowconcat(bnf_get_gen(bnf), bid_get_gen(bid));
     413       45633 :   GEN nf = bnf_get_nf(bnf);
     414       45633 :   long i, l = lg(El);
     415       62888 :   for (i = 1; i < l; i++)
     416             :   {
     417       17255 :     GEN e = gel(El,i);
     418       17255 :     if (!isint1(e)) gel(Gen,i) = idealmul(nf, gel(El,i), gel(Gen,i));
     419             :   }
     420       45633 :   return Gen;
     421             : }
     422             : 
     423             : static GEN
     424       44226 : Buchraymod_i(GEN bnf, GEN module, long flag, GEN MOD)
     425             : {
     426             :   GEN nf, cyc0, cyc, gen, Cyc, clg, h, logU, U, Ui, vu;
     427             :   GEN bid, cycbid, H, El;
     428             :   long RU, Ri, j, ngen;
     429       44226 :   const long add_gen = flag & nf_GEN;
     430       44226 :   const long do_init = flag & nf_INIT;
     431             : 
     432       44226 :   if (MOD && typ(MOD) != t_INT)
     433           0 :     pari_err_TYPE("bnrinit [incorrect cycmod]", MOD);
     434       44226 :   bnf = checkbnf(bnf);
     435       44219 :   nf = bnf_get_nf(bnf);
     436       44219 :   RU = lg(nf_get_roots(nf))-1; /* #K.futu */
     437       44219 :   El = NULL; /* gcc -Wall */
     438       44219 :   cyc = cyc0 = bnf_get_cyc(bnf);
     439       44219 :   gen = bnf_get_gen(bnf); ngen = lg(cyc)-1;
     440             : 
     441       44219 :   bid = checkbid_i(module);
     442       44219 :   if (!bid) bid = Idealstarmod(nf,module,nf_GEN|nf_INIT, MOD);
     443       44219 :   cycbid = bid_get_cyc(bid);
     444       44219 :   if (MOD)
     445             :   {
     446        2506 :     cyc = ZV_snf_gcd(cyc, MOD);
     447        2506 :     cycbid = ZV_snf_gcd(cycbid, MOD);
     448             :   }
     449       44219 :   Ri = lg(cycbid)-1;
     450       44219 :   if (Ri || add_gen || do_init)
     451             :   {
     452       44219 :     GEN fx = bid_get_fact(bid);
     453       44219 :     El = cgetg(ngen+1,t_VEC);
     454       65632 :     for (j=1; j<=ngen; j++)
     455             :     {
     456       21413 :       GEN c = idealcoprimefact(nf, gel(gen,j), fx);
     457       21413 :       gel(El,j) = nf_to_scalar_or_basis(nf,c);
     458             :     }
     459             :   }
     460       44219 :   if (!Ri)
     461             :   {
     462        7021 :     GEN hK = bnf_get_no(bnf);
     463        7021 :     clg = add_gen? mkvec3(hK, cyc, get_Gen(bnf, bid, El)): mkvec2(hK, cyc);
     464        7021 :     if (!do_init) return clg;
     465        7021 :     U = matid(ngen);
     466        7021 :     U = mkvec3(U, cgetg(1,t_MAT), U);
     467        7021 :     vu = mkvec3(cgetg(1,t_MAT), matid(RU), gen_1);
     468        7021 :     return mkvecn(6, bnf, bid, El, U, clg, vu);
     469             :   }
     470             : 
     471       37198 :   logU = ideallog_units0(bnf, bid, MOD);
     472       37198 :   if (do_init)
     473             :   { /* (log(Units)|D) * u = (0 | H) */
     474       37198 :     GEN c1,c2, u,u1,u2, Hi, D = shallowconcat(logU, diagonal_shallow(cycbid));
     475       37198 :     H = ZM_hnfall_i(D, &u, 1);
     476       37198 :     u1 = matslice(u, 1,RU, 1,RU);
     477       37198 :     u2 = matslice(u, 1,RU, RU+1,lg(u)-1);
     478             :     /* log(Units) (u1|u2) = (0|H) (mod D), H HNF */
     479             : 
     480       37198 :     u1 = ZM_lll(u1, 0.99, LLL_INPLACE);
     481       37198 :     Hi = Q_primitive_part(RgM_inv_upper(H), &c1);
     482       37198 :     u2 = ZM_mul(ZM_reducemodmatrix(u2,u1), Hi);
     483       37198 :     u2 = Q_primitive_part(u2, &c2);
     484       37198 :     u2 = ZM_content_mul(u2, mul_content(c1,c2), &c2);
     485       37198 :     vu = mkvec3(u2,u1,c2); /* u2/c2 = H^(-1) (mod Im u1) */
     486             :   }
     487             :   else
     488             :   {
     489           0 :     H = ZM_hnfmodid(logU, cycbid);
     490           0 :     vu = NULL; /* -Wall */
     491             :   }
     492       37198 :   if (!ngen)
     493       24122 :     h = H;
     494             :   else
     495             :   {
     496       13076 :     GEN L = cgetg(ngen+1, t_MAT), cycgen = bnf_build_cycgen(bnf);
     497       26796 :     for (j=1; j<=ngen; j++)
     498             :     {
     499       13720 :       GEN c = gel(cycgen,j), e = gel(El,j);
     500       13720 :       if (!equali1(e)) c = famat_mulpow_shallow(c, e, gel(cyc0,j));
     501       13720 :       gel(L,j) = ideallogmod(nf, c, bid, MOD); /* = log(Gen[j]^cyc[j]) */
     502             :     }
     503             :     /* [cyc0, 0; -L, H] = relation matrix for generators Gen of Cl_f */
     504       13076 :     h = shallowconcat(vconcat(diagonal_shallow(cyc0), ZM_neg(L)),
     505             :                       vconcat(zeromat(ngen, Ri), H));
     506       13076 :     h = MOD? ZM_hnfmodid(h, MOD): ZM_hnf(h);
     507             :   }
     508       37198 :   Cyc = ZM_snf_group(h, &U, &Ui);
     509             :   /* Gen = clg.gen*U, clg.gen = Gen*Ui */
     510       32592 :   clg = add_gen? bnr_grp(nf, Ui, get_Gen(bnf, bid, El), Cyc, bid)
     511       37198 :                : mkvec2(ZV_prod(Cyc), Cyc);
     512       37198 :   if (!do_init) return clg;
     513       37198 :   U = mkvec3(vecslice(U, 1,ngen), vecslice(U,ngen+1,lg(U)-1), Ui);
     514       37198 :   return mkvecn(6, bnf, bid, El, U, clg, vu);
     515             : }
     516             : GEN
     517       40509 : Buchray(GEN bnf, GEN f, long flag)
     518       40509 : { return Buchraymod(bnf, f, flag, NULL); }
     519             : GEN
     520       42210 : Buchraymod(GEN bnf, GEN f, long flag, GEN MOD)
     521             : {
     522       42210 :   pari_sp av = avma;
     523       42210 :   return gerepilecopy(av, Buchraymod_i(bnf, f, flag, MOD));
     524             : }
     525             : GEN
     526        1169 : bnrinitmod(GEN bnf, GEN f, long flag, GEN MOD)
     527             : {
     528        1169 :   switch(flag)
     529             :   {
     530        1127 :     case 0: flag = nf_INIT; break;
     531          42 :     case 1: flag = nf_INIT | nf_GEN; break;
     532           0 :     default: pari_err_FLAG("bnrinit");
     533             :   }
     534        1169 :   return Buchraymod(bnf, f, flag, MOD);
     535             : }
     536             : GEN
     537           0 : bnrinit0(GEN bnf, GEN ideal, long flag)
     538           0 : { return bnrinitmod(bnf, ideal, flag, NULL); }
     539             : 
     540             : GEN
     541         112 : bnrclassno(GEN bnf,GEN ideal)
     542             : {
     543             :   GEN h, D, bid, cycbid;
     544         112 :   pari_sp av = avma;
     545             : 
     546         112 :   bnf = checkbnf(bnf);
     547         112 :   h = bnf_get_no(bnf);
     548         112 :   bid = checkbid_i(ideal);
     549         112 :   if (!bid) bid = Idealstar(bnf, ideal, nf_INIT);
     550         105 :   cycbid = bid_get_cyc(bid);
     551         105 :   if (lg(cycbid) == 1) { set_avma(av); return icopy(h); }
     552          84 :   D = ideallog_units(bnf, bid); /* (Z_K/f)^* / units ~ Z^n / D */
     553          84 :   D = ZM_hnfmodid(D,cycbid);
     554          84 :   return gerepileuptoint(av, mulii(h, ZM_det_triangular(D)));
     555             : }
     556             : GEN
     557         105 : bnrclassno0(GEN A, GEN B, GEN C)
     558             : {
     559         105 :   pari_sp av = avma;
     560         105 :   GEN h, H = NULL;
     561             :   /* adapted from ABC_to_bnr, avoid costly bnrinit if possible */
     562         105 :   if (typ(A) == t_VEC)
     563         105 :     switch(lg(A))
     564             :     {
     565          14 :       case 7: /* bnr */
     566          14 :         checkbnr(A); H = B;
     567          14 :         break;
     568          91 :       case 11: /* bnf */
     569          91 :         if (!B) pari_err_TYPE("bnrclassno [bnf+missing conductor]",A);
     570          91 :         if (!C) return bnrclassno(A, B);
     571           7 :         A = Buchray(A, B, nf_INIT); H = C;
     572           7 :         break;
     573           0 :       default: checkbnf(A);/*error*/
     574             :     }
     575           0 :   else checkbnf(A);/*error*/
     576             : 
     577          21 :   H = bnr_subgroup_check(A, H, &h);
     578          21 :   if (!H) { set_avma(av); return icopy(h); }
     579          14 :   return gerepileuptoint(av, h);
     580             : }
     581             : 
     582             : /* ZMV_ZCV_mul for two matrices U = [Ux,Uy], it may have more components
     583             :  * (ignored) and vectors x,y */
     584             : static GEN
     585      271983 : ZM2_ZC2_mul(GEN U, GEN x, GEN y)
     586             : {
     587      271983 :   GEN Ux = gel(U,1), Uy = gel(U,2);
     588      271983 :   if (lg(Ux) == 1) return ZM_ZC_mul(Uy,y);
     589      154907 :   if (lg(Uy) == 1) return ZM_ZC_mul(Ux,x);
     590      154907 :   return ZC_add(ZM_ZC_mul(Ux,x), ZM_ZC_mul(Uy,y));
     591             : }
     592             : 
     593             : GEN
     594      385559 : bnrisprincipalmod(GEN bnr, GEN x, GEN MOD, long flag)
     595             : {
     596      385559 :   pari_sp av = avma;
     597             :   GEN E, G, clgp, bnf, nf, bid, ex, cycray, alpha, El;
     598             :   int trivialbid;
     599             : 
     600      385559 :   checkbnr(bnr);
     601      385559 :   El = bnr_get_El(bnr);
     602      385559 :   cycray = bnr_get_cyc(bnr);
     603      385559 :   if (MOD && flag) pari_err_FLAG("bnrisprincipalmod [MOD!=NULL and flag!=0]");
     604      385559 :   if (lg(cycray) == 1 && !(flag & nf_GEN)) return cgetg(1,t_COL);
     605      385559 :   if (MOD) cycray = ZV_snf_gcd(cycray, MOD);
     606             : 
     607      385559 :   bnf = bnr_get_bnf(bnr); nf = bnf_get_nf(bnf);
     608      385559 :   bid = bnr_get_bid(bnr);
     609      385559 :   trivialbid = lg(bid_get_cyc(bid)) == 1;
     610      385559 :   if (trivialbid) ex = isprincipal(bnf, x);
     611             :   else
     612             :   {
     613      271983 :     GEN v = bnfisprincipal0(bnf, x, nf_FORCE|nf_GENMAT);
     614      271983 :     GEN e = gel(v,1), b = gel(v,2);
     615      271983 :     long i, j = lg(e);
     616      432007 :     for (i = 1; i < j; i++) /* modify b as if bnf.gen were El*bnf.gen */
     617      160024 :       if (typ(gel(El,i)) != t_INT && signe(gel(e,i))) /* <==> != 1 */
     618       29078 :         b = famat_mulpow_shallow(b, gel(El,i), negi(gel(e,i)));
     619      271983 :     if (!MOD && !(flag & nf_GEN)) MOD = gel(cycray,1);
     620      271983 :     ex = ZM2_ZC2_mul(bnr_get_U(bnr), e, ideallogmod(nf, b, bid, MOD));
     621             :   }
     622      385559 :   ex = vecmodii(ex, cycray);
     623      385559 :   if (!(flag & (nf_GEN|nf_GENMAT))) return gerepileupto(av, ex);
     624             : 
     625             :   /* compute generator */
     626        6552 :   E = ZC_neg(ex);
     627        6552 :   clgp = bnr_get_clgp(bnr);
     628        6552 :   if (lg(clgp) == 4)
     629          21 :     G = abgrp_get_gen(clgp);
     630             :   else
     631             :   {
     632        6531 :     G = get_Gen(bnf, bid, El);
     633        6531 :     E = ZM_ZC_mul(bnr_get_Ui(bnr), E);
     634             :   }
     635        6552 :   alpha = isprincipalfact(bnf, x, G, E, nf_GENMAT|nf_GEN_IF_PRINCIPAL|nf_FORCE);
     636        6552 :   if (alpha == gen_0) pari_err_BUG("isprincipalray");
     637        6552 :   if (!trivialbid)
     638             :   {
     639        6552 :     GEN v = gel(bnr,6), u2 = gel(v,1), u1 = gel(v,2), du2 = gel(v,3);
     640        6552 :     GEN y = ZM_ZC_mul(u2, ideallog(nf, alpha, bid));
     641        6552 :     if (!is_pm1(du2)) y = ZC_Z_divexact(y,du2);
     642        6552 :     y = ZC_reducemodmatrix(y, u1);
     643        6552 :     if (!ZV_equal0(y))
     644             :     {
     645        4501 :       GEN U = shallowcopy(bnf_build_units(bnf));
     646        4501 :       settyp(U, t_COL);
     647        4501 :       alpha = famat_div_shallow(alpha, mkmat2(U,y));
     648             :     }
     649             :   }
     650        6552 :   alpha = famat_reduce(alpha);
     651        6552 :   if (!(flag & nf_GENMAT)) alpha = nffactorback(nf, alpha, NULL);
     652        6552 :   return gerepilecopy(av, mkvec2(ex,alpha));
     653             : }
     654             : 
     655             : GEN
     656      364523 : bnrisprincipal(GEN bnr, GEN x, long flag)
     657      364523 : { return bnrisprincipalmod(bnr, x, NULL, flag); }
     658             : 
     659             : GEN
     660      357943 : isprincipalray(GEN bnr, GEN x) { return bnrisprincipal(bnr,x,0); }
     661             : GEN
     662           0 : isprincipalraygen(GEN bnr, GEN x) { return bnrisprincipal(bnr,x,nf_GEN); }
     663             : 
     664             : /* N! / N^N * (4/pi)^r2 * sqrt(|D|) */
     665             : GEN
     666           0 : minkowski_bound(GEN D, long N, long r2, long prec)
     667             : {
     668           0 :   pari_sp av = avma;
     669           0 :   GEN c = divri(mpfactr(N,prec), powuu(N,N));
     670           0 :   if (r2) c = mulrr(c, powru(divur(4,mppi(prec)), r2));
     671           0 :   c = mulrr(c, gsqrt(absi_shallow(D),prec));
     672           0 :   return gerepileuptoleaf(av, c);
     673             : }
     674             : 
     675             : /* N = [K:Q] > 1, D = disc(K) */
     676             : static GEN
     677          56 : zimmertbound(GEN D, long N, long R2)
     678             : {
     679          56 :   pari_sp av = avma;
     680             :   GEN w;
     681             : 
     682          56 :   if (N > 20) w = minkowski_bound(D, N, R2, DEFAULTPREC);
     683             :   else
     684             :   {
     685          56 :     const double c[19][11] = {
     686             : {/*2*/  0.6931,     0.45158},
     687             : {/*3*/  1.71733859, 1.37420604},
     688             : {/*4*/  2.91799837, 2.50091538, 2.11943331},
     689             : {/*5*/  4.22701425, 3.75471588, 3.31196660},
     690             : {/*6*/  5.61209925, 5.09730381, 4.60693851, 4.14303665},
     691             : {/*7*/  7.05406203, 6.50550021, 5.97735406, 5.47145968},
     692             : {/*8*/  8.54052636, 7.96438858, 7.40555445, 6.86558259, 6.34608077},
     693             : {/*9*/ 10.0630022,  9.46382812, 8.87952524, 8.31139202, 7.76081149},
     694             : {/*10*/11.6153797, 10.9966020, 10.3907654,  9.79895170, 9.22232770, 8.66213267},
     695             : {/*11*/13.1930961, 12.5573772, 11.9330458, 11.3210061, 10.7222412, 10.1378082},
     696             : {/*12*/14.7926394, 14.1420915, 13.5016616, 12.8721114, 12.2542699, 11.6490374,
     697             :        11.0573775},
     698             : {/*13*/16.4112395, 15.7475710, 15.0929680, 14.4480777, 13.8136054, 13.1903162,
     699             :        12.5790381},
     700             : {/*14*/18.0466672, 17.3712806, 16.7040780, 16.0456127, 15.3964878, 14.7573587,
     701             :        14.1289364, 13.5119848},
     702             : {/*15*/19.6970961, 19.0111606, 18.3326615, 17.6620757, 16.9999233, 16.3467686,
     703             :        15.7032228, 15.0699480},
     704             : {/*16*/21.3610081, 20.6655103, 19.9768082, 19.2953176, 18.6214885, 17.9558093,
     705             :        17.2988108, 16.6510652, 16.0131906},
     706             : 
     707             : {/*17*/23.0371259, 22.3329066, 21.6349299, 20.9435607, 20.2591899, 19.5822454,
     708             :        18.9131878, 18.2525157, 17.6007672},
     709             : 
     710             : {/*18*/24.7243611, 24.0121449, 23.3056902, 22.6053167, 21.9113705, 21.2242247,
     711             :        20.5442836, 19.8719830, 19.2077941, 18.5522234},
     712             : 
     713             : {/*19*/26.4217792, 25.7021950, 24.9879497, 24.2793271, 23.5766321, 22.8801952,
     714             :        22.1903709, 21.5075437, 20.8321263, 20.1645647},
     715             : {/*20*/28.1285704, 27.4021674, 26.6807314, 25.9645140, 25.2537867, 24.5488420,
     716             :        23.8499943, 23.1575823, 22.4719720, 21.7935548, 21.1227537}
     717             :     };
     718          56 :     w = mulrr(dbltor(exp(-c[N-2][R2])), gsqrt(absi_shallow(D),DEFAULTPREC));
     719             :   }
     720          56 :   return gerepileuptoint(av, ceil_safe(w));
     721             : }
     722             : 
     723             : /* return \gamma_n^n if known, an upper bound otherwise */
     724             : static GEN
     725          56 : hermiteconstant(long n)
     726             : {
     727             :   GEN h,h1;
     728             :   pari_sp av;
     729             : 
     730          56 :   switch(n)
     731             :   {
     732          28 :     case 1: return gen_1;
     733          14 :     case 2: return mkfrac(utoipos(4), utoipos(3));
     734           7 :     case 3: return gen_2;
     735           7 :     case 4: return utoipos(4);
     736           0 :     case 5: return utoipos(8);
     737           0 :     case 6: return mkfrac(utoipos(64), utoipos(3));
     738           0 :     case 7: return utoipos(64);
     739           0 :     case 8: return utoipos(256);
     740             :   }
     741           0 :   av = avma;
     742           0 :   h  = powru(divur(2,mppi(DEFAULTPREC)), n);
     743           0 :   h1 = sqrr(ggamma(sstoQ(n+4,2),DEFAULTPREC));
     744           0 :   return gerepileuptoleaf(av, mulrr(h,h1));
     745             : }
     746             : 
     747             : /* 1 if L (= nf != Q) primitive for sure, 0 if MAYBE imprimitive (may have a
     748             :  * subfield K) */
     749             : static long
     750          28 : isprimitive(GEN nf)
     751             : {
     752          28 :   long p, i, l, ep, N = nf_get_degree(nf);
     753             :   GEN D, fa;
     754             : 
     755          28 :   p = ucoeff(factoru(N), 1,1); /* smallest prime | N */
     756          28 :   if (p == N) return 1; /* prime degree */
     757             : 
     758             :   /* N = [L:Q] = product of primes >= p, same is true for [L:K]
     759             :    * d_L = t d_K^[L:K] --> check that some q^p divides d_L */
     760           0 :   D = nf_get_disc(nf);
     761           0 :   fa = gel(absZ_factor_limit(D,0),2); /* list of v_q(d_L). Don't check large primes */
     762           0 :   if (mod2(D)) i = 1;
     763             :   else
     764             :   { /* q = 2 */
     765           0 :     ep = itos(gel(fa,1));
     766           0 :     if ((ep>>1) >= p) return 0; /* 2 | d_K ==> 4 | d_K */
     767           0 :     i = 2;
     768             :   }
     769           0 :   l = lg(fa);
     770           0 :   for ( ; i < l; i++)
     771             :   {
     772           0 :     ep = itos(gel(fa,i));
     773           0 :     if (ep >= p) return 0;
     774             :   }
     775           0 :   return 1;
     776             : }
     777             : 
     778             : static GEN
     779           0 : dft_bound(void)
     780             : {
     781           0 :   if (DEBUGLEVEL>1) err_printf("Default bound for regulator: 0.2\n");
     782           0 :   return dbltor(0.2);
     783             : }
     784             : 
     785             : static GEN
     786          28 : regulatorbound(GEN bnf)
     787             : {
     788             :   long N, R1, R2, R;
     789             :   GEN nf, dK, p1, c1;
     790             : 
     791          28 :   nf = bnf_get_nf(bnf); N = nf_get_degree(nf);
     792          28 :   if (!isprimitive(nf)) return dft_bound();
     793             : 
     794          28 :   dK = absi_shallow(nf_get_disc(nf));
     795          28 :   nf_get_sign(nf, &R1, &R2); R = R1+R2-1;
     796          28 :   c1 = (!R2 && N<12)? int2n(N & (~1UL)): powuu(N,N);
     797          28 :   if (cmpii(dK,c1) <= 0) return dft_bound();
     798             : 
     799          28 :   p1 = sqrr(glog(gdiv(dK,c1),DEFAULTPREC));
     800          28 :   p1 = divru(gmul2n(powru(divru(mulru(p1,3),N*(N*N-1)-6*R2),R),R2), N);
     801          28 :   p1 = sqrtr(gdiv(p1, hermiteconstant(R)));
     802          28 :   if (DEBUGLEVEL>1) err_printf("Mahler bound for regulator: %Ps\n",p1);
     803          28 :   return gmax_shallow(p1, dbltor(0.2));
     804             : }
     805             : 
     806             : static int
     807       70553 : is_unit(GEN M, long r1, GEN x)
     808             : {
     809       70553 :   pari_sp av = avma;
     810       70553 :   GEN Nx = ground( embed_norm(RgM_zc_mul(M,x), r1) );
     811       70553 :   return gc_bool(av, is_pm1(Nx));
     812             : }
     813             : 
     814             : /* FIXME: should use smallvectors */
     815             : static double
     816          42 : minimforunits(GEN nf, long BORNE, ulong w)
     817             : {
     818          42 :   const long prec = MEDDEFAULTPREC;
     819          42 :   long n, r1, i, j, k, *x, cnt = 0;
     820          42 :   pari_sp av = avma;
     821             :   GEN r, M;
     822             :   double p, norme, normin;
     823             :   double **q,*v,*y,*z;
     824          42 :   double eps=0.000001, BOUND = BORNE * 1.00001;
     825             : 
     826          42 :   if (DEBUGLEVEL>=2)
     827             :   {
     828           0 :     err_printf("Searching minimum of T2-form on units:\n");
     829           0 :     if (DEBUGLEVEL>2) err_printf("   BOUND = %ld\n",BORNE);
     830             :   }
     831          42 :   n = nf_get_degree(nf); r1 = nf_get_r1(nf);
     832          42 :   minim_alloc(n+1, &q, &x, &y, &z, &v);
     833          42 :   M = gprec_w(nf_get_M(nf), prec);
     834          42 :   r = gaussred_from_QR(nf_get_G(nf), prec);
     835         231 :   for (j=1; j<=n; j++)
     836             :   {
     837         189 :     v[j] = gtodouble(gcoeff(r,j,j));
     838         651 :     for (i=1; i<j; i++) q[i][j] = gtodouble(gcoeff(r,i,j));
     839             :   }
     840          42 :   normin = (double)BORNE*(1-eps);
     841          42 :   k=n; y[n]=z[n]=0;
     842          42 :   x[n] = (long)(sqrt(BOUND/v[n]));
     843             : 
     844       70553 :   for(;;x[1]--)
     845             :   {
     846             :     do
     847             :     {
     848       71953 :       if (k>1)
     849             :       {
     850        1400 :         long l = k-1;
     851        1400 :         z[l] = 0;
     852        5334 :         for (j=k; j<=n; j++) z[l] += q[l][j]*x[j];
     853        1400 :         p = (double)x[k] + z[k];
     854        1400 :         y[l] = y[k] + p*p*v[k];
     855        1400 :         x[l] = (long)floor(sqrt((BOUND-y[l])/v[l])-z[l]);
     856        1400 :         k = l;
     857             :       }
     858             :       for(;;)
     859             :       {
     860       74459 :         p = (double)x[k] + z[k];
     861       73206 :         if (y[k] + p*p*v[k] <= BOUND) break;
     862        1253 :         k++; x[k]--;
     863             :       }
     864             :     }
     865       71953 :     while (k>1);
     866       70595 :     if (!x[1] && y[1]<=eps) break;
     867             : 
     868       70567 :     if (DEBUGLEVEL>8) err_printf(".");
     869       70567 :     if (++cnt == 5000) return -1.; /* too expensive */
     870             : 
     871       70553 :     p = (double)x[1] + z[1]; norme = y[1] + p*p*v[1];
     872       70553 :     if (is_unit(M, r1, x) && norme < normin)
     873             :     {
     874             :       /* exclude roots of unity */
     875          56 :       if (norme < 2*n)
     876             :       {
     877          42 :         GEN t = nfpow_u(nf, zc_to_ZC(x), w);
     878          42 :         if (typ(t) != t_COL || ZV_isscalar(t)) continue;
     879             :       }
     880          21 :       normin = norme*(1-eps);
     881          21 :       if (DEBUGLEVEL>=2) err_printf("*");
     882             :     }
     883             :   }
     884          28 :   if (DEBUGLEVEL>=2) err_printf("\n");
     885          28 :   set_avma(av);
     886          28 :   return normin;
     887             : }
     888             : 
     889             : #undef NBMAX
     890             : static int
     891        1804 : is_zero(GEN x, long bitprec) { return (gexpo(x) < -bitprec); }
     892             : 
     893             : static int
     894        1228 : is_complex(GEN x, long bitprec) { return !is_zero(imag_i(x), bitprec); }
     895             : 
     896             : /* assume M_star t_REAL
     897             :  * FIXME: what does this do ? To be rewritten */
     898             : static GEN
     899          28 : compute_M0(GEN M_star,long N)
     900             : {
     901             :   long m1,m2,n1,n2,n3,lr,lr1,lr2,i,j,l,vx,vy,vz,vM;
     902             :   GEN pol,p1,p2,p3,p4,p5,p6,p7,p8,p9,u,v,w,r,r1,r2,M0,M0_pro,S,P,M;
     903             :   GEN f1,f2,f3,g1,g2,g3,pg1,pg2,pg3,pf1,pf2,pf3,X,Y,Z;
     904          28 :   long bitprec = 24;
     905             : 
     906          28 :   if (N == 2) return gmul2n(sqrr(gacosh(gmul2n(M_star,-1),0)), -1);
     907          21 :   vx = fetch_var(); X = pol_x(vx);
     908          21 :   vy = fetch_var(); Y = pol_x(vy);
     909          21 :   vz = fetch_var(); Z = pol_x(vz);
     910          21 :   vM = fetch_var(); M = pol_x(vM);
     911             : 
     912          21 :   M0 = NULL; m1 = N/3;
     913          56 :   for (n1=1; n1<=m1; n1++) /* 1 <= n1 <= n2 <= n3 < N */
     914             :   {
     915          35 :     m2 = (N-n1)>>1;
     916         112 :     for (n2=n1; n2<=m2; n2++)
     917             :     {
     918          77 :       pari_sp av = avma; n3=N-n1-n2;
     919          77 :       if (n1==n2 && n1==n3) /* n1 = n2 = n3 = m1 = N/3 */
     920             :       {
     921           7 :         p1 = divru(M_star, m1);
     922           7 :         p4 = sqrtr_abs( mulrr(addsr(1,p1),subrs(p1,3)) );
     923           7 :         p5 = subrs(p1,1);
     924           7 :         u = gen_1;
     925           7 :         v = gmul2n(addrr(p5,p4),-1);
     926           7 :         w = gmul2n(subrr(p5,p4),-1);
     927           7 :         M0_pro=gmul2n(mulur(m1,addrr(sqrr(logr_abs(v)),sqrr(logr_abs(w)))), -2);
     928           7 :         if (DEBUGLEVEL>2)
     929           0 :           err_printf("[ %ld, %ld, %ld ]: %.28Pg\n",n1,n2,n3,M0_pro);
     930           7 :         if (!M0 || gcmp(M0_pro,M0) < 0) M0 = M0_pro;
     931             :       }
     932          70 :       else if (n1==n2 || n2==n3)
     933          42 :       { /* n3 > N/3 >= n1 */
     934          42 :         long k = N - 2*n2;
     935          42 :         p2 = deg1pol_shallow(stoi(-n2), M_star, vx); /* M* - n2 X */
     936          42 :         p3 = gmul(powuu(k,k),
     937             :                   gpowgs(gsubgs(RgX_Rg_mul(p2, M_star), k*k), n2));
     938          42 :         pol = gsub(p3, RgX_mul(monomial(powuu(n2,n2), n2, vx),
     939             :                                gpowgs(p2, N-n2)));
     940          42 :         r = roots(pol, DEFAULTPREC); lr = lg(r);
     941         378 :         for (i=1; i<lr; i++)
     942             :         {
     943             :           GEN n2S;
     944         336 :           S = real_i(gel(r,i));
     945         336 :           if (is_complex(gel(r,i), bitprec) || signe(S) <= 0) continue;
     946             : 
     947         182 :           n2S = mulur(n2,S);
     948         182 :           p4 = subrr(M_star, n2S);
     949         182 :           P = divrr(mulrr(n2S,p4), subrs(mulrr(M_star,p4),k*k));
     950         182 :           p5 = subrr(sqrr(S), gmul2n(P,2));
     951         182 :           if (gsigne(p5) < 0) continue;
     952             : 
     953         140 :           p6 = sqrtr(p5);
     954         140 :           v = gmul2n(subrr(S,p6),-1);
     955         140 :           if (gsigne(v) <= 0) continue;
     956             : 
     957         126 :           u = gmul2n(addrr(S,p6),-1);
     958         126 :           w = gpow(P, sstoQ(-n2,k), 0);
     959         126 :           p6 = mulur(n2, addrr(sqrr(logr_abs(u)), sqrr(logr_abs(v))));
     960         126 :           M0_pro = gmul2n(addrr(p6, mulur(k, sqrr(logr_abs(w)))),-2);
     961         126 :           if (DEBUGLEVEL>2)
     962           0 :             err_printf("[ %ld, %ld, %ld ]: %.28Pg\n",n1,n2,n3,M0_pro);
     963         126 :           if (!M0 || gcmp(M0_pro,M0) < 0) M0 = M0_pro;
     964             :         }
     965             :       }
     966             :       else
     967             :       {
     968          28 :         f1 = gsub(gadd(gmulsg(n1,X),gadd(gmulsg(n2,Y),gmulsg(n3,Z))), M);
     969          28 :         f2 =         gmulsg(n1,gmul(Y,Z));
     970          28 :         f2 = gadd(f2,gmulsg(n2,gmul(X,Z)));
     971          28 :         f2 = gadd(f2,gmulsg(n3,gmul(X,Y)));
     972          28 :         f2 = gsub(f2,gmul(M,gmul(X,gmul(Y,Z))));
     973          28 :         f3 = gsub(gmul(gpowgs(X,n1),gmul(gpowgs(Y,n2),gpowgs(Z,n3))), gen_1);
     974             :         /* f1 = n1 X + n2 Y + n3 Z - M */
     975             :         /* f2 = n1 YZ + n2 XZ + n3 XY */
     976             :         /* f3 = X^n1 Y^n2 Z^n3 - 1*/
     977          28 :         g1=resultant(f1,f2); g1=primpart(g1);
     978          28 :         g2=resultant(f1,f3); g2=primpart(g2);
     979          28 :         g3=resultant(g1,g2); g3=primpart(g3);
     980          28 :         pf1=gsubst(f1,vM,M_star); pg1=gsubst(g1,vM,M_star);
     981          28 :         pf2=gsubst(f2,vM,M_star); pg2=gsubst(g2,vM,M_star);
     982          28 :         pf3=gsubst(f3,vM,M_star); pg3=gsubst(g3,vM,M_star);
     983             :         /* g3 = Res_Y,Z(f1,f2,f3) */
     984          28 :         r = roots(pg3,DEFAULTPREC); lr = lg(r);
     985         476 :         for (i=1; i<lr; i++)
     986             :         {
     987         448 :           w = real_i(gel(r,i));
     988         448 :           if (is_complex(gel(r,i), bitprec) || signe(w) <= 0) continue;
     989         140 :           p1=gsubst(pg1,vz,w);
     990         140 :           p2=gsubst(pg2,vz,w);
     991         140 :           p3=gsubst(pf1,vz,w);
     992         140 :           p4=gsubst(pf2,vz,w);
     993         140 :           p5=gsubst(pf3,vz,w);
     994         140 :           r1 = roots(p1, DEFAULTPREC); lr1 = lg(r1);
     995         420 :           for (j=1; j<lr1; j++)
     996             :           {
     997         280 :             v = real_i(gel(r1,j));
     998         280 :             if (is_complex(gel(r1,j), bitprec) || signe(v) <= 0
     999         280 :              || !is_zero(gsubst(p2,vy,v), bitprec)) continue;
    1000             : 
    1001         164 :             p7=gsubst(p3,vy,v);
    1002         164 :             p8=gsubst(p4,vy,v);
    1003         164 :             p9=gsubst(p5,vy,v);
    1004         164 :             r2 = roots(p7, DEFAULTPREC); lr2 = lg(r2);
    1005         328 :             for (l=1; l<lr2; l++)
    1006             :             {
    1007         164 :               u = real_i(gel(r2,l));
    1008         164 :               if (is_complex(gel(r2,l), bitprec) || signe(u) <= 0
    1009         164 :                || !is_zero(gsubst(p8,vx,u), bitprec)
    1010         164 :                || !is_zero(gsubst(p9,vx,u), bitprec)) continue;
    1011             : 
    1012         164 :               M0_pro =              mulur(n1, sqrr(logr_abs(u)));
    1013         164 :               M0_pro = gadd(M0_pro, mulur(n2, sqrr(logr_abs(v))));
    1014         164 :               M0_pro = gadd(M0_pro, mulur(n3, sqrr(logr_abs(w))));
    1015         164 :               M0_pro = gmul2n(M0_pro,-2);
    1016         164 :               if (DEBUGLEVEL>2)
    1017           0 :                 err_printf("[ %ld, %ld, %ld ]: %.28Pg\n",n1,n2,n3,M0_pro);
    1018         164 :               if (!M0 || gcmp(M0_pro,M0) < 0) M0 = M0_pro;
    1019             :             }
    1020             :           }
    1021             :         }
    1022             :       }
    1023          77 :       if (!M0) set_avma(av); else M0 = gerepilecopy(av, M0);
    1024             :     }
    1025             :   }
    1026         105 :   for (i=1;i<=4;i++) (void)delete_var();
    1027          21 :   return M0? M0: gen_0;
    1028             : }
    1029             : 
    1030             : static GEN
    1031          56 : lowerboundforregulator(GEN bnf, GEN units)
    1032             : {
    1033          56 :   long i, N, R2, RU = lg(units)-1;
    1034             :   GEN nf, M0, M, G, minunit;
    1035             :   double bound;
    1036             : 
    1037          56 :   if (!RU) return gen_1;
    1038          56 :   nf = bnf_get_nf(bnf);
    1039          56 :   N = nf_get_degree(nf);
    1040          56 :   R2 = nf_get_r2(nf);
    1041             : 
    1042          56 :   G = nf_get_G(nf);
    1043          56 :   minunit = gnorml2(RgM_RgC_mul(G, gel(units,1))); /* T2(units[1]) */
    1044         105 :   for (i=2; i<=RU; i++)
    1045             :   {
    1046          49 :     GEN t = gnorml2(RgM_RgC_mul(G, gel(units,i)));
    1047          49 :     if (gcmp(t,minunit) < 0) minunit = t;
    1048             :   }
    1049          56 :   if (gexpo(minunit) > 30) return NULL;
    1050             : 
    1051          42 :   bound = minimforunits(nf, itos(gceil(minunit)), bnf_get_tuN(bnf));
    1052          42 :   if (bound < 0) return NULL;
    1053          28 :   if (DEBUGLEVEL>1) err_printf("M* = %Ps\n", dbltor(bound));
    1054          28 :   M0 = compute_M0(dbltor(bound), N);
    1055          28 :   if (DEBUGLEVEL>1) err_printf("M0 = %.28Pg\n",M0);
    1056          28 :   M = gmul2n(divru(gdiv(powrs(M0,RU),hermiteconstant(RU)),N),R2);
    1057          28 :   if (cmprr(M, dbltor(0.04)) < 0) return NULL;
    1058          28 :   M = sqrtr(M);
    1059          28 :   if (DEBUGLEVEL>1)
    1060           0 :     err_printf("(lower bound for regulator) M = %.28Pg\n",M);
    1061          28 :   return M;
    1062             : }
    1063             : 
    1064             : /* upper bound for the index of bnf.fu in the full unit group */
    1065             : static GEN
    1066          56 : bound_unit_index(GEN bnf, GEN units)
    1067             : {
    1068          56 :   pari_sp av = avma;
    1069          56 :   GEN x = lowerboundforregulator(bnf, units);
    1070          56 :   if (!x) { set_avma(av); x = regulatorbound(bnf); }
    1071          56 :   return gerepileuptoint(av, ground(gdiv(bnf_get_reg(bnf), x)));
    1072             : }
    1073             : 
    1074             : /* Compute a square matrix of rank #beta attached to a family
    1075             :  * (P_i), 1<=i<=#beta, of primes s.t. N(P_i) = 1 mod p, and
    1076             :  * (P_i,beta[j]) = 1 for all i,j. nf = true nf */
    1077             : static void
    1078        1554 : primecertify(GEN nf, GEN beta, ulong p, GEN bad)
    1079             : {
    1080        1554 :   long lb = lg(beta), rmax = lb - 1;
    1081             :   GEN M, vQ, L;
    1082             :   ulong q;
    1083             :   forprime_t T;
    1084             : 
    1085        1554 :   if (p == 2)
    1086          42 :     L = cgetg(1,t_VECSMALL);
    1087             :   else
    1088        1512 :     L = mkvecsmall(p);
    1089        1554 :   (void)u_forprime_arith_init(&T, 1, ULONG_MAX, 1, p);
    1090        1554 :   M = cgetg(lb,t_MAT); setlg(M,1);
    1091        3220 :   while ((q = u_forprime_next(&T)))
    1092             :   {
    1093             :     GEN qq, gg, og;
    1094             :     long lQ, i, j;
    1095             :     ulong g, m;
    1096        3220 :     if (!umodiu(bad,q)) continue;
    1097             : 
    1098        2940 :     qq = utoipos(q);
    1099        2940 :     vQ = idealprimedec_limit_f(nf,qq,1);
    1100        2940 :     lQ = lg(vQ); if (lQ == 1) continue;
    1101             : 
    1102             :     /* cf rootsof1_Fl */
    1103        1946 :     g = pgener_Fl_local(q, L);
    1104        1946 :     m = (q-1) / p;
    1105        1946 :     gg = utoipos( Fl_powu(g, m, q) ); /* order p in (Z/q)^* */
    1106        1946 :     og = mkmat2(mkcol(utoi(p)), mkcol(gen_1)); /* order of g */
    1107             : 
    1108        1946 :     if (DEBUGLEVEL>3) err_printf("       generator of (Zk/Q)^*: %lu\n", g);
    1109        2611 :     for (i = 1; i < lQ; i++)
    1110             :     {
    1111        2219 :       GEN C = cgetg(lb, t_VECSMALL);
    1112        2219 :       GEN Q = gel(vQ,i); /* degree 1 */
    1113        2219 :       GEN modpr = zkmodprinit(nf, Q);
    1114             :       long r;
    1115             : 
    1116        6643 :       for (j = 1; j < lb; j++)
    1117             :       {
    1118        4424 :         GEN t = nf_to_Fp_coprime(nf, gel(beta,j), modpr);
    1119        4424 :         t = utoipos( Fl_powu(t[2], m, q) );
    1120        4424 :         C[j] = itou( Fp_log(t, gg, og, qq) ) % p;
    1121             :       }
    1122        2219 :       r = lg(M);
    1123        2219 :       gel(M,r) = C; setlg(M, r+1);
    1124        2219 :       if (Flm_rank(M, p) != r) { setlg(M,r); continue; }
    1125             : 
    1126        2016 :       if (DEBUGLEVEL>2)
    1127             :       {
    1128           0 :         if (DEBUGLEVEL>3)
    1129             :         {
    1130           0 :           err_printf("       prime ideal Q: %Ps\n",Q);
    1131           0 :           err_printf("       matrix log(b_j mod Q_i): %Ps\n", M);
    1132             :         }
    1133           0 :         err_printf("       new rank: %ld\n",r);
    1134             :       }
    1135        2016 :       if (r == rmax) return;
    1136             :     }
    1137             :   }
    1138           0 :   pari_err_BUG("primecertify");
    1139             : }
    1140             : 
    1141             : struct check_pr {
    1142             :   long w; /* #mu(K) */
    1143             :   GEN mu; /* generator of mu(K) */
    1144             :   GEN fu;
    1145             :   GEN cyc;
    1146             :   GEN cycgen;
    1147             :   GEN bad; /* p | bad <--> p | some element occurring in cycgen */
    1148             : };
    1149             : 
    1150             : static void
    1151        1554 : check_prime(ulong p, GEN nf, struct check_pr *S)
    1152             : {
    1153        1554 :   pari_sp av = avma;
    1154        1554 :   long i,b, lc = lg(S->cyc), lf = lg(S->fu);
    1155        1554 :   GEN beta = cgetg(lf+lc, t_VEC);
    1156             : 
    1157        1554 :   if (DEBUGLEVEL>1) err_printf("  *** testing p = %lu\n",p);
    1158        1617 :   for (b=1; b<lc; b++)
    1159             :   {
    1160        1323 :     if (umodiu(gel(S->cyc,b), p)) break; /* p \nmid cyc[b] */
    1161          63 :     if (b==1 && DEBUGLEVEL>2) err_printf("     p divides h(K)\n");
    1162          63 :     gel(beta,b) = gel(S->cycgen,b);
    1163             :   }
    1164        1554 :   if (S->w % p == 0)
    1165             :   {
    1166          42 :     if (DEBUGLEVEL>2) err_printf("     p divides w(K)\n");
    1167          42 :     gel(beta,b++) = S->mu;
    1168             :   }
    1169        3465 :   for (i=1; i<lf; i++) gel(beta,b++) = gel(S->fu,i);
    1170        1554 :   setlg(beta, b); /* beta = [cycgen[i] if p|cyc[i], tu if p|w, fu] */
    1171        1554 :   if (DEBUGLEVEL>3) err_printf("     Beta list = %Ps\n",beta);
    1172        1554 :   primecertify(nf, beta, p, S->bad); set_avma(av);
    1173        1554 : }
    1174             : 
    1175             : static void
    1176          56 : init_bad(struct check_pr *S, GEN nf, GEN gen)
    1177             : {
    1178          56 :   long i, l = lg(gen);
    1179          56 :   GEN bad = gen_1;
    1180             : 
    1181         112 :   for (i=1; i < l; i++)
    1182          56 :     bad = lcmii(bad, gcoeff(gel(gen,i),1,1));
    1183         112 :   for (i = 1; i < l; i++)
    1184             :   {
    1185          56 :     GEN c = gel(S->cycgen,i);
    1186             :     long j;
    1187          56 :     if (typ(c) == t_MAT)
    1188             :     {
    1189          56 :       GEN g = gel(c,1);
    1190         392 :       for (j = 1; j < lg(g); j++)
    1191             :       {
    1192         336 :         GEN h = idealhnf_shallow(nf, gel(g,j));
    1193         336 :         bad = lcmii(bad, gcoeff(h,1,1));
    1194             :       }
    1195             :     }
    1196             :   }
    1197          56 :   S->bad = bad;
    1198          56 : }
    1199             : 
    1200             : long
    1201          56 : bnfcertify0(GEN bnf, long flag)
    1202             : {
    1203          56 :   pari_sp av = avma;
    1204             :   long N;
    1205             :   GEN nf, cyc, B, U;
    1206             :   ulong bound, p;
    1207             :   struct check_pr S;
    1208             :   forprime_t T;
    1209             : 
    1210          56 :   bnf = checkbnf(bnf);
    1211          56 :   nf = bnf_get_nf(bnf);
    1212          56 :   N = nf_get_degree(nf); if (N==1) return 1;
    1213          56 :   B = zimmertbound(nf_get_disc(nf), N, nf_get_r2(nf));
    1214          56 :   if (is_bigint(B))
    1215           0 :     pari_warn(warner,"Zimmert's bound is large (%Ps), certification will take a long time", B);
    1216          56 :   if (!is_pm1(nf_get_index(nf)))
    1217             :   {
    1218          35 :     GEN D = nf_get_diff(nf), L;
    1219          35 :     if (DEBUGLEVEL>1) err_printf("**** Testing Different = %Ps\n",D);
    1220          35 :     L = bnfisprincipal0(bnf, D, nf_FORCE);
    1221          35 :     if (DEBUGLEVEL>1) err_printf("     is %Ps\n", L);
    1222             :   }
    1223          56 :   if (DEBUGLEVEL)
    1224             :   {
    1225           0 :     err_printf("PHASE 1 [CLASS GROUP]: are all primes good ?\n");
    1226           0 :     err_printf("  Testing primes <= %Ps\n", B);
    1227             :   }
    1228          56 :   bnftestprimes(bnf, B);
    1229          56 :   if (flag) return 1;
    1230             : 
    1231          56 :   U = bnf_build_units(bnf);
    1232          56 :   cyc = bnf_get_cyc(bnf);
    1233          56 :   S.w = bnf_get_tuN(bnf);
    1234          56 :   S.mu = gel(U,1);
    1235          56 :   S.fu = vecslice(U,2,lg(U)-1);
    1236          56 :   S.cyc = cyc;
    1237          56 :   S.cycgen = bnf_build_cycgen(bnf);
    1238          56 :   init_bad(&S, nf, bnf_get_gen(bnf));
    1239             : 
    1240          56 :   B = bound_unit_index(bnf, S.fu);
    1241          56 :   if (DEBUGLEVEL)
    1242             :   {
    1243           0 :     err_printf("PHASE 2 [UNITS/RELATIONS]: are all primes good ?\n");
    1244           0 :     err_printf("  Testing primes <= %Ps\n", B);
    1245             :   }
    1246          56 :   bound = itou_or_0(B);
    1247          56 :   if (!bound) pari_err_OVERFLOW("bnfcertify [too many primes to check]");
    1248          56 :   if (u_forprime_init(&T, 2, bound))
    1249        1589 :     while ( (p = u_forprime_next(&T)) ) check_prime(p, nf, &S);
    1250          56 :   if (lg(cyc) > 1)
    1251             :   {
    1252          21 :     GEN f = Z_factor(cyc_get_expo(cyc)), P = gel(f,1);
    1253             :     long i;
    1254          21 :     if (DEBUGLEVEL>1) err_printf("  Primes dividing h(K)\n\n");
    1255          28 :     for (i = lg(P)-1; i; i--)
    1256             :     {
    1257          21 :       p = itou(gel(P,i)); if (p <= bound) break;
    1258           7 :       check_prime(p, nf, &S);
    1259             :     }
    1260             :   }
    1261          56 :   return gc_long(av,1);
    1262             : }
    1263             : long
    1264          28 : bnfcertify(GEN bnf) { return bnfcertify0(bnf, 0); }
    1265             : 
    1266             : /*******************************************************************/
    1267             : /*                                                                 */
    1268             : /*        RAY CLASS FIELDS: CONDUCTORS AND DISCRIMINANTS           */
    1269             : /*                                                                 */
    1270             : /*******************************************************************/
    1271             : /* \chi(gen[i]) = zeta_D^chic[i])
    1272             :  * denormalize: express chi(gen[i]) in terms of zeta_{cyc[i]} */
    1273             : GEN
    1274      162813 : char_denormalize(GEN cyc, GEN D, GEN chic)
    1275             : {
    1276      162813 :   long i, l = lg(chic);
    1277      162813 :   GEN chi = cgetg(l, t_VEC);
    1278             :   /* \chi(gen[i]) = e(chic[i] / D) = e(chi[i] / cyc[i])
    1279             :    * hence chi[i] = chic[i]cyc[i]/ D  mod cyc[i] */
    1280      615468 :   for (i = 1; i < l; ++i)
    1281             :   {
    1282      452655 :     GEN di = gel(cyc, i), t = diviiexact(mulii(di, gel(chic,i)), D);
    1283      452655 :     gel(chi, i) = modii(t, di);
    1284             :   }
    1285      162813 :   return chi;
    1286             : }
    1287             : static GEN
    1288         462 : bnrchar_i(GEN bnr, GEN g, GEN v)
    1289             : {
    1290         462 :   long i, h, l = lg(g);
    1291             :   GEN CH, D, U, U2, H, cyc, cycD, dv, dchi;
    1292         462 :   checkbnr(bnr);
    1293         462 :   switch(typ(g))
    1294             :   {
    1295             :     GEN G;
    1296          14 :     case t_VEC:
    1297          14 :       G = cgetg(l, t_MAT);
    1298          49 :       for (i = 1; i < l; i++) gel(G,i) = isprincipalray(bnr, gel(g,i));
    1299          14 :       g = G; break;
    1300         448 :     case t_MAT:
    1301         448 :       if (RgM_is_ZM(g)) break;
    1302             :     default:
    1303           0 :       pari_err_TYPE("bnrchar",g);
    1304             :   }
    1305         462 :   cyc = bnr_get_cyc(bnr);
    1306         462 :   H = ZM_hnfall_i(shallowconcat(g,diagonal_shallow(cyc)), v? &U: NULL, 1);
    1307         462 :   dv = NULL;
    1308         462 :   if (v)
    1309             :   {
    1310          28 :     GEN w = Q_remove_denom(v, &dv);
    1311          28 :     if (typ(v)!=t_VEC || lg(v)!=l || !RgV_is_ZV(w)) pari_err_TYPE("bnrchar",v);
    1312          28 :     if (!dv) v = NULL;
    1313             :     else
    1314             :     {
    1315          28 :       U = rowslice(U, 1, l-1);
    1316          28 :       w = FpV_red(ZV_ZM_mul(w, U), dv);
    1317         105 :       for (i = 1; i < l; i++)
    1318          84 :         if (signe(gel(w,i))) pari_err_TYPE("bnrchar [inconsistent values]",v);
    1319          21 :       v = vecslice(w,l,lg(w)-1);
    1320             :     }
    1321             :   }
    1322             :   /* chi defined on subgroup H, chi(H[i]) = e(v[i] / dv)
    1323             :    * unless v = NULL: chi|H = 1*/
    1324         455 :   h = itos( ZM_det_triangular(H) ); /* #(clgp/H) = number of chars */
    1325         455 :   if (h == 1) /* unique character, H = Id */
    1326             :   {
    1327           7 :     if (v)
    1328           7 :       v = char_denormalize(cyc,dv,v);
    1329             :     else
    1330           0 :       v = zerovec(lg(cyc)-1); /* trivial char */
    1331           7 :     return mkvec(v);
    1332             :   }
    1333             : 
    1334             :   /* chi defined on a subgroup of index h > 1; U H V = D diagonal,
    1335             :    * Z^#H / (H) = Z^#H / (D) ~ \oplus (Z/diZ) */
    1336         448 :   D = ZM_snfall_i(H, &U, NULL, 1);
    1337         448 :   cycD = cyc_normalize(D); gel(cycD,1) = gen_1; /* cycD[i] = d1/di */
    1338         448 :   dchi = gel(D,1);
    1339         448 :   U2 = ZM_diag_mul(cycD, U);
    1340         448 :   if (v)
    1341             :   {
    1342          14 :     GEN Ui = ZM_inv(U, NULL);
    1343          14 :     GEN Z = hnf_solve(H, ZM_mul_diag(Ui, D));
    1344          14 :     v = ZV_ZM_mul(ZV_ZM_mul(v, Z), U2);
    1345          14 :     dchi = mulii(dchi, dv);
    1346          14 :     U2 = ZM_Z_mul(U2, dv);
    1347             :   }
    1348         448 :   CH = cyc2elts(D);
    1349        1820 :   for (i = 1; i <= h; i++)
    1350             :   {
    1351        1372 :     GEN c = zv_ZM_mul(gel(CH,i), U2);
    1352        1372 :     if (v) c = ZC_add(c, v);
    1353        1372 :     gel(CH,i) = char_denormalize(cyc, dchi, c);
    1354             :   }
    1355         448 :   return CH;
    1356             : }
    1357             : GEN
    1358         462 : bnrchar(GEN bnr, GEN g, GEN v)
    1359             : {
    1360         462 :   pari_sp av = avma;
    1361         462 :   return gerepilecopy(av, bnrchar_i(bnr,g,v));
    1362             : }
    1363             : 
    1364             : /* Let bnr1, bnr2 be such that mod(bnr2) | mod(bnr1), compute surjective map
    1365             :  *   p: Cl(bnr1) ->> Cl(bnr2).
    1366             :  * Write (bnr gens) for the concatenation of the bnf [corrected by El] and bid
    1367             :  * generators; and bnr.gen for the SNF generators. Then
    1368             :  *   bnr.gen = (bnf.gen*bnr.El | bid.gen) bnr.Ui
    1369             :  *  (bnf.gen*bnr.El | bid.gen) = bnr.gen * bnr.U */
    1370             : GEN
    1371        2114 : bnrsurjection(GEN bnr1, GEN bnr2)
    1372             : {
    1373        2114 :   GEN bnf = bnr_get_bnf(bnr2), nf = bnf_get_nf(bnf);
    1374        2114 :   GEN M, U = bnr_get_U(bnr2), bid2 = bnr_get_bid(bnr2);
    1375        2114 :   GEN gen1 = bid_get_gen(bnr_get_bid(bnr1));
    1376        2114 :   GEN cyc2 = bnr_get_cyc(bnr2), e2 = cyc_get_expo(cyc2);
    1377        2114 :   long i, l = lg(bnf_get_cyc(bnf)), lb = lg(gen1);
    1378             :   /* p(bnr1.gen) = p(bnr1 gens) * bnr1.Ui
    1379             :    *             = (bnr2 gens) * P * bnr1.Ui
    1380             :    *             = bnr2.gen * (bnr2.U * P * bnr1.Ui) */
    1381             : 
    1382             :   /* p(bid1.gen) on bid2.gen */
    1383        2114 :   M = cgetg(lb, t_MAT);
    1384        9744 :   for (i = 1; i < lb; i++) gel(M,i) = ideallogmod(nf, gel(gen1,i), bid2, e2);
    1385             :   /* [U[1], U[2]] * [Id, 0; N, M] = [U[1] + U[2]*N, U[2]*M] */
    1386        2114 :   M = ZM_mul(gel(U,2), M);
    1387        2114 :   if (l > 1)
    1388             :   { /* non trivial class group */
    1389             :     /* p(bnf.gen * bnr1.El) in terms of bnf.gen * bnr2.El and bid2.gen */
    1390         805 :     GEN El2 = bnr_get_El(bnr2), El1 = bnr_get_El(bnr1);
    1391         805 :     long ngen2 = lg(bid_get_gen(bid2))-1;
    1392         805 :     if (!ngen2)
    1393         546 :       M = gel(U,1);
    1394             :     else
    1395             :     {
    1396         259 :       GEN U1 = gel(U,1), U2 = gel(U,2), T = cgetg(l, t_MAT);
    1397             :       /* T = U1 + U2 log(El2/El1) */
    1398         539 :       for (i = 1; i < l; i++)
    1399             :       { /* bnf gen in bnr1 is bnf.gen * El1 = bnf gen in bnr 2 * El1/El2 */
    1400         280 :         GEN c = gel(U1,i);
    1401         280 :         if (typ(gel(El1,i)) != t_INT) /* else El1[i] = 1 => El2[i] = 1 */
    1402             :         {
    1403         112 :           GEN z = nfdiv(nf,gel(El1,i),gel(El2,i));
    1404         112 :           c = ZC_add(c, ZM_ZC_mul(U2, ideallogmod(nf, z, bid2, e2)));
    1405             :         }
    1406         280 :         gel(T,i) = c;
    1407             :       }
    1408         259 :       M = shallowconcat(T, M);
    1409             :     }
    1410             :   }
    1411             :   /* could reduce the matrix mod cyc2 */
    1412        2114 :   return mkvec3(ZM_mul(M, bnr_get_Ui(bnr1)), bnr_get_cyc(bnr1), cyc2);
    1413             : }
    1414             : 
    1415             : /* nchi a normalized character, S a surjective map ; return S(nchi)
    1416             :  * still normalized wrt the original cyclic structure (S[2]) */
    1417             : static GEN
    1418         721 : ag_nchar_image(GEN S, GEN nchi)
    1419             : {
    1420         721 :   GEN U, M = gel(S,1), Mc = diagonal_shallow(gel(S,3));
    1421         721 :   long l = lg(M);
    1422             : 
    1423         721 :   (void)ZM_hnfall_i(shallowconcat(M, Mc), &U, 1); /* identity */
    1424         721 :   U = matslice(U,1,l-1, l,lg(U)-1);
    1425         721 :   return char_simplify(gel(nchi,1), ZV_ZM_mul(gel(nchi,2), U));
    1426             : }
    1427             : static GEN
    1428         504 : ag_char_image(GEN S, GEN chi)
    1429             : {
    1430         504 :   GEN nchi = char_normalize(chi, cyc_normalize(gel(S,2)));
    1431         504 :   GEN DC = ag_nchar_image(S, nchi);
    1432         504 :   return char_denormalize(gel(S,3), gel(DC,1), gel(DC,2));
    1433             : }
    1434             : 
    1435             : GEN
    1436         294 : bnrmap(GEN A, GEN B)
    1437             : {
    1438         294 :   pari_sp av = avma;
    1439             :   GEN KA, KB, M, c, C;
    1440         294 :   if ((KA = checkbnf_i(A)))
    1441             :   {
    1442           7 :     checkbnr(A); checkbnr(B); KB = bnr_get_bnf(B);
    1443           7 :     if (!gidentical(KA, KB))
    1444           0 :       pari_err_TYPE("bnrmap [different fields]", mkvec2(KA,KB));
    1445           7 :     return gerepilecopy(av, bnrsurjection(A,B));
    1446             :   }
    1447         287 :   if (lg(A) != 4 || typ(A) != t_VEC) pari_err_TYPE("bnrmap [not a map]", A);
    1448         280 :   M = gel(A,1); c = gel(A,2); C = gel(A,3);
    1449         280 :   if (typ(M) != t_MAT || !RgM_is_ZM(M) || typ(c) != t_VEC ||
    1450         280 :       typ(C) != t_VEC || lg(c) != lg(M) || (lg(M) > 1 && lgcols(M) != lg(C)))
    1451           0 :         pari_err_TYPE("bnrmap [not a map]", A);
    1452         280 :   switch(typ(B))
    1453             :   {
    1454           7 :     case t_INT: /* subgroup */
    1455           7 :       B = scalarmat_shallow(B, lg(C)-1);
    1456           7 :       B = ZM_hnfmodid(B, C); break;
    1457         231 :     case t_MAT: /* subgroup */
    1458         231 :       if (!RgM_is_ZM(B)) pari_err_TYPE("bnrmap [not a subgroup]", B);
    1459         224 :       B = ZM_hnfmodid(B, c); B = ag_subgroup_image(A, B); break;
    1460          21 :     case t_VEC: /* character */
    1461          21 :       if (!char_check(c, B))
    1462          14 :         pari_err_TYPE("bnrmap [not a character mod mA]", B);
    1463           7 :       B = ag_char_image(A, B); break;
    1464          21 :     case t_COL: /* discrete log mod mA */
    1465          21 :       if (lg(B) != lg(c) || !RgV_is_ZV(B))
    1466          14 :         pari_err_TYPE("bnrmap [not a discrete log]", B);
    1467           7 :       B = vecmodii(ZM_ZC_mul(M, B), C);
    1468           7 :       return gerepileupto(av, B);
    1469             :   }
    1470         231 :   return gerepilecopy(av, B);
    1471             : }
    1472             : 
    1473             : /* Given normalized chi on bnr.clgp of conductor bnrc.mod,
    1474             :  * compute primitive character chic on bnrc.clgp equivalent to chi,
    1475             :  * still normalized wrt. bnr:
    1476             :  *   chic(genc[i]) = zeta_C^chic[i]), C = cyc_normalize(bnr.cyc)[1] */
    1477             : GEN
    1478         217 : bnrchar_primitive(GEN bnr, GEN nchi, GEN bnrc)
    1479         217 : { return ag_nchar_image(bnrsurjection(bnr, bnrc), nchi); }
    1480             : 
    1481             : /* s: <gen> = Cl_f -> Cl_f2 -> 0, H subgroup of Cl_f (generators given as
    1482             :  * HNF on [gen]). Return subgroup s(H) in Cl_f2 */
    1483             : static GEN
    1484        1043 : imageofgroup(GEN bnr, GEN bnr2, GEN H)
    1485             : {
    1486        1043 :   if (!H) return diagonal_shallow(bnr_get_cyc(bnr2));
    1487         868 :   return ag_subgroup_image(bnrsurjection(bnr, bnr2), H);
    1488             : }
    1489             : GEN
    1490         497 : bnrchar_primitive_raw(GEN bnr, GEN bnrc, GEN chi)
    1491             : {
    1492         497 :   GEN S = bnrsurjection(bnr, bnrc);
    1493         497 :   return ag_char_image(S, chi);
    1494             : }
    1495             : 
    1496             : /* convert A,B,C to [bnr, H] */
    1497             : GEN
    1498         273 : ABC_to_bnr(GEN A, GEN B, GEN C, GEN *H, int gen)
    1499             : {
    1500         273 :   if (typ(A) == t_VEC)
    1501         273 :     switch(lg(A))
    1502             :     {
    1503         119 :       case 7: /* bnr */
    1504         119 :         *H = B; return A;
    1505         154 :       case 11: /* bnf */
    1506         154 :         if (!B) pari_err_TYPE("ABC_to_bnr [bnf+missing conductor]",A);
    1507         154 :         *H = C; return Buchray(A,B, gen? nf_INIT | nf_GEN: nf_INIT);
    1508             :     }
    1509           0 :   pari_err_TYPE("ABC_to_bnr",A);
    1510             :   *H = NULL; return NULL; /* LCOV_EXCL_LINE */
    1511             : }
    1512             : 
    1513             : /* OBSOLETE */
    1514             : GEN
    1515          63 : bnrconductor0(GEN A, GEN B, GEN C, long flag)
    1516             : {
    1517          63 :   pari_sp av = avma;
    1518          63 :   GEN H, bnr = ABC_to_bnr(A,B,C,&H, 0);
    1519          63 :   return gerepilecopy(av, bnrconductor_i(bnr, H, flag));
    1520             : }
    1521             : 
    1522             : long
    1523          35 : bnrisconductor0(GEN A,GEN B,GEN C)
    1524             : {
    1525          35 :   GEN H, bnr = ABC_to_bnr(A,B,C,&H, 0);
    1526          35 :   return bnrisconductor(bnr, H);
    1527             : }
    1528             : 
    1529             : static GEN
    1530      115360 : ideallog_to_bnr_i(GEN Ubid, GEN cyc, GEN z)
    1531      115360 : { return (lg(Ubid)==1)? zerocol(lg(cyc)-1): vecmodii(ZM_ZC_mul(Ubid,z), cyc); }
    1532             : /* return bnrisprincipal(bnr, (x)), assuming z = ideallog(x); allow a
    1533             :  * t_MAT for z, understood as a collection of ideallog(x_i) */
    1534             : static GEN
    1535      112007 : ideallog_to_bnr(GEN bnr, GEN z)
    1536             : {
    1537      112007 :   GEN U = gel(bnr_get_U(bnr), 2); /* bid part */
    1538      112007 :   GEN y, cyc = bnr_get_cyc(bnr);
    1539             :   long i, l;
    1540      112007 :   if (typ(z) == t_COL) return ideallog_to_bnr_i(U, cyc, z);
    1541       95445 :   y = cgetg_copy(z, &l);
    1542      194243 :   for (i = 1; i < l; i++) gel(y,i) = ideallog_to_bnr_i(U, cyc, gel(z,i));
    1543       95445 :   return y;
    1544             : }
    1545             : static GEN
    1546       95445 : bnr_log_gen_pr(GEN bnr, zlog_S *S, long e, long index)
    1547       95445 : { return ideallog_to_bnr(bnr, log_gen_pr(S, index, bnr_get_nf(bnr), e)); }
    1548             : static GEN
    1549       16562 : bnr_log_gen_arch(GEN bnr, zlog_S *S, long index)
    1550       16562 : { return ideallog_to_bnr(bnr, log_gen_arch(S, index)); }
    1551             : 
    1552             : /* A \subset H ? Allow H = NULL = trivial subgroup */
    1553             : static int
    1554      110684 : contains(GEN H, GEN A)
    1555      110684 : { return H? (hnf_solve(H, A) != NULL): gequal0(A); }
    1556             : 
    1557             : /* finite part of the conductor of H is S.P^e2*/
    1558             : static GEN
    1559        7049 : cond0_e(GEN bnr, GEN H, zlog_S *S)
    1560             : {
    1561        7049 :   long j, k, l = lg(S->k), iscond0 = S->no2;
    1562        7049 :   GEN e = S->k, e2 = cgetg(l, t_COL);
    1563       15771 :   for (k = 1; k < l; k++)
    1564             :   {
    1565       14679 :     for (j = itos(gel(e,k)); j > 0; j--)
    1566             :     {
    1567       11305 :       if (!contains(H, bnr_log_gen_pr(bnr, S, j, k))) break;
    1568        5957 :       iscond0 = 0;
    1569             :     }
    1570        8722 :     gel(e2,k) = utoi(j);
    1571             :   }
    1572        7049 :   return iscond0? NULL: e2;
    1573             : }
    1574             : /* infinite part of the conductor of H in archp form */
    1575             : static GEN
    1576        7049 : condoo_archp(GEN bnr, GEN H, zlog_S *S)
    1577             : {
    1578        7049 :   GEN archp = S->archp, archp2 = leafcopy(archp);
    1579        7049 :   long j, k, l = lg(archp);
    1580       13251 :   for (k = j = 1; k < l; k++)
    1581             :   {
    1582        6202 :     if (!contains(H, bnr_log_gen_arch(bnr, S, k)))
    1583             :     {
    1584        2982 :       archp2[j++] = archp[k];
    1585        2982 :       continue;
    1586             :     }
    1587             :   }
    1588        7049 :   if (j == l) return S->archp;
    1589        2478 :   setlg(archp2, j); return archp2;
    1590             : }
    1591             : /* MOD useless in this function */
    1592             : static GEN
    1593        3108 : bnrconductor_factored_i(GEN bnr, GEN H, long raw)
    1594             : {
    1595        3108 :   GEN nf, bid, ideal, arch, archp, e, fa, cond = NULL;
    1596             :   zlog_S S;
    1597             : 
    1598        3108 :   checkbnr(bnr);
    1599        3108 :   bid = bnr_get_bid(bnr); init_zlog(&S, bid);
    1600        3108 :   nf = bnr_get_nf(bnr);
    1601        3108 :   H = bnr_subgroup_check(bnr, H, NULL);
    1602        3108 :   e = cond0_e(bnr, H, &S); /* in terms of S.P */
    1603        3108 :   archp = condoo_archp(bnr, H, &S);
    1604        3108 :   ideal = e? factorbackprime(nf, S.P, e): bid_get_ideal(bid);
    1605        3108 :   if (archp == S.archp)
    1606             :   {
    1607        1491 :     if (!e) cond = bnr_get_mod(bnr);
    1608        1491 :     arch = bid_get_arch(bid);
    1609             :   }
    1610             :   else
    1611        1617 :     arch = indices_to_vec01(archp, nf_get_r1(nf));
    1612        3108 :   if (!cond) cond = mkvec2(ideal, arch);
    1613        3108 :   if (raw) return cond;
    1614         133 :   fa = e? famat_remove_trivial(mkmat2(S.P, e)): bid_get_fact(bid);
    1615         133 :   return mkvec2(cond, fa);
    1616             : }
    1617             : GEN
    1618         133 : bnrconductor_factored(GEN bnr, GEN H)
    1619         133 : { return bnrconductor_factored_i(bnr, H, 0); }
    1620             : GEN
    1621        2975 : bnrconductor_raw(GEN bnr, GEN H)
    1622        2975 : { return bnrconductor_factored_i(bnr, H, 1); }
    1623             : 
    1624             : /* (see bnrdisc_i). Given a bnr, and a subgroup
    1625             :  * H0 (possibly given as a character chi, in which case H0 = ker chi) of the
    1626             :  * ray class group, compute the conductor of H if flag=0. If flag > 0, compute
    1627             :  * also the corresponding H' and output
    1628             :  * if flag = 1: [[ideal,arch],[hm,cyc,gen],H']
    1629             :  * if flag = 2: [[ideal,arch],newbnr,H'] */
    1630             : GEN
    1631        3941 : bnrconductormod(GEN bnr, GEN H0, GEN MOD)
    1632             : {
    1633        3941 :   GEN nf, bid, arch, archp, bnrc, e, H, cond = NULL;
    1634             :   int ischi;
    1635             :   zlog_S S;
    1636             : 
    1637        3941 :   checkbnr(bnr);
    1638        3941 :   bid = bnr_get_bid(bnr); init_zlog(&S, bid);
    1639        3941 :   nf = bnr_get_nf(bnr);
    1640        3941 :   H = bnr_subgroup_check(bnr, H0, NULL);
    1641        3941 :   e = cond0_e(bnr, H, &S);
    1642        3941 :   archp = condoo_archp(bnr, H, &S);
    1643        3941 :   if (archp == S.archp)
    1644             :   {
    1645        3080 :     if (!e) cond = bnr_get_mod(bnr);
    1646        3080 :     arch = gel(bnr_get_mod(bnr), 2);
    1647             :   }
    1648             :   else
    1649         861 :     arch = indices_to_vec01(archp, nf_get_r1(nf));
    1650             : 
    1651             :   /* character or subgroup ? */
    1652        3941 :   ischi = H0 && typ(H0) == t_VEC;
    1653        3941 :   if (cond)
    1654             :   { /* same conductor */
    1655        2401 :     bnrc = bnr;
    1656        2401 :     if (ischi)
    1657         644 :       H = H0;
    1658        1757 :     else if (!H)
    1659        1029 :       H = diagonal_shallow(bnr_get_cyc(bnr));
    1660             :   }
    1661             :   else
    1662             :   {
    1663        1540 :     long fl = lg(bnr_get_clgp(bnr)) == 4? nf_INIT | nf_GEN: nf_INIT;
    1664        1540 :     GEN fa = famat_remove_trivial(mkmat2(S.P, e? e: S.k)), bid;
    1665        1540 :     bid = Idealstarmod(nf, mkvec2(fa, arch), nf_INIT | nf_GEN, MOD);
    1666        1540 :     bnrc = Buchraymod_i(bnr, bid, fl, MOD);
    1667        1540 :     cond = bnr_get_mod(bnrc);
    1668        1540 :     if (ischi)
    1669         497 :       H = bnrchar_primitive_raw(bnr, bnrc, H0);
    1670             :     else
    1671        1043 :       H = imageofgroup(bnr, bnrc, H);
    1672             :   }
    1673        3941 :   return mkvec3(cond, bnrc, H);
    1674             : }
    1675             : /* OBSOLETE */
    1676             : GEN
    1677          63 : bnrconductor_i(GEN bnr, GEN H, long flag)
    1678             : {
    1679             :   GEN v;
    1680          63 :   if (flag == 0) return bnrconductor_raw(bnr, H);
    1681           0 :   v = bnrconductormod(bnr, H, NULL);
    1682           0 :   if (flag == 1) gel(v,2) = bnr_get_clgp(gel(v,2));
    1683           0 :   return v;
    1684             : }
    1685             : /* OBSOLETE */
    1686             : GEN
    1687           0 : bnrconductor(GEN bnr, GEN H, long flag)
    1688             : {
    1689           0 :   pari_sp av = avma;
    1690           0 :   if (flag > 2 || flag < 0) pari_err_FLAG("bnrconductor");
    1691           0 :   return gerepilecopy(av, bnrconductor_i(bnr, H, flag));
    1692             : }
    1693             : 
    1694             : long
    1695       92001 : bnrisconductor(GEN bnr, GEN H0)
    1696             : {
    1697       92001 :   pari_sp av = avma;
    1698             :   long j, k, l;
    1699             :   GEN archp, e, H;
    1700             :   zlog_S S;
    1701             : 
    1702       92001 :   checkbnr(bnr);
    1703       92001 :   init_zlog(&S, bnr_get_bid(bnr));
    1704       92001 :   if (!S.no2) return 0;
    1705       76223 :   H = bnr_subgroup_check(bnr, H0, NULL);
    1706             : 
    1707       76223 :   archp = S.archp;
    1708       76223 :   e     = S.k; l = lg(e);
    1709      123466 :   for (k = 1; k < l; k++)
    1710             :   {
    1711       83405 :     j = itos(gel(e,k));
    1712       83405 :     if (contains(H, bnr_log_gen_pr(bnr, &S, j, k))) return gc_long(av,0);
    1713             :   }
    1714       40061 :   l = lg(archp);
    1715       44989 :   for (k = 1; k < l; k++)
    1716        9625 :     if (contains(H, bnr_log_gen_arch(bnr, &S, k))) return gc_long(av,0);
    1717       35364 :   return gc_long(av,1);
    1718             : }
    1719             : 
    1720             : /* return the norm group corresponding to the relative extension given by
    1721             :  * polrel over bnr.bnf, assuming it is abelian and the modulus of bnr is a
    1722             :  * multiple of the conductor */
    1723             : static GEN
    1724         490 : rnfnormgroup_i(GEN bnr, GEN polrel)
    1725             : {
    1726             :   long i, j, degrel, degnf, k;
    1727             :   GEN bnf, index, discnf, nf, G, detG, fa, gdegrel;
    1728             :   GEN fac, col, cnd;
    1729             :   forprime_t S;
    1730             :   ulong p;
    1731             : 
    1732         490 :   checkbnr(bnr); bnf = bnr_get_bnf(bnr);
    1733         490 :   nf = bnf_get_nf(bnf);
    1734         490 :   cnd = gel(bnr_get_mod(bnr), 1);
    1735         490 :   polrel = RgX_nffix("rnfnormgroup", nf_get_pol(nf),polrel,1);
    1736         490 :   if (!gequal1(leading_coeff(polrel)))
    1737           0 :     pari_err_IMPL("rnfnormgroup for non-monic polynomials");
    1738             : 
    1739         490 :   degrel = degpol(polrel);
    1740         490 :   if (umodiu(bnr_get_no(bnr), degrel)) return NULL;
    1741             :   /* degrel-th powers are in norm group */
    1742         483 :   gdegrel = utoipos(degrel);
    1743         483 :   G = ZV_snf_gcd(bnr_get_cyc(bnr), gdegrel);
    1744         483 :   detG = ZV_prod(G);
    1745         483 :   k = abscmpiu(detG,degrel);
    1746         483 :   if (k < 0) return NULL;
    1747         483 :   if (!k) return diagonal(G);
    1748             : 
    1749         210 :   G = diagonal_shallow(G);
    1750         210 :   discnf = nf_get_disc(nf);
    1751         210 :   index  = nf_get_index(nf);
    1752         210 :   degnf = nf_get_degree(nf);
    1753         210 :   u_forprime_init(&S, 2, ULONG_MAX);
    1754        1057 :   while ( (p = u_forprime_next(&S)) )
    1755             :   {
    1756             :     long oldf, nfa;
    1757             :     /* If all pr are unramified and have the same residue degree, p =prod pr
    1758             :      * and including last pr^f or p^f is the same, but the last isprincipal
    1759             :      * is much easier! oldf is used to track this */
    1760             : 
    1761        1057 :     if (!umodiu(index, p)) continue; /* can't be treated efficiently */
    1762             : 
    1763             :     /* primes of degree 1 are enough, and simpler */
    1764        1057 :     fa = idealprimedec_limit_f(nf, utoipos(p), 1);
    1765        1057 :     nfa = lg(fa)-1;
    1766        1057 :     if (!nfa) continue;
    1767             :     /* all primes above p included ? */
    1768        1029 :     oldf = (nfa == degnf)? -1: 0;
    1769        1862 :     for (i=1; i<=nfa; i++)
    1770             :     {
    1771        1043 :       GEN pr = gel(fa,i), pp, T, polr, modpr;
    1772             :       long f, nfac;
    1773             :       /* if pr (probably) ramified, we have to use all (non-ram) P | pr */
    1774        1372 :       if (idealval(nf,cnd,pr)) { oldf = 0; continue; }
    1775         721 :       modpr = zk_to_Fq_init(nf, &pr, &T, &pp); /* T = NULL, pp ignored */
    1776         721 :       polr = nfX_to_FqX(polrel, nf, modpr); /* in Fp[X] */
    1777         721 :       polr = ZX_to_Flx(polr, p);
    1778         721 :       if (!Flx_is_squarefree(polr, p)) { oldf = 0; continue; }
    1779             : 
    1780         700 :       fac = gel(Flx_factor(polr, p), 1);
    1781         700 :       f = degpol(gel(fac,1));
    1782         700 :       if (f == degrel) continue; /* degrel-th powers already included */
    1783         392 :       nfac = lg(fac)-1;
    1784             :       /* check decomposition of pr has Galois type */
    1785        1008 :       for (j=2; j<=nfac; j++)
    1786         826 :         if (degpol(gel(fac,j)) != f) return NULL;
    1787         385 :       if (oldf < 0) oldf = f; else if (oldf != f) oldf = 0;
    1788             : 
    1789             :       /* last prime & all pr^f, pr | p, included. Include p^f instead */
    1790         385 :       if (oldf && i == nfa && degrel == nfa*f && !umodiu(discnf, p))
    1791           0 :         pr = utoipos(p);
    1792             : 
    1793             :       /* pr^f = N P, P | pr, hence is in norm group */
    1794         385 :       col = bnrisprincipalmod(bnr,pr,gdegrel,0);
    1795         385 :       if (f > 1) col = ZC_z_mul(col, f);
    1796         385 :       G = ZM_hnf(shallowconcat(G, col));
    1797         385 :       detG = ZM_det_triangular(G);
    1798         385 :       k = abscmpiu(detG,degrel);
    1799         385 :       if (k < 0) return NULL;
    1800         385 :       if (!k) { cgiv(detG); return G; }
    1801             :     }
    1802             :   }
    1803           0 :   return NULL;
    1804             : }
    1805             : GEN
    1806          14 : rnfnormgroup(GEN bnr, GEN polrel)
    1807             : {
    1808          14 :   pari_sp av = avma;
    1809          14 :   GEN G = rnfnormgroup_i(bnr, polrel);
    1810          14 :   if (!G) { set_avma(av); return cgetg(1,t_MAT); }
    1811           7 :   return gerepileupto(av, G);
    1812             : }
    1813             : 
    1814             : GEN
    1815          21 : nf_deg1_prime(GEN nf)
    1816             : {
    1817          21 :   GEN z, T = nf_get_pol(nf), D = nf_get_disc(nf), f = nf_get_index(nf);
    1818          21 :   long degnf = degpol(T);
    1819             :   forprime_t S;
    1820             :   pari_sp av;
    1821             :   ulong p;
    1822          21 :   u_forprime_init(&S, degnf, ULONG_MAX);
    1823          21 :   av = avma;
    1824         749 :   while ( (p = u_forprime_next(&S)) )
    1825             :   {
    1826             :     ulong r;
    1827         749 :     if (!umodiu(D, p) || !umodiu(f, p)) continue;
    1828         686 :     r = Flx_oneroot(ZX_to_Flx(T,p), p);
    1829         686 :     if (r != p)
    1830             :     {
    1831          21 :       z = utoi(Fl_neg(r, p));
    1832          21 :       z = deg1pol_shallow(gen_1, z, varn(T));
    1833          21 :       return idealprimedec_kummer(nf, z, 1, utoipos(p));
    1834             :     }
    1835         665 :     set_avma(av);
    1836             :   }
    1837           0 :   return NULL;
    1838             : }
    1839             : 
    1840             : static long
    1841          70 : rnfisabelian_i(GEN nf, GEN pol)
    1842             : {
    1843             :   GEN modpr, pr, T, Tnf, pp, ro, nfL, C, a, sig, eq;
    1844             :   long i, j, l, v;
    1845             :   ulong p, k, ka;
    1846             : 
    1847          70 :   if (typ(nf) == t_POL)
    1848          63 :     Tnf = nf;
    1849             :   else {
    1850           7 :     nf = checknf(nf);
    1851           7 :     Tnf = nf_get_pol(nf);
    1852             :   }
    1853          70 :   v = varn(Tnf);
    1854          70 :   if (degpol(Tnf) != 1 && typ(pol) == t_POL && RgX_is_QX(pol)
    1855          21 :                        && rnfisabelian_i(pol_x(v), pol)) return 1;
    1856          63 :   pol = RgX_nffix("rnfisabelian",Tnf,pol,1);
    1857          63 :   eq = nf_rnfeq(nf,pol); /* init L := K[x]/(pol), nf attached to K */
    1858          63 :   C = gel(eq,1); setvarn(C, v); /* L = Q[t]/(C) */
    1859          63 :   a = gel(eq,2); setvarn(a, v); /* root of K.pol in L */
    1860          63 :   nfL = C;
    1861          63 :   ro = nfroots_if_split(&nfL, QXX_QXQ_eval(pol, a, C));
    1862          63 :   if (!ro) return 0;
    1863          42 :   l = lg(ro)-1;
    1864             :   /* small groups are abelian, as are groups of prime order */
    1865          42 :   if (l < 6 || uisprime(l)) return 1;
    1866             : 
    1867          21 :   pr = nf_deg1_prime(nfL);
    1868          21 :   modpr = nf_to_Fq_init(nfL, &pr, &T, &pp);
    1869          21 :   p = itou(pp);
    1870          21 :   k = umodiu(gel(eq,3), p);
    1871          21 :   ka = (k * itou(nf_to_Fq(nfL, a, modpr))) % p;
    1872          21 :   sig= cgetg(l+1, t_VECSMALL);
    1873             :   /* image of c = ro[1] + k a [distinguished root of C] by the l automorphisms
    1874             :    * sig[i]: ro[1] -> ro[i] */
    1875         147 :   for (i = 1; i <= l; i++)
    1876         126 :     sig[i] = Fl_add(ka, itou(nf_to_Fq(nfL, gel(ro,i), modpr)), p);
    1877          21 :   ro = Q_primpart(ro);
    1878         126 :   for (i=2; i<=l; i++) { /* start at 2, since sig[1] = identity */
    1879         105 :     gel(ro,i) = ZX_to_Flx(gel(ro,i), p);
    1880         315 :     for (j=2; j<i; j++)
    1881         420 :       if (Flx_eval(gel(ro,j), sig[i], p)
    1882         210 :        != Flx_eval(gel(ro,i), sig[j], p)) return 0;
    1883             :   }
    1884          21 :   return 1;
    1885             : }
    1886             : long
    1887          49 : rnfisabelian(GEN nf, GEN pol)
    1888          49 : { pari_sp av = avma; return gc_long(av, rnfisabelian_i(nf, pol)); }
    1889             : 
    1890             : /* Given bnf and T defining an abelian relative extension, compute the
    1891             :  * corresponding conductor and congruence subgroup. Return
    1892             :  * [cond,bnr(cond),H] where cond=[ideal,arch] is the conductor. */
    1893             : GEN
    1894         476 : rnfconductor0(GEN bnf, GEN T, long flag)
    1895             : {
    1896         476 :   pari_sp av = avma;
    1897             :   GEN D, nf, module, bnr, H, lim, Tr, MOD;
    1898             : 
    1899         476 :   bnf = checkbnf(bnf); nf = bnf_get_nf(bnf);
    1900         476 :   Tr = rnfdisc_get_T(nf, T, &lim);
    1901         476 :   T = nfX_to_monic(nf, Tr, NULL);
    1902         476 :   if (!lim)
    1903         462 :     D = rnfdisc_factored(nf, T, NULL);
    1904             :   else
    1905             :   {
    1906             :     GEN P, E, Ez;
    1907          14 :     long i, l, degT = degpol(T);
    1908          14 :     D = idealfactor_partial(nf, nfX_disc(nf, Q_primpart(Tr)), lim);
    1909          14 :     P = gel(D,1); l = lg(P);
    1910          14 :     E = gel(D,2); Ez = ZV_to_zv(E);
    1911          14 :     if (l > 1 && vecsmall_max(Ez) > 1)
    1912             :     { /* cheaply update tame primes */
    1913          77 :       for (i = 1; i < l; i++)
    1914             :       { /* v_pr(f) = 1 + \sum_{0 < i < l} g_i/g_0
    1915             :                    <= 1 + max_{i>0} g_i/(g_i-1) \sum_{0 < i < l} g_i -1
    1916             :                    <= 1 + (p/(p-1)) * v_P(e(L/K, pr)), P | pr | p */
    1917          63 :         GEN pr = gel(P,i), p = pr_get_p(pr), e = gen_1;
    1918          63 :         long q, v = z_pvalrem(degT, p, &q);
    1919          63 :         if (v)
    1920             :         { /* e = e_tame * e_wild, e_wild | p^v */
    1921          14 :           long ee, pp = itou(p);
    1922          14 :           long t = ugcd(umodiu(subiu(pr_norm(pr),1), q), q); /* e_tame | t */
    1923             :           /* upper bound for 1 + p/(p-1) * v * e(L/Q,p) */
    1924          14 :           ee = 1 + (pp * v * pr_get_e(pr) * upowuu(pp,v) * t) / (pp-1);
    1925          14 :           e = utoi(minss(ee, Ez[i]));
    1926             :         }
    1927          63 :         gel(E,i) = e;
    1928             :       }
    1929             :     }
    1930             :   }
    1931         476 :   module = mkvec2(D, identity_perm(nf_get_r1(nf)));
    1932         476 :   MOD = flag? utoipos(degpol(T)): NULL;
    1933         476 :   bnr = Buchraymod_i(bnf, module, nf_INIT|nf_GEN, MOD);
    1934         476 :   H = rnfnormgroup_i(bnr,T); if (!H) return gc_const(av,gen_0);
    1935         469 :   return gerepilecopy(av, bnrconductormod(bnr, H, MOD));
    1936             : }
    1937             : GEN
    1938           0 : rnfconductor(GEN bnf, GEN T) { return rnfconductor0(bnf, T, 0); }
    1939             : 
    1940             : static GEN
    1941          98 : prV_norms(GEN v)
    1942             : {
    1943             :   long i, l;
    1944          98 :   GEN w = cgetg_copy(v, &l);
    1945         189 :   for (i = 1; i < l; i++) gel(w,i) = pr_norm(gel(v,i));
    1946          98 :   return w;
    1947             : }
    1948             : 
    1949             : /* Given a number field bnf=bnr[1], a ray class group structure bnr, and a
    1950             :  * subgroup H (HNF form) of the ray class group, compute [n, r1, dk]
    1951             :  * attached to H. If flag & rnf_COND, abort (return NULL) if module is not the
    1952             :  * conductor. If flag & rnf_REL, return relative data, else absolute */
    1953             : static GEN
    1954         175 : bnrdisc_i(GEN bnr, GEN H, long flag)
    1955             : {
    1956         175 :   const long flcond = flag & rnf_COND;
    1957             :   GEN nf, clhray, E, ED, dk;
    1958             :   long k, d, l, n, r1;
    1959             :   zlog_S S;
    1960             : 
    1961         175 :   checkbnr(bnr);
    1962         175 :   init_zlog(&S, bnr_get_bid(bnr));
    1963         175 :   nf = bnr_get_nf(bnr);
    1964         175 :   H = bnr_subgroup_check(bnr, H, &clhray);
    1965         175 :   d = itos(clhray);
    1966         175 :   if (!H) H = diagonal_shallow(bnr_get_cyc(bnr));
    1967         175 :   E = S.k; ED = cgetg_copy(E, &l);
    1968         308 :   for (k = 1; k < l; k++)
    1969             :   {
    1970         147 :     long j, e = itos(gel(E,k)), eD = e*d;
    1971         147 :     GEN H2 = H;
    1972         266 :     for (j = e; j > 0; j--)
    1973             :     {
    1974         182 :       GEN z = bnr_log_gen_pr(bnr, &S, j, k);
    1975             :       long d2;
    1976         182 :       H2 = ZM_hnf(shallowconcat(H2, z));
    1977         182 :       d2 = itos( ZM_det_triangular(H2) );
    1978         182 :       if (flcond && j==e && d2 == d) return NULL;
    1979         168 :       if (d2 == 1) { eD -= j; break; }
    1980         119 :       eD -= d2;
    1981             :     }
    1982         133 :     gel(ED,k) = utoi(eD); /* v_{P[k]}(relative discriminant) */
    1983             :   }
    1984         161 :   l = lg(S.archp); r1 = nf_get_r1(nf);
    1985         280 :   for (k = 1; k < l; k++)
    1986             :   {
    1987         147 :     if (!contains(H, bnr_log_gen_arch(bnr, &S, k))) { r1--; continue; }
    1988          98 :     if (flcond) return NULL;
    1989             :   }
    1990             :   /* d = relative degree
    1991             :    * r1 = number of unramified real places;
    1992             :    * [P,ED] = factorization of relative discriminant */
    1993         133 :   if (flag & rnf_REL)
    1994             :   {
    1995          35 :     n  = d;
    1996          35 :     dk = factorbackprime(nf, S.P, ED);
    1997             :   }
    1998             :   else
    1999             :   {
    2000          98 :     n = d * nf_get_degree(nf);
    2001          98 :     r1= d * r1;
    2002          98 :     dk = factorback2(prV_norms(S.P), ED);
    2003          98 :     if (((n-r1)&3) == 2) dk = negi(dk); /* (2r2) mod 4 = 2: r2(relext) is odd */
    2004          98 :     dk = mulii(dk, powiu(absi_shallow(nf_get_disc(nf)), d));
    2005             :   }
    2006         133 :   return mkvec3(utoipos(n), utoi(r1), dk);
    2007             : }
    2008             : GEN
    2009         175 : bnrdisc(GEN bnr, GEN H, long flag)
    2010             : {
    2011         175 :   pari_sp av = avma;
    2012         175 :   GEN D = bnrdisc_i(bnr, H, flag);
    2013         175 :   return D? gerepilecopy(av, D): gc_const(av, gen_0);
    2014             : }
    2015             : GEN
    2016         175 : bnrdisc0(GEN A, GEN B, GEN C, long flag)
    2017             : {
    2018         175 :   GEN H, bnr = ABC_to_bnr(A,B,C,&H, 0);
    2019         175 :   return bnrdisc(bnr,H,flag);
    2020             : }
    2021             : 
    2022             : /* Given a number field bnf=bnr[1], a ray class group structure bnr and a
    2023             :  * vector chi representing a character on the generators bnr[2][3], compute
    2024             :  * the conductor of chi. */
    2025             : GEN
    2026           7 : bnrconductorofchar(GEN bnr, GEN chi)
    2027             : {
    2028           7 :   pari_sp av = avma;
    2029           7 :   return gerepilecopy(av, bnrconductor_raw(bnr, chi));
    2030             : }
    2031             : 
    2032             : /* \sum U[i]*y[i], U[i],y[i] ZM, we allow lg(y) > lg(U). */
    2033             : static GEN
    2034         910 : ZMV_mul(GEN U, GEN y)
    2035             : {
    2036         910 :   long i, l = lg(U);
    2037         910 :   GEN z = NULL;
    2038         910 :   if (l == 1) return cgetg(1,t_MAT);
    2039        2324 :   for (i = 1; i < l; i++)
    2040             :   {
    2041        1442 :     GEN u = ZM_mul(gel(U,i), gel(y,i));
    2042        1442 :     z = z? ZM_add(z, u): u;
    2043             :   }
    2044         882 :   return z;
    2045             : }
    2046             : 
    2047             : /* t = [bid,U], h = #Cl(K) */
    2048             : static GEN
    2049         910 : get_classno(GEN t, GEN h)
    2050             : {
    2051         910 :   GEN bid = gel(t,1), m = gel(t,2), cyc = bid_get_cyc(bid), U = bid_get_U(bid);
    2052         910 :   return mulii(h, ZM_det_triangular(ZM_hnfmodid(ZMV_mul(U,m), cyc)));
    2053             : }
    2054             : 
    2055             : static void
    2056          28 : chk_listBU(GEN L, const char *s) {
    2057          28 :   if (typ(L) != t_VEC) pari_err_TYPE(s,L);
    2058          28 :   if (lg(L) > 1) {
    2059          28 :     GEN z = gel(L,1);
    2060          28 :     if (typ(z) != t_VEC) pari_err_TYPE(s,z);
    2061          28 :     if (lg(z) == 1) return;
    2062          28 :     z = gel(z,1); /* [bid,U] */
    2063          28 :     if (typ(z) != t_VEC || lg(z) != 3) pari_err_TYPE(s,z);
    2064          28 :     checkbid(gel(z,1));
    2065             :   }
    2066             : }
    2067             : 
    2068             : /* Given lists of [bid, unit ideallogs], return lists of ray class numbers */
    2069             : GEN
    2070           7 : bnrclassnolist(GEN bnf,GEN L)
    2071             : {
    2072           7 :   pari_sp av = avma;
    2073           7 :   long i, l = lg(L);
    2074             :   GEN V, h;
    2075             : 
    2076           7 :   chk_listBU(L, "bnrclassnolist");
    2077           7 :   if (l == 1) return cgetg(1, t_VEC);
    2078           7 :   bnf = checkbnf(bnf);
    2079           7 :   h = bnf_get_no(bnf);
    2080           7 :   V = cgetg(l,t_VEC);
    2081         392 :   for (i = 1; i < l; i++)
    2082             :   {
    2083         385 :     GEN v, z = gel(L,i);
    2084         385 :     long j, lz = lg(z);
    2085         385 :     gel(V,i) = v = cgetg(lz,t_VEC);
    2086         826 :     for (j=1; j<lz; j++) gel(v,j) = get_classno(gel(z,j), h);
    2087             :   }
    2088           7 :   return gerepilecopy(av, V);
    2089             : }
    2090             : 
    2091             : static GEN
    2092        1484 : Lbnrclassno(GEN L, GEN fac)
    2093             : {
    2094        1484 :   long i, l = lg(L);
    2095        2184 :   for (i=1; i<l; i++)
    2096        2184 :     if (gequal(gmael(L,i,1),fac)) return gmael(L,i,2);
    2097           0 :   pari_err_BUG("Lbnrclassno");
    2098             :   return NULL; /* LCOV_EXCL_LINE */
    2099             : }
    2100             : 
    2101             : static GEN
    2102         406 : factordivexact(GEN fa1,GEN fa2)
    2103             : {
    2104             :   long i, j, k, c, l;
    2105             :   GEN P, E, P1, E1, P2, E2, p1;
    2106             : 
    2107         406 :   P1 = gel(fa1,1); E1 = gel(fa1,2); l = lg(P1);
    2108         406 :   P2 = gel(fa2,1); E2 = gel(fa2,2);
    2109         406 :   P = cgetg(l,t_COL);
    2110         406 :   E = cgetg(l,t_COL);
    2111         903 :   for (c = i = 1; i < l; i++)
    2112             :   {
    2113         497 :     j = RgV_isin(P2,gel(P1,i));
    2114         497 :     if (!j) { gel(P,c) = gel(P1,i); gel(E,c) = gel(E1,i); c++; }
    2115             :     else
    2116             :     {
    2117         497 :       p1 = subii(gel(E1,i), gel(E2,j)); k = signe(p1);
    2118         497 :       if (k < 0) pari_err_BUG("factordivexact [not exact]");
    2119         497 :       if (k > 0) { gel(P,c) = gel(P1,i); gel(E,c) = p1; c++; }
    2120             :     }
    2121             :   }
    2122         406 :   setlg(P, c);
    2123         406 :   setlg(E, c); return mkmat2(P, E);
    2124             : }
    2125             : /* remove index k */
    2126             : static GEN
    2127        1169 : factorsplice(GEN fa, long k)
    2128             : {
    2129        1169 :   GEN p = gel(fa,1), e = gel(fa,2), P, E;
    2130        1169 :   long i, l = lg(p) - 1;
    2131        1169 :   P = cgetg(l, typ(p));
    2132        1169 :   E = cgetg(l, typ(e));
    2133        1344 :   for (i=1; i<k; i++) { P[i] = p[i]; E[i] = e[i]; }
    2134        1169 :   p++; e++;
    2135        1694 :   for (   ; i<l; i++) { P[i] = p[i]; E[i] = e[i]; }
    2136        1169 :   return mkvec2(P,E);
    2137             : }
    2138             : static GEN
    2139         812 : factorpow(GEN fa, long n)
    2140             : {
    2141         812 :   if (!n) return trivial_fact();
    2142         812 :   return mkmat2(gel(fa,1), gmulsg(n, gel(fa,2)));
    2143             : }
    2144             : static GEN
    2145        1043 : factormul(GEN fa1,GEN fa2)
    2146             : {
    2147        1043 :   GEN p, pnew, e, enew, v, P, y = famat_mul_shallow(fa1,fa2);
    2148             :   long i, c, lx;
    2149             : 
    2150        1043 :   p = gel(y,1); v = indexsort(p); lx = lg(p);
    2151        1043 :   e = gel(y,2);
    2152        1043 :   pnew = vecpermute(p, v);
    2153        1043 :   enew = vecpermute(e, v);
    2154        1043 :   P = gen_0; c = 0;
    2155        2933 :   for (i=1; i<lx; i++)
    2156             :   {
    2157        1890 :     if (gequal(gel(pnew,i),P))
    2158          49 :       gel(e,c) = addii(gel(e,c),gel(enew,i));
    2159             :     else
    2160             :     {
    2161        1841 :       c++; P = gel(pnew,i);
    2162        1841 :       gel(p,c) = P;
    2163        1841 :       gel(e,c) = gel(enew,i);
    2164             :     }
    2165             :   }
    2166        1043 :   setlg(p, c+1);
    2167        1043 :   setlg(e, c+1); return y;
    2168             : }
    2169             : 
    2170             : static long
    2171         168 : get_nz(GEN bnf, GEN ideal, GEN arch, long clhray)
    2172             : {
    2173             :   GEN arch2, mod;
    2174         168 :   long nz = 0, l = lg(arch), k, clhss;
    2175         168 :   if (typ(arch) == t_VECSMALL)
    2176          14 :     arch2 = indices_to_vec01(arch,nf_get_r1(bnf_get_nf(bnf)));
    2177             :   else
    2178         154 :     arch2 = leafcopy(arch);
    2179         168 :   mod = mkvec2(ideal, arch2);
    2180         448 :   for (k = 1; k < l; k++)
    2181             :   { /* FIXME: this is wasteful. Use the same algorithm as bnrconductor */
    2182         301 :     if (signe(gel(arch2,k)))
    2183             :     {
    2184          28 :       gel(arch2,k) = gen_0; clhss = itos(bnrclassno(bnf,mod));
    2185          28 :       gel(arch2,k) = gen_1;
    2186          28 :       if (clhss == clhray) return -1;
    2187             :     }
    2188         273 :     else nz++;
    2189             :   }
    2190         147 :   return nz;
    2191             : }
    2192             : 
    2193             : static GEN
    2194         427 : get_NR1D(long Nf, long clhray, long degk, long nz, GEN fadkabs, GEN idealrel)
    2195             : {
    2196             :   long n, R1;
    2197             :   GEN dlk;
    2198         427 :   if (nz < 0) return mkvec3(gen_0,gen_0,gen_0); /*EMPTY*/
    2199         406 :   n  = clhray * degk;
    2200         406 :   R1 = clhray * nz;
    2201         406 :   dlk = factordivexact(factorpow(Z_factor(utoipos(Nf)),clhray), idealrel);
    2202             :   /* r2 odd, set dlk = -dlk */
    2203         406 :   if (((n-R1)&3)==2) dlk = factormul(to_famat_shallow(gen_m1,gen_1), dlk);
    2204         406 :   return mkvec3(utoipos(n),
    2205             :                 stoi(R1),
    2206             :                 factormul(dlk,factorpow(fadkabs,clhray)));
    2207             : }
    2208             : 
    2209             : /* t = [bid,U], h = #Cl(K) */
    2210             : static GEN
    2211         469 : get_discdata(GEN t, GEN h)
    2212             : {
    2213         469 :   GEN bid = gel(t,1), fa = bid_get_fact(bid);
    2214         469 :   GEN P = gel(fa,1), E = vec_to_vecsmall(gel(fa,2));
    2215         469 :   return mkvec3(mkvec2(P, E), (GEN)itou(get_classno(t, h)), bid_get_mod(bid));
    2216             : }
    2217             : typedef struct _disc_data {
    2218             :   long degk;
    2219             :   GEN bnf, fadk, idealrelinit, V;
    2220             : } disc_data;
    2221             : 
    2222             : static GEN
    2223         469 : get_discray(disc_data *D, GEN V, GEN z, long N)
    2224             : {
    2225         469 :   GEN idealrel = D->idealrelinit;
    2226         469 :   GEN mod = gel(z,3), Fa = gel(z,1);
    2227         469 :   GEN P = gel(Fa,1), E = gel(Fa,2);
    2228         469 :   long k, nz, clhray = z[2], lP = lg(P);
    2229         700 :   for (k=1; k<lP; k++)
    2230             :   {
    2231         546 :     GEN pr = gel(P,k), p = pr_get_p(pr);
    2232         546 :     long e, ep = E[k], f = pr_get_f(pr);
    2233         546 :     long S = 0, norm = N, Npr = upowuu(p[2],f), clhss;
    2234         798 :     for (e=1; e<=ep; e++)
    2235             :     {
    2236             :       GEN fad;
    2237         574 :       if (e < ep) { E[k] = ep-e; fad = Fa; }
    2238         462 :       else fad = factorsplice(Fa, k);
    2239         574 :       norm /= Npr;
    2240         574 :       clhss = (long)Lbnrclassno(gel(V,norm), fad);
    2241         574 :       if (e==1 && clhss==clhray) { E[k] = ep; return cgetg(1, t_VEC); }
    2242         259 :       if (clhss == 1) { S += ep-e+1; break; }
    2243         252 :       S += clhss;
    2244             :     }
    2245         231 :     E[k] = ep;
    2246         231 :     idealrel = factormul(idealrel, to_famat_shallow(p, utoi(f * S)));
    2247             :   }
    2248         154 :   nz = get_nz(D->bnf, gel(mod,1), gel(mod,2), clhray);
    2249         154 :   return get_NR1D(N, clhray, D->degk, nz, D->fadk, idealrel);
    2250             : }
    2251             : 
    2252             : /* Given a list of bids and attached unit log matrices, return the
    2253             :  * list of discrayabs. Only keep moduli which are conductors. */
    2254             : GEN
    2255          21 : discrayabslist(GEN bnf, GEN L)
    2256             : {
    2257          21 :   pari_sp av = avma;
    2258          21 :   long i, l = lg(L);
    2259             :   GEN nf, V, D, h;
    2260             :   disc_data ID;
    2261             : 
    2262          21 :   chk_listBU(L, "discrayabslist");
    2263          21 :   if (l == 1) return cgetg(1, t_VEC);
    2264          21 :   ID.bnf = bnf = checkbnf(bnf);
    2265          21 :   nf = bnf_get_nf(bnf);
    2266          21 :   h = bnf_get_no(bnf);
    2267          21 :   ID.degk = nf_get_degree(nf);
    2268          21 :   ID.fadk = absZ_factor(nf_get_disc(nf));
    2269          21 :   ID.idealrelinit = trivial_fact();
    2270          21 :   V = cgetg(l, t_VEC);
    2271          21 :   D = cgetg(l, t_VEC);
    2272         448 :   for (i = 1; i < l; i++)
    2273             :   {
    2274         427 :     GEN z = gel(L,i), v, d;
    2275         427 :     long j, lz = lg(z);
    2276         427 :     gel(V,i) = v = cgetg(lz,t_VEC);
    2277         427 :     gel(D,i) = d = cgetg(lz,t_VEC);
    2278         896 :     for (j=1; j<lz; j++) {
    2279         469 :       gel(d,j) = get_discdata(gel(z,j), h);
    2280         469 :       gel(v,j) = get_discray(&ID, D, gel(d,j), i);
    2281             :     }
    2282             :   }
    2283          21 :   return gerepilecopy(av, V);
    2284             : }
    2285             : 
    2286             : /* a zsimp is [fa, cyc, v]
    2287             :  * fa: vecsmall factorisation,
    2288             :  * cyc: ZV (concatenation of (Z_K/pr^k)^* SNFs), the generators
    2289             :  * are positive at all real places [defined implicitly by weak approximation]
    2290             :  * v: ZC (log of units on (Z_K/pr^k)^* components) */
    2291             : static GEN
    2292          28 : zsimp(void)
    2293             : {
    2294          28 :   GEN empty = cgetg(1, t_VECSMALL);
    2295          28 :   return mkvec3(mkvec2(empty,empty), cgetg(1,t_VEC), cgetg(1,t_MAT));
    2296             : }
    2297             : 
    2298             : /* fa a vecsmall factorization, append p^e */
    2299             : static GEN
    2300         175 : fasmall_append(GEN fa, long p, long e)
    2301             : {
    2302         175 :   GEN P = gel(fa,1), E = gel(fa,2);
    2303         175 :   retmkvec2(vecsmall_append(P,p), vecsmall_append(E,e));
    2304             : }
    2305             : 
    2306             : static GEN
    2307         518 : sprk_get_cyc(GEN s) { return gel(s,1); }
    2308             : 
    2309             : /* sprk = sprkinit(pr,k), b zsimp with modulus coprime to pr */
    2310             : static GEN
    2311         518 : zsimpjoin(GEN b, GEN sprk, GEN U_pr, long prcode, long e)
    2312             : {
    2313         518 :   GEN fa, cyc = sprk_get_cyc(sprk);
    2314         518 :   if (lg(gel(b,2)) == 1) /* trivial group */
    2315         343 :     fa = mkvec2(mkvecsmall(prcode),mkvecsmall(e));
    2316             :   else
    2317             :   {
    2318         175 :     fa = fasmall_append(gel(b,1), prcode, e);
    2319         175 :     cyc = shallowconcat(gel(b,2), cyc); /* no SNF ! */
    2320         175 :     U_pr = vconcat(gel(b,3),U_pr);
    2321             :   }
    2322         518 :   return mkvec3(fa, cyc, U_pr);
    2323             : }
    2324             : /* B a zsimp, sgnU = [cyc[f_oo], sgn_{f_oo}(units)] */
    2325             : static GEN
    2326          28 : bnrclassno_1(GEN B, ulong h, GEN sgnU)
    2327             : {
    2328          28 :   long lx = lg(B), j;
    2329          28 :   GEN L = cgetg(lx,t_VEC);
    2330          56 :   for (j=1; j<lx; j++)
    2331             :   {
    2332          28 :     pari_sp av = avma;
    2333          28 :     GEN b = gel(B,j), cyc = gel(b,2), qm = gel(b,3);
    2334             :     ulong z;
    2335          28 :     cyc = shallowconcat(cyc, gel(sgnU,1));
    2336          28 :     qm = vconcat(qm, gel(sgnU,2));
    2337          28 :     z = itou( mului(h, ZM_det_triangular(ZM_hnfmodid(qm, cyc))) );
    2338          28 :     set_avma(av);
    2339          28 :     gel(L,j) = mkvec2(gel(b,1), mkvecsmall(z));
    2340             :   }
    2341          28 :   return L;
    2342             : }
    2343             : 
    2344             : static void
    2345        1344 : vecselect_p(GEN A, GEN B, GEN p, long init, long lB)
    2346             : {
    2347        1344 :   long i; setlg(B, lB);
    2348        2688 :   for (i=init; i<lB; i++) B[i] = A[p[i]];
    2349        1344 : }
    2350             : /* B := p . A = row selection according to permutation p. Treat only lower
    2351             :  * right corner init x init */
    2352             : static void
    2353        1022 : rowselect_p(GEN A, GEN B, GEN p, long init)
    2354             : {
    2355        1022 :   long i, lB = lg(A), lp = lg(p);
    2356        2436 :   for (i=1; i<init; i++) setlg(B[i],lp);
    2357        2366 :   for (   ; i<lB;   i++) vecselect_p(gel(A,i),gel(B,i),p,init,lp);
    2358        1022 : }
    2359             : static ulong
    2360        1022 : hdet(ulong h, GEN m)
    2361             : {
    2362        1022 :   pari_sp av = avma;
    2363        1022 :   GEN z = mului(h, ZM_det_triangular(ZM_hnf(m)));
    2364        1022 :   return gc_ulong(av, itou(z));
    2365             : }
    2366             : static GEN
    2367        1106 : bnrclassno_all(GEN B, ulong h, GEN sgnU)
    2368             : {
    2369             :   long lx, k, kk, j, r1, jj, nba, nbarch;
    2370             :   GEN _2, L, m, H, mm, rowsel;
    2371             : 
    2372        1106 :   if (typ(sgnU) == t_VEC) return bnrclassno_1(B,h,sgnU);
    2373        1078 :   lx = lg(B); if (lx == 1) return B;
    2374             : 
    2375         371 :   r1 = nbrows(sgnU); _2 = const_vec(r1, gen_2);
    2376         371 :   L = cgetg(lx,t_VEC); nbarch = 1L<<r1;
    2377         889 :   for (j=1; j<lx; j++)
    2378             :   {
    2379         518 :     pari_sp av = avma;
    2380         518 :     GEN b = gel(B,j), cyc = gel(b,2), qm = gel(b,3);
    2381         518 :     long nc = lg(cyc)-1;
    2382             :     /* [ qm   cyc 0 ]
    2383             :      * [ sgnU  0  2 ] */
    2384         518 :     m = ZM_hnfmodid(vconcat(qm, sgnU), shallowconcat(cyc,_2));
    2385         518 :     mm = RgM_shallowcopy(m);
    2386         518 :     rowsel = cgetg(nc+r1+1,t_VECSMALL);
    2387         518 :     H = cgetg(nbarch+1,t_VECSMALL);
    2388        1540 :     for (k = 0; k < nbarch; k++)
    2389             :     {
    2390        1022 :       nba = nc+1;
    2391        2366 :       for (kk=k,jj=1; jj<=r1; jj++,kk>>=1)
    2392        1344 :         if (kk&1) rowsel[nba++] = nc + jj;
    2393        1022 :       setlg(rowsel, nba);
    2394        1022 :       rowselect_p(m, mm, rowsel, nc+1);
    2395        1022 :       H[k+1] = hdet(h, mm);
    2396             :     }
    2397         518 :     H = gerepileuptoleaf(av, H);
    2398         518 :     gel(L,j) = mkvec2(gel(b,1), H);
    2399             :   }
    2400         371 :   return L;
    2401             : }
    2402             : 
    2403             : static int
    2404          21 : is_module(GEN v)
    2405             : {
    2406          21 :   if (lg(v) != 3 || (typ(v) != t_MAT && typ(v) != t_VEC)) return 0;
    2407          21 :   return typ(gel(v,1)) == t_VECSMALL && typ(gel(v,2)) == t_VECSMALL;
    2408             : }
    2409             : GEN
    2410          21 : decodemodule(GEN nf, GEN fa)
    2411             : {
    2412             :   long n, nn, k;
    2413          21 :   pari_sp av = avma;
    2414             :   GEN G, E, id, pr;
    2415             : 
    2416          21 :   nf = checknf(nf);
    2417          21 :   if (!is_module(fa)) pari_err_TYPE("decodemodule [not a factorization]", fa);
    2418          21 :   n = nf_get_degree(nf); nn = n*n; id = NULL;
    2419          21 :   G = gel(fa,1);
    2420          21 :   E = gel(fa,2);
    2421          35 :   for (k=1; k<lg(G); k++)
    2422             :   {
    2423          14 :     long code = G[k], p = code / nn, j = (code%n)+1;
    2424          14 :     GEN P = idealprimedec(nf, utoipos(p)), e = stoi(E[k]);
    2425          14 :     if (lg(P) <= j) pari_err_BUG("decodemodule [incorrect hash code]");
    2426          14 :     pr = gel(P,j);
    2427          14 :     id = id? idealmulpowprime(nf,id, pr,e)
    2428          14 :            : idealpow(nf, pr,e);
    2429             :   }
    2430          21 :   if (!id) { set_avma(av); return matid(n); }
    2431          14 :   return gerepileupto(av,id);
    2432             : }
    2433             : 
    2434             : /* List of ray class fields. Do all from scratch, bound < 2^30. No subgroups.
    2435             :  *
    2436             :  * Output: a vector V, V[k] contains the ideals of norm k. Given such an ideal
    2437             :  * m, the component is as follows:
    2438             :  *
    2439             :  * + if arch = NULL, run through all possible archimedean parts; archs are
    2440             :  * ordered using inverse lexicographic order, [0,..,0], [1,0,..,0], [0,1,..,0],
    2441             :  * Component is [m,V] where V is a vector with 2^r1 entries, giving for each
    2442             :  * arch the triple [N,R1,D], with N, R1, D as in discrayabs; D is in factored
    2443             :  * form.
    2444             :  *
    2445             :  * + otherwise [m,N,R1,D] */
    2446             : GEN
    2447          28 : discrayabslistarch(GEN bnf, GEN arch, ulong bound)
    2448             : {
    2449          28 :   int allarch = (arch==NULL), flbou = 0;
    2450             :   long degk, j, k, l, nba, nbarch, r1, c, sqbou;
    2451          28 :   pari_sp av0 = avma,  av,  av1;
    2452             :   GEN nf, p, Z, fa, Disc, U, sgnU, EMPTY, empty, archp;
    2453             :   GEN res, Ray, discall, idealrel, idealrelinit, fadkabs, BOUND;
    2454             :   ulong i, h;
    2455             :   forprime_t S;
    2456             : 
    2457          28 :   if (bound == 0)
    2458           0 :     pari_err_DOMAIN("discrayabslistarch","bound","==",gen_0,utoi(bound));
    2459          28 :   res = discall = NULL; /* -Wall */
    2460             : 
    2461          28 :   bnf = checkbnf(bnf);
    2462          28 :   nf = bnf_get_nf(bnf);
    2463          28 :   r1 = nf_get_r1(nf);
    2464          28 :   degk = nf_get_degree(nf);
    2465          28 :   fadkabs = absZ_factor(nf_get_disc(nf));
    2466          28 :   h = itou(bnf_get_no(bnf));
    2467             : 
    2468          28 :   if (allarch)
    2469             :   {
    2470          21 :     if (r1>15) pari_err_IMPL("r1>15 in discrayabslistarch");
    2471          21 :     arch = const_vec(r1, gen_1);
    2472             :   }
    2473           7 :   else if (lg(arch)-1 != r1)
    2474           0 :     pari_err_TYPE("Idealstar [incorrect archimedean component]",arch);
    2475          28 :   U = log_prk_units_init(bnf);
    2476          28 :   archp = vec01_to_indices(arch);
    2477          28 :   nba = lg(archp)-1;
    2478          28 :   sgnU = zm_to_ZM( nfsign_units(bnf, archp, 1) );
    2479          28 :   if (!allarch) sgnU = mkvec2(const_vec(nba,gen_2), sgnU);
    2480             : 
    2481          28 :   empty = cgetg(1,t_VEC);
    2482             :   /* what follows was rewritten from Ideallist */
    2483          28 :   BOUND = utoipos(bound);
    2484          28 :   p = cgetipos(3);
    2485          28 :   u_forprime_init(&S, 2, bound);
    2486          28 :   av = avma;
    2487          28 :   sqbou = (long)sqrt((double)bound) + 1;
    2488          28 :   Z = const_vec(bound, empty);
    2489          28 :   gel(Z,1) = mkvec(zsimp());
    2490          28 :   if (DEBUGLEVEL>1) err_printf("Starting zidealstarunits computations\n");
    2491             :   /* The goal is to compute Ray (lists of bnrclassno). Z contains "zsimps",
    2492             :    * simplified bid, from which bnrclassno is easy to compute.
    2493             :    * Once p > sqbou, delete Z[i] for i > sqbou and compute directly Ray */
    2494          28 :   Ray = Z;
    2495         294 :   while ((p[2] = u_forprime_next(&S)))
    2496             :   {
    2497         266 :     if (!flbou && p[2] > sqbou)
    2498             :     {
    2499          21 :       flbou = 1;
    2500          21 :       if (DEBUGLEVEL>1) err_printf("\nStarting bnrclassno computations\n");
    2501          21 :       Z = gerepilecopy(av,Z);
    2502          21 :       Ray = cgetg(bound+1, t_VEC);
    2503         889 :       for (i=1; i<=bound; i++) gel(Ray,i) = bnrclassno_all(gel(Z,i),h,sgnU);
    2504          21 :       Z = vecslice(Z, 1, sqbou);
    2505             :     }
    2506         266 :     fa = idealprimedec_limit_norm(nf,p,BOUND);
    2507         504 :     for (j=1; j<lg(fa); j++)
    2508             :     {
    2509         238 :       GEN pr = gel(fa,j);
    2510         238 :       long prcode, f = pr_get_f(pr);
    2511         238 :       ulong q, Q = upowuu(p[2], f);
    2512             : 
    2513             :       /* p, f-1, j-1 as a single integer in "base degk" (f,j <= degk)*/
    2514         238 :       prcode = (p[2]*degk + f-1)*degk + j-1;
    2515         238 :       q = Q;
    2516             :       /* FIXME: if Q = 2, should start at l = 2 */
    2517         238 :       for (l = 1;; l++) /* Q <= bound */
    2518         105 :       {
    2519             :         ulong iQ;
    2520         343 :         GEN sprk = log_prk_init(nf, pr, l, NULL);
    2521         343 :         GEN U_pr = log_prk_units(nf, U, sprk);
    2522        1582 :         for (iQ = Q, i = 1; iQ <= bound; iQ += Q, i++)
    2523             :         {
    2524        1239 :           GEN pz, p2, p1 = gel(Z,i);
    2525        1239 :           long lz = lg(p1);
    2526        1239 :           if (lz == 1) continue;
    2527             : 
    2528         595 :           p2 = cgetg(lz,t_VEC); c = 0;
    2529        1113 :           for (k=1; k<lz; k++)
    2530             :           {
    2531         658 :             GEN z = gel(p1,k), v = gmael(z,1,1); /* primes in zsimp's fact. */
    2532         658 :             long lv = lg(v);
    2533             :             /* If z has a power of pr in its modulus, skip it */
    2534         658 :             if (i != 1 && lv > 1 && v[lv-1] == prcode) break;
    2535         518 :             gel(p2,++c) = zsimpjoin(z,sprk,U_pr,prcode,l);
    2536             :           }
    2537         595 :           setlg(p2, c+1);
    2538         595 :           pz = gel(Ray,iQ);
    2539         595 :           if (flbou) p2 = bnrclassno_all(p2,h,sgnU);
    2540         595 :           if (lg(pz) > 1) p2 = shallowconcat(pz,p2);
    2541         595 :           gel(Ray,iQ) = p2;
    2542             :         }
    2543         343 :         Q = itou_or_0( muluu(Q, q) );
    2544         343 :         if (!Q || Q > bound) break;
    2545             :       }
    2546             :     }
    2547         266 :     if (gc_needed(av,1))
    2548             :     {
    2549           0 :       if(DEBUGMEM>1) pari_warn(warnmem,"[1]: discrayabslistarch");
    2550           0 :       gerepileall(av, flbou? 2: 1, &Z, &Ray);
    2551             :     }
    2552             :   }
    2553          28 :   if (!flbou) /* occurs iff bound = 1,2,4 */
    2554             :   {
    2555           7 :     if (DEBUGLEVEL>1) err_printf("\nStarting bnrclassno computations\n");
    2556           7 :     Ray = cgetg(bound+1, t_VEC);
    2557          35 :     for (i=1; i<=bound; i++) gel(Ray,i) = bnrclassno_all(gel(Z,i),h,sgnU);
    2558             :   }
    2559          28 :   Ray = gerepilecopy(av, Ray);
    2560             : 
    2561          28 :   if (DEBUGLEVEL>1) err_printf("Starting discrayabs computations\n");
    2562          28 :   if (allarch) nbarch = 1L<<r1;
    2563             :   else
    2564             :   {
    2565           7 :     nbarch = 1;
    2566           7 :     discall = cgetg(2,t_VEC);
    2567             :   }
    2568          28 :   EMPTY = mkvec3(gen_0,gen_0,gen_0);
    2569          28 :   idealrelinit = trivial_fact();
    2570          28 :   av1 = avma;
    2571          28 :   Disc = const_vec(bound, empty);
    2572         924 :   for (i=1; i<=bound; i++)
    2573             :   {
    2574         896 :     GEN sousdisc, sous = gel(Ray,i);
    2575         896 :     long ls = lg(sous);
    2576         896 :     gel(Disc,i) = sousdisc = cgetg(ls,t_VEC);
    2577        1442 :     for (j=1; j<ls; j++)
    2578             :     {
    2579         546 :       GEN b = gel(sous,j), clhrayall = gel(b,2), Fa = gel(b,1);
    2580         546 :       GEN P = gel(Fa,1), E = gel(Fa,2);
    2581         546 :       long lP = lg(P), karch;
    2582             : 
    2583         546 :       if (allarch) discall = cgetg(nbarch+1,t_VEC);
    2584        1596 :       for (karch=0; karch<nbarch; karch++)
    2585             :       {
    2586        1050 :         long nz, clhray = clhrayall[karch+1];
    2587        1050 :         if (allarch)
    2588             :         {
    2589             :           long ka, k2;
    2590        1022 :           nba = 0;
    2591        2366 :           for (ka=karch,k=1; k<=r1; k++,ka>>=1)
    2592        1344 :             if (ka & 1) nba++;
    2593        1918 :           for (k2=1,k=1; k<=r1; k++,k2<<=1)
    2594        1190 :             if (karch&k2 && clhrayall[karch-k2+1] == clhray)
    2595         294 :               { res = EMPTY; goto STORE; }
    2596             :         }
    2597         756 :         idealrel = idealrelinit;
    2598        1078 :         for (k=1; k<lP; k++) /* cf get_discray */
    2599             :         {
    2600         805 :           long e, ep = E[k], pf = P[k] / degk, f = (pf%degk) + 1, S = 0;
    2601         805 :           ulong normi = i, Npr;
    2602         805 :           p = utoipos(pf / degk);
    2603         805 :           Npr = upowuu(p[2],f);
    2604        1204 :           for (e=1; e<=ep; e++)
    2605             :           {
    2606             :             long clhss;
    2607             :             GEN fad;
    2608         910 :             if (e < ep) { E[k] = ep-e; fad = Fa; }
    2609         707 :             else fad = factorsplice(Fa, k);
    2610         910 :             normi /= Npr;
    2611         910 :             clhss = Lbnrclassno(gel(Ray,normi),fad)[karch+1];
    2612         910 :             if (e==1 && clhss==clhray) { E[k] = ep; res = EMPTY; goto STORE; }
    2613         427 :             if (clhss == 1) { S += ep-e+1; break; }
    2614         399 :             S += clhss;
    2615             :           }
    2616         322 :           E[k] = ep;
    2617         322 :           idealrel = factormul(idealrel, to_famat_shallow(p, utoi(f * S)));
    2618             :         }
    2619         273 :         if (!allarch && nba)
    2620          14 :           nz = get_nz(bnf, decodemodule(nf,Fa), arch, clhray);
    2621             :         else
    2622         259 :           nz = r1 - nba;
    2623         273 :         res = get_NR1D(i, clhray, degk, nz, fadkabs, idealrel);
    2624        1050 : STORE:  gel(discall,karch+1) = res;
    2625             :       }
    2626         518 :       res = allarch? mkvec2(Fa, discall)
    2627         546 :                    : mkvec4(Fa, gel(res,1), gel(res,2), gel(res,3));
    2628         546 :       gel(sousdisc,j) = res;
    2629         546 :       if (gc_needed(av1,1))
    2630             :       {
    2631             :         long jj;
    2632           0 :         if(DEBUGMEM>1) pari_warn(warnmem,"[2]: discrayabslistarch");
    2633           0 :         for (jj=j+1; jj<ls; jj++) gel(sousdisc,jj) = gen_0; /* dummy */
    2634           0 :         Disc = gerepilecopy(av1, Disc);
    2635           0 :         sousdisc = gel(Disc,i);
    2636             :       }
    2637             :     }
    2638             :   }
    2639          28 :   return gerepilecopy(av0, Disc);
    2640             : }
    2641             : 
    2642             : int
    2643        1134 : subgroup_conductor_ok(GEN H, GEN L)
    2644             : { /* test conductor */
    2645        1134 :   long i, l = lg(L);
    2646        3185 :   for (i = 1; i < l; i++)
    2647        2464 :     if ( hnf_solve(H, gel(L,i)) ) return 0;
    2648         721 :   return 1;
    2649             : }
    2650             : static GEN
    2651         497 : conductor_elts(GEN bnr)
    2652             : {
    2653             :   long le, la, i, k;
    2654             :   GEN e, L;
    2655             :   zlog_S S;
    2656             : 
    2657         497 :   if (!bnrisconductor(bnr, NULL)) return NULL;
    2658         483 :   init_zlog(&S, bnr_get_bid(bnr));
    2659         483 :   e = S.k; le = lg(e); la = lg(S.archp);
    2660         483 :   L = cgetg(le + la - 1, t_VEC);
    2661         483 :   i = 1;
    2662        1036 :   for (k = 1; k < le; k++)
    2663         553 :     gel(L,i++) = bnr_log_gen_pr(bnr, &S, itos(gel(e,k)), k);
    2664        1071 :   for (k = 1; k < la; k++)
    2665         588 :     gel(L,i++) = bnr_log_gen_arch(bnr, &S, k);
    2666         483 :   return L;
    2667             : }
    2668             : 
    2669             : /* Let C a congruence group in bnr, compute its subgroups whose index is
    2670             :  * described by bound (see subgrouplist) as subgroups of Clk(bnr).
    2671             :  * Restrict to subgroups having the same conductor as bnr */
    2672             : GEN
    2673         448 : subgrouplist_cond_sub(GEN bnr, GEN C, GEN bound)
    2674             : {
    2675         448 :   pari_sp av = avma;
    2676             :   long l, i, j;
    2677         448 :   GEN D, Mr, U, T, subgrp, L, cyc = bnr_get_cyc(bnr);
    2678             : 
    2679         448 :   L = conductor_elts(bnr); if (!L) return cgetg(1,t_VEC);
    2680         448 :   Mr = diagonal_shallow(cyc);
    2681         448 :   D = ZM_snfall_i(hnf_solve(C, Mr), &U, NULL, 1);
    2682         448 :   T = ZM_mul(C, RgM_inv(U));
    2683         448 :   subgrp  = subgrouplist(D, bound);
    2684         448 :   l = lg(subgrp);
    2685         952 :   for (i = j = 1; i < l; i++)
    2686             :   {
    2687         504 :     GEN H = ZM_hnfmodid(ZM_mul(T, gel(subgrp,i)), cyc);
    2688         504 :     if (subgroup_conductor_ok(H, L)) gel(subgrp, j++) = H;
    2689             :   }
    2690         448 :   setlg(subgrp, j);
    2691         448 :   return gerepilecopy(av, subgrp);
    2692             : }
    2693             : 
    2694             : static GEN
    2695          49 : subgroupcond(GEN bnr, GEN indexbound)
    2696             : {
    2697          49 :   pari_sp av = avma;
    2698          49 :   GEN L = conductor_elts(bnr);
    2699             : 
    2700          49 :   if (!L) return cgetg(1, t_VEC);
    2701          35 :   L = subgroupcondlist(bnr_get_cyc(bnr), indexbound, L);
    2702          35 :   if (indexbound && typ(indexbound) != t_VEC)
    2703             :   { /* sort by increasing index if not single value */
    2704          14 :     long i, l = lg(L);
    2705          14 :     GEN D = cgetg(l,t_VEC);
    2706         245 :     for (i=1; i<l; i++) gel(D,i) = ZM_det_triangular(gel(L,i));
    2707          14 :     L = vecreverse( vecpermute(L, indexsort(D)) );
    2708             :   }
    2709          35 :   return gerepilecopy(av, L);
    2710             : }
    2711             : 
    2712             : GEN
    2713         119 : subgrouplist0(GEN cyc, GEN indexbound, long all)
    2714             : {
    2715         119 :   if (!all && checkbnr_i(cyc)) return subgroupcond(cyc,indexbound);
    2716          70 :   if (typ(cyc) != t_VEC || !RgV_is_ZV(cyc)) cyc = member_cyc(cyc);
    2717          63 :   return subgrouplist(cyc,indexbound);
    2718             : }
    2719             : 
    2720             : GEN
    2721          49 : bnrdisclist0(GEN bnf, GEN L, GEN arch)
    2722             : {
    2723          49 :   if (typ(L)!=t_INT) return discrayabslist(bnf,L);
    2724          28 :   return discrayabslistarch(bnf,arch,itos(L));
    2725             : }
    2726             : 
    2727             : /****************************************************************************/
    2728             : /*                                Galois action on a BNR                    */
    2729             : /****************************************************************************/
    2730             : GEN
    2731         462 : bnrautmatrix(GEN bnr, GEN aut)
    2732             : {
    2733         462 :   pari_sp av = avma;
    2734         462 :   GEN bnf = bnr_get_bnf(bnr), nf = bnf_get_nf(bnf), bid = bnr_get_bid(bnr);
    2735         462 :   GEN M, Gen = get_Gen(bnf, bid, bnr_get_El(bnr)), cyc = bnr_get_cyc(bnr);
    2736         462 :   long i, l = lg(Gen);
    2737             : 
    2738         462 :   M = cgetg(l, t_MAT); aut = nfgaloismatrix(nf, aut);
    2739             :   /* Gen = clg.gen*U, clg.gen = Gen*Ui */
    2740        2772 :   for (i = 1; i < l; i++)
    2741        2310 :     gel(M,i) = isprincipalray(bnr, nfgaloismatrixapply(nf, aut, gel(Gen,i)));
    2742         462 :   M = ZM_mul(M, bnr_get_Ui(bnr));
    2743         462 :   l = lg(M);
    2744        2310 :   for (i = 1; i < l; i++) gel(M,i) = vecmodii(gel(M,i), cyc);
    2745         462 :   return gerepilecopy(av, M);
    2746             : }
    2747             : 
    2748             : GEN
    2749         231 : bnrgaloismatrix(GEN bnr, GEN aut)
    2750             : {
    2751         231 :   checkbnr(bnr);
    2752         231 :   switch (typ(aut))
    2753             :   {
    2754           0 :     case t_POL:
    2755             :     case t_COL:
    2756           0 :       return bnrautmatrix(bnr, aut);
    2757         231 :     case t_VEC:
    2758             :     {
    2759         231 :       pari_sp av = avma;
    2760         231 :       long i, l = lg(aut);
    2761             :       GEN v;
    2762         231 :       if (l == 9)
    2763             :       {
    2764           7 :         GEN g = gal_get_gen(aut);
    2765           7 :         if (typ(g) == t_VEC) { aut = galoispermtopol(aut, g); l = lg(aut); }
    2766             :       }
    2767         231 :       v = cgetg(l, t_VEC);
    2768         693 :       for(i = 1; i < l; i++) gel(v,i) = bnrautmatrix(bnr, gel(aut,i));
    2769         231 :       return gerepileupto(av, v);
    2770             :     }
    2771           0 :     default:
    2772           0 :       pari_err_TYPE("bnrgaloismatrix", aut);
    2773             :       return NULL; /*LCOV_EXCL_LINE*/
    2774             :   }
    2775             : }
    2776             : 
    2777             : GEN
    2778        1008 : bnrgaloisapply(GEN bnr, GEN mat, GEN x)
    2779             : {
    2780        1008 :   pari_sp av=avma;
    2781             :   GEN cyc;
    2782        1008 :   checkbnr(bnr);
    2783        1008 :   cyc = bnr_get_cyc(bnr);
    2784        1008 :   if (typ(mat)!=t_MAT || !RgM_is_ZM(mat))
    2785           0 :     pari_err_TYPE("bnrgaloisapply",mat);
    2786        1008 :   if (typ(x)!=t_MAT || !RgM_is_ZM(x))
    2787           0 :     pari_err_TYPE("bnrgaloisapply",x);
    2788        1008 :   return gerepileupto(av, ZM_hnfmodid(ZM_mul(mat, x), cyc));
    2789             : }
    2790             : 
    2791             : static GEN
    2792         448 : check_bnrgal(GEN bnr, GEN M)
    2793             : {
    2794         448 :   checkbnr(bnr);
    2795         448 :   if (typ(M)==t_MAT)
    2796           0 :     return mkvec(M);
    2797         448 :   else if (typ(M)==t_VEC && lg(M)==9 && typ(gal_get_gen(M))==t_VEC)
    2798             :   {
    2799         224 :     pari_sp av = avma;
    2800         224 :     GEN V = galoispermtopol(M, gal_get_gen(M));
    2801         224 :     return gerepileupto(av, bnrgaloismatrix(bnr, V));
    2802             :   }
    2803         224 :   else if (!is_vec_t(typ(M)))
    2804           0 :     pari_err_TYPE("bnrisgalois",M);
    2805         224 :   return M;
    2806             : }
    2807             : 
    2808             : long
    2809         448 : bnrisgalois(GEN bnr, GEN M, GEN H)
    2810             : {
    2811         448 :   pari_sp av = avma;
    2812             :   long i, l;
    2813         448 :   if (typ(H)!=t_MAT || !RgM_is_ZM(H))
    2814           0 :     pari_err_TYPE("bnrisgalois",H);
    2815         448 :   M = check_bnrgal(bnr, M); l = lg(M);
    2816         616 :   for (i=1; i<l; i++)
    2817             :   {
    2818         560 :     long res = ZM_equal(bnrgaloisapply(bnr,gel(M,i), H), H);
    2819         560 :     if (!res) return gc_long(av,0);
    2820             :   }
    2821          56 :   return gc_long(av,1);
    2822             : }

Generated by: LCOV version 1.13