Testing latest pari + WASM + node.js... and it works?! Wow.
License: GPL3
ubuntu2004
/* Copyright (C) 2000 The PARI group.12This file is part of the PARI/GP package.34PARI/GP is free software; you can redistribute it and/or modify it under the5terms of the GNU General Public License as published by the Free Software6Foundation; either version 2 of the License, or (at your option) any later7version. It is distributed in the hope that it will be useful, but WITHOUT8ANY WARRANTY WHATSOEVER.910Check the License for details. You should have received a copy of it, along11with the package; see the file 'COPYING'. If not, write to the Free Software12Foundation, Inc., 51 Franklin Street, Fifth Floor, Boston, MA 02110-1301 USA. */1314/*******************************************************************/15/* */16/* RAY CLASS FIELDS */17/* */18/*******************************************************************/19#include "pari.h"20#include "paripriv.h"2122#define DEBUGLEVEL DEBUGLEVEL_bnr2324static GEN25bnr_get_El(GEN bnr) { return gel(bnr,3); }26static GEN27bnr_get_U(GEN bnr) { return gel(bnr,4); }28static GEN29bnr_get_Ui(GEN bnr) { return gmael(bnr,4,3); }3031/* faster than Buchray */32GEN33bnfnarrow(GEN bnf)34{35GEN nf, cyc, gen, Cyc, Gen, A, GD, v, w, H, invpi, L, R, u, U0, Uoo, archp, sarch;36long r1, j, l, t, RU;37pari_sp av;3839bnf = checkbnf(bnf);40nf = bnf_get_nf(bnf);41r1 = nf_get_r1(nf); if (!r1) return gcopy( bnf_get_clgp(bnf) );4243/* simplified version of nfsign_units; r1 > 0 so bnf.tu = -1 */44av = avma; archp = identity_perm(r1);45A = bnf_get_logfu(bnf); RU = lg(A)+1;46invpi = invr( mppi(nf_get_prec(nf)) );47v = cgetg(RU,t_MAT); gel(v, 1) = const_vecsmall(r1, 1); /* nfsign(-1) */48for (j=2; j<RU; j++) gel(v,j) = nfsign_from_logarch(gel(A,j-1), invpi, archp);49/* up to here */5051v = Flm_image(v, 2); t = lg(v)-1;52if (t == r1) { set_avma(av); return gcopy( bnf_get_clgp(bnf) ); }5354v = Flm_suppl(v,2); /* v = (sgn(U)|H) in GL_r1(F_2) */55H = zm_to_ZM( vecslice(v, t+1, r1) ); /* supplement H of sgn(U) */56w = rowslice(Flm_inv(v,2), t+1, r1); /* H*w*z = proj of z on H // sgn(U) */5758sarch = nfarchstar(nf, NULL, archp);59cyc = bnf_get_cyc(bnf);60gen = bnf_get_gen(bnf); l = lg(gen);61L = cgetg(l,t_MAT); GD = gmael(bnf,9,3);62for (j=1; j<l; j++)63{64GEN z = nfsign_from_logarch(gel(GD,j), invpi, archp);65gel(L,j) = zc_to_ZC( Flm_Flc_mul(w, z, 2) );66}67/* [cyc, 0; L, 2] = relation matrix for Cl_f */68R = shallowconcat(69vconcat(diagonal_shallow(cyc), L),70vconcat(zeromat(l-1, r1-t), scalarmat_shallow(gen_2,r1-t)));71Cyc = ZM_snf_group(R, NULL, &u);72U0 = rowslice(u, 1, l-1);73Uoo = ZM_mul(H, rowslice(u, l, nbrows(u)));74l = lg(Cyc); Gen = cgetg(l,t_VEC);75for (j = 1; j < l; j++)76{77GEN g = gel(U0,j), s = gel(Uoo,j);78g = (lg(g) == 1)? gen_1: Q_primpart( idealfactorback(nf,gen,g,0) );79if (!ZV_equal0(s))80{81GEN a = set_sign_mod_divisor(nf, ZV_to_Flv(s,2), gen_1, sarch);82g = is_pm1(g)? a: idealmul(nf, a, g);83}84gel(Gen,j) = g;85}86return gerepilecopy(av, mkvec3(shifti(bnf_get_no(bnf),r1-t), Cyc, Gen));87}8889/********************************************************************/90/** **/91/** REDUCTION MOD IDELE **/92/** **/93/********************************************************************/9495static GEN96compute_fact(GEN nf, GEN U, GEN gen)97{98long i, j, l = lg(U), h = lgcols(U); /* l > 1 */99GEN basecl = cgetg(l,t_VEC), G;100101G = mkvec2(NULL, trivial_fact());102for (j = 1; j < l; j++)103{104GEN z = NULL;105for (i = 1; i < h; i++)106{107GEN g, e = gcoeff(U,i,j); if (!signe(e)) continue;108109g = gel(gen,i);110if (typ(g) != t_MAT)111{112if (z)113gel(z,2) = famat_mulpow_shallow(gel(z,2), g, e);114else115z = mkvec2(NULL, to_famat_shallow(g, e));116continue;117}118gel(G,1) = g;119g = idealpowred(nf,G,e);120z = z? idealmulred(nf,z,g): g;121}122gel(z,2) = famat_reduce(gel(z,2));123gel(basecl,j) = z;124}125return basecl;126}127128static int129too_big(GEN nf, GEN bet)130{131GEN x = nfnorm(nf,bet);132switch (typ(x))133{134case t_INT: return abscmpii(x, gen_1);135case t_FRAC: return abscmpii(gel(x,1), gel(x,2));136}137pari_err_BUG("wrong type in too_big");138return 0; /* LCOV_EXCL_LINE */139}140141/* true nf; GTM 193: Algo 4.3.4. Reduce x mod divisor */142static GEN143idealmoddivisor_aux(GEN nf, GEN x, GEN f, GEN sarch)144{145pari_sp av = avma;146GEN a, A;147148if ( is_pm1(gcoeff(f,1,1)) ) /* f = 1 */149{150A = idealred(nf, mkvec2(x, gen_1));151A = nfinv(nf, gel(A,2));152}153else154{/* given coprime integral ideals x and f (f HNF), compute "small"155* G in x, such that G = 1 mod (f). GTM 193: Algo 4.3.3 */156GEN G = idealaddtoone_raw(nf, x, f);157GEN D = idealaddtoone_i(nf, idealdiv(nf,G,x), f);158A = nfdiv(nf,D,G);159}160if (too_big(nf,A) > 0) return gc_const(av, x);161a = set_sign_mod_divisor(nf, NULL, A, sarch);162if (a != A && too_big(nf,A) > 0) return gc_const(av, x);163return idealmul(nf, a, x);164}165166GEN167idealmoddivisor(GEN bnr, GEN x)168{169GEN nf = bnr_get_nf(bnr), bid = bnr_get_bid(bnr);170return idealmoddivisor_aux(nf, x, bid_get_ideal(bid), bid_get_sarch(bid));171}172173/* v_pr(L0 * cx) */174static long175fast_val(GEN L0, GEN cx, GEN pr)176{177pari_sp av = avma;178long v = typ(L0) == t_INT? 0: ZC_nfval(L0,pr);179if (cx)180{181long w = Q_pval(cx, pr_get_p(pr));182if (w) v += w * pr_get_e(pr);183}184return gc_long(av,v);185}186187/* x coprime to fZ, return y = x mod fZ, y integral */188static GEN189make_integral_Z(GEN x, GEN fZ)190{191GEN d, y = Q_remove_denom(x, &d);192if (d) y = FpC_Fp_mul(y, Fp_inv(d, fZ), fZ);193return y;194}195196/* p pi^(-1) mod f */197static GEN198get_pinvpi(GEN nf, GEN fZ, GEN p, GEN pi, GEN *v)199{200if (!*v) {201GEN invpi = nfinv(nf, pi);202*v = make_integral_Z(RgC_Rg_mul(invpi, p), mulii(p, fZ));203}204return *v;205}206/* uniformizer pi for pr, coprime to F/p */207static GEN208get_pi(GEN F, GEN pr, GEN *v)209{210if (!*v) *v = pr_uniformizer(pr, F);211return *v;212}213214/* true nf */215static GEN216bnr_grp(GEN nf, GEN U, GEN gen, GEN cyc, GEN bid)217{218GEN h = ZV_prod(cyc);219GEN f, fZ, basecl, fa, pr, t, EX, sarch, F, P, vecpi, vecpinvpi;220long i,j,l,lp;221222if (lg(U) == 1) return mkvec3(h, cyc, cgetg(1, t_VEC));223basecl = compute_fact(nf, U, gen); /* generators in factored form */224EX = gel(bid_get_cyc(bid),1); /* exponent of (O/f)^* */225f = bid_get_ideal(bid); fZ = gcoeff(f,1,1);226fa = bid_get_fact(bid);227sarch = bid_get_sarch(bid);228P = gel(fa,1); F = prV_lcm_capZ(P);229230lp = lg(P);231vecpinvpi = cgetg(lp, t_VEC);232vecpi = cgetg(lp, t_VEC);233for (i=1; i<lp; i++)234{235pr = gel(P,i);236gel(vecpi,i) = NULL; /* to be computed if needed */237gel(vecpinvpi,i) = NULL; /* to be computed if needed */238}239240l = lg(basecl);241for (i=1; i<l; i++)242{243GEN p, pi, pinvpi, dmulI, mulI, G, I, A, e, L, newL;244long la, v, k;245pari_sp av;246/* G = [I, A=famat(L,e)] is a generator, I integral */247G = gel(basecl,i);248I = gel(G,1);249A = gel(G,2); L = gel(A,1); e = gel(A,2);250/* if no reduction took place in compute_fact, everybody is still coprime251* to f + no denominators */252if (!I) { gel(basecl,i) = famat_to_nf_moddivisor(nf, L, e, bid); continue; }253if (lg(A) == 1) { gel(basecl,i) = I; continue; }254255/* compute mulI so that mulI * I coprime to f256* FIXME: use idealcoprime ??? (Less efficient. Fix idealcoprime!) */257dmulI = mulI = NULL;258for (j=1; j<lp; j++)259{260pr = gel(P,j);261v = idealval(nf, I, pr);262if (!v) continue;263p = pr_get_p(pr);264pi = get_pi(F, pr, &gel(vecpi,j));265pinvpi = get_pinvpi(nf, fZ, p, pi, &gel(vecpinvpi,j));266t = nfpow_u(nf, pinvpi, (ulong)v);267mulI = mulI? nfmuli(nf, mulI, t): t;268t = powiu(p, v);269dmulI = dmulI? mulii(dmulI, t): t;270}271272/* make all components of L coprime to f.273* Assuming (L^e * I, f) = 1, then newL^e * mulI = L^e */274la = lg(e); newL = cgetg(la, t_VEC);275for (k=1; k<la; k++)276{277GEN cx, LL = nf_to_scalar_or_basis(nf, gel(L,k));278GEN L0 = Q_primitive_part(LL, &cx); /* LL = L0*cx (faster nfval) */279for (j=1; j<lp; j++)280{281pr = gel(P,j);282v = fast_val(L0,cx, pr); /* = val_pr(LL) */283if (!v) continue;284p = pr_get_p(pr);285pi = get_pi(F, pr, &gel(vecpi,j));286if (v > 0)287{288pinvpi = get_pinvpi(nf, fZ, p, pi, &gel(vecpinvpi,j));289t = nfpow_u(nf,pinvpi, (ulong)v);290LL = nfmul(nf, LL, t);291LL = gdiv(LL, powiu(p, v));292}293else294{295t = nfpow_u(nf,pi,(ulong)(-v));296LL = nfmul(nf, LL, t);297}298}299LL = make_integral(nf,LL,f,P);300gel(newL,k) = typ(LL) == t_INT? LL: FpC_red(LL, fZ);301}302303av = avma;304/* G in nf, = L^e mod f */305G = famat_to_nf_modideal_coprime(nf, newL, e, f, EX);306if (mulI)307{308G = nfmuli(nf, G, mulI);309G = typ(G) == t_COL? ZC_hnfrem(G, ZM_Z_mul(f, dmulI))310: modii(G, mulii(fZ,dmulI));311G = RgC_Rg_div(G, dmulI);312}313G = set_sign_mod_divisor(nf,A,G,sarch);314I = idealmul(nf,I,G);315/* more or less useless, but cheap at this point */316I = idealmoddivisor_aux(nf,I,f,sarch);317gel(basecl,i) = gerepilecopy(av, I);318}319return mkvec3(h, cyc, basecl);320}321322/********************************************************************/323/** **/324/** INIT RAY CLASS GROUP **/325/** **/326/********************************************************************/327GEN328bnr_subgroup_check(GEN bnr, GEN H, GEN *pdeg)329{330GEN no = bnr_get_no(bnr);331if (H && isintzero(H)) H = NULL;332if (H)333{334GEN h, cyc = bnr_get_cyc(bnr);335switch(typ(H))336{337case t_INT:338H = scalarmat_shallow(H, lg(cyc)-1);339/* fall through */340case t_MAT:341RgM_check_ZM(H, "bnr_subgroup_check");342H = ZM_hnfmodid(H, cyc);343break;344case t_VEC:345if (char_check(cyc, H)) { H = charker(cyc, H); break; }346default: pari_err_TYPE("bnr_subgroup_check", H);347}348h = ZM_det_triangular(H);349if (equalii(h, no)) H = NULL; else no = h;350}351if (pdeg) *pdeg = no;352return H;353}354355void356bnr_subgroup_sanitize(GEN *pbnr, GEN *pH)357{358GEN D, cnd, mod, cyc, bnr = *pbnr, H = *pH;359360if (nftyp(bnr)==typ_BNF) bnr = Buchray(bnr, gen_1, nf_INIT);361else checkbnr(bnr);362cyc = bnr_get_cyc(bnr);363if (!H) mod = cyc_get_expo(cyc);364else switch(typ(H))365{366case t_INT: mod = H; break;367case t_VEC:368if (!char_check(cyc, H))369pari_err_TYPE("bnr_subgroup_sanitize [character]", H);370H = charker(cyc, H); /* character -> subgroup */371case t_MAT:372H = hnfmodid(H, cyc); /* make sure H is a left divisor of Mat(cyc) */373D = ZM_snf(H); /* structure of Cl_f / H */374mod = lg(D) == 1? gen_1: gel(D,1);375break;376default: pari_err_TYPE("bnr_subroup_sanitize [subgroup]", H);377mod = NULL;378}379cnd = bnrconductormod(bnr, H, mod);380*pbnr = gel(cnd,2); *pH = gel(cnd,3);381}382void383bnr_char_sanitize(GEN *pbnr, GEN *pchi)384{385GEN cnd, cyc, bnr = *pbnr, chi = *pchi;386387if (nftyp(bnr)==typ_BNF) bnr = Buchray(bnr, gen_1, nf_INIT);388else checkbnr(bnr);389cyc = bnr_get_cyc(bnr);390if (typ(chi) != t_VEC || !char_check(cyc, chi))391pari_err_TYPE("bnr_char_sanitize [character]", chi);392cnd = bnrconductormod(bnr, chi, charorder(cyc, chi));393*pbnr = gel(cnd,2); *pchi = gel(cnd,3);394}395396397/* c a rational content (NULL or t_INT or t_FRAC), return u*c as a ZM/d */398static GEN399ZM_content_mul(GEN u, GEN c, GEN *pd)400{401*pd = gen_1;402if (c)403{404if (typ(c) == t_FRAC) { *pd = gel(c,2); c = gel(c,1); }405if (!is_pm1(c)) u = ZM_Z_mul(u, c);406}407return u;408}409410/* bnr natural generators: bnf gens made coprime to modulus + bid gens */411static GEN412get_Gen(GEN bnf, GEN bid, GEN El)413{414GEN Gen = shallowconcat(bnf_get_gen(bnf), bid_get_gen(bid));415GEN nf = bnf_get_nf(bnf);416long i, l = lg(El);417for (i = 1; i < l; i++)418{419GEN e = gel(El,i);420if (!isint1(e)) gel(Gen,i) = idealmul(nf, gel(El,i), gel(Gen,i));421}422return Gen;423}424425static GEN426Buchraymod_i(GEN bnf, GEN module, long flag, GEN MOD)427{428GEN nf, cyc0, cyc, gen, Cyc, clg, h, logU, U, Ui, vu;429GEN bid, cycbid, H, El;430long RU, Ri, j, ngen;431const long add_gen = flag & nf_GEN;432const long do_init = flag & nf_INIT;433434if (MOD && typ(MOD) != t_INT)435pari_err_TYPE("bnrinit [incorrect cycmod]", MOD);436bnf = checkbnf(bnf);437nf = bnf_get_nf(bnf);438RU = lg(nf_get_roots(nf))-1; /* #K.futu */439El = NULL; /* gcc -Wall */440cyc = cyc0 = bnf_get_cyc(bnf);441gen = bnf_get_gen(bnf); ngen = lg(cyc)-1;442443bid = checkbid_i(module);444if (!bid) bid = Idealstarmod(nf,module,nf_GEN|nf_INIT, MOD);445cycbid = bid_get_cyc(bid);446if (MOD)447{448cyc = ZV_snf_gcd(cyc, MOD);449cycbid = ZV_snf_gcd(cycbid, MOD);450}451Ri = lg(cycbid)-1;452if (Ri || add_gen || do_init)453{454GEN fx = bid_get_fact(bid);455El = cgetg(ngen+1,t_VEC);456for (j=1; j<=ngen; j++)457{458GEN c = idealcoprimefact(nf, gel(gen,j), fx);459gel(El,j) = nf_to_scalar_or_basis(nf,c);460}461}462if (!Ri)463{464GEN hK = bnf_get_no(bnf);465clg = add_gen? mkvec3(hK, cyc, get_Gen(bnf, bid, El)): mkvec2(hK, cyc);466if (!do_init) return clg;467U = matid(ngen);468U = mkvec3(U, cgetg(1,t_MAT), U);469vu = mkvec3(cgetg(1,t_MAT), matid(RU), gen_1);470return mkvecn(6, bnf, bid, El, U, clg, vu);471}472473logU = ideallog_units0(bnf, bid, MOD);474if (do_init)475{ /* (log(Units)|D) * u = (0 | H) */476GEN c1,c2, u,u1,u2, Hi, D = shallowconcat(logU, diagonal_shallow(cycbid));477H = ZM_hnfall_i(D, &u, 1);478u1 = matslice(u, 1,RU, 1,RU);479u2 = matslice(u, 1,RU, RU+1,lg(u)-1);480/* log(Units) (u1|u2) = (0|H) (mod D), H HNF */481482u1 = ZM_lll(u1, 0.99, LLL_INPLACE);483Hi = Q_primitive_part(RgM_inv_upper(H), &c1);484u2 = ZM_mul(ZM_reducemodmatrix(u2,u1), Hi);485u2 = Q_primitive_part(u2, &c2);486u2 = ZM_content_mul(u2, mul_content(c1,c2), &c2);487vu = mkvec3(u2,u1,c2); /* u2/c2 = H^(-1) (mod Im u1) */488}489else490{491H = ZM_hnfmodid(logU, cycbid);492vu = NULL; /* -Wall */493}494if (!ngen)495h = H;496else497{498GEN L = cgetg(ngen+1, t_MAT), cycgen = bnf_build_cycgen(bnf);499for (j=1; j<=ngen; j++)500{501GEN c = gel(cycgen,j), e = gel(El,j);502if (!equali1(e)) c = famat_mulpow_shallow(c, e, gel(cyc0,j));503gel(L,j) = ideallogmod(nf, c, bid, MOD); /* = log(Gen[j]^cyc[j]) */504}505/* [cyc0, 0; -L, H] = relation matrix for generators Gen of Cl_f */506h = shallowconcat(vconcat(diagonal_shallow(cyc0), ZM_neg(L)),507vconcat(zeromat(ngen, Ri), H));508h = MOD? ZM_hnfmodid(h, MOD): ZM_hnf(h);509}510Cyc = ZM_snf_group(h, &U, &Ui);511/* Gen = clg.gen*U, clg.gen = Gen*Ui */512clg = add_gen? bnr_grp(nf, Ui, get_Gen(bnf, bid, El), Cyc, bid)513: mkvec2(ZV_prod(Cyc), Cyc);514if (!do_init) return clg;515U = mkvec3(vecslice(U, 1,ngen), vecslice(U,ngen+1,lg(U)-1), Ui);516return mkvecn(6, bnf, bid, El, U, clg, vu);517}518GEN519Buchray(GEN bnf, GEN f, long flag)520{ return Buchraymod(bnf, f, flag, NULL); }521GEN522Buchraymod(GEN bnf, GEN f, long flag, GEN MOD)523{524pari_sp av = avma;525return gerepilecopy(av, Buchraymod_i(bnf, f, flag, MOD));526}527GEN528bnrinitmod(GEN bnf, GEN f, long flag, GEN MOD)529{530switch(flag)531{532case 0: flag = nf_INIT; break;533case 1: flag = nf_INIT | nf_GEN; break;534default: pari_err_FLAG("bnrinit");535}536return Buchraymod(bnf, f, flag, MOD);537}538GEN539bnrinit0(GEN bnf, GEN ideal, long flag)540{ return bnrinitmod(bnf, ideal, flag, NULL); }541542GEN543bnrclassno(GEN bnf,GEN ideal)544{545GEN h, D, bid, cycbid;546pari_sp av = avma;547548bnf = checkbnf(bnf);549h = bnf_get_no(bnf);550bid = checkbid_i(ideal);551if (!bid) bid = Idealstar(bnf, ideal, nf_INIT);552cycbid = bid_get_cyc(bid);553if (lg(cycbid) == 1) { set_avma(av); return icopy(h); }554D = ideallog_units(bnf, bid); /* (Z_K/f)^* / units ~ Z^n / D */555D = ZM_hnfmodid(D,cycbid);556return gerepileuptoint(av, mulii(h, ZM_det_triangular(D)));557}558GEN559bnrclassno0(GEN A, GEN B, GEN C)560{561pari_sp av = avma;562GEN h, H = NULL;563/* adapted from ABC_to_bnr, avoid costly bnrinit if possible */564if (typ(A) == t_VEC)565switch(lg(A))566{567case 7: /* bnr */568checkbnr(A); H = B;569break;570case 11: /* bnf */571if (!B) pari_err_TYPE("bnrclassno [bnf+missing conductor]",A);572if (!C) return bnrclassno(A, B);573A = Buchray(A, B, nf_INIT); H = C;574break;575default: checkbnf(A);/*error*/576}577else checkbnf(A);/*error*/578579H = bnr_subgroup_check(A, H, &h);580if (!H) { set_avma(av); return icopy(h); }581return gerepileuptoint(av, h);582}583584/* ZMV_ZCV_mul for two matrices U = [Ux,Uy], it may have more components585* (ignored) and vectors x,y */586static GEN587ZM2_ZC2_mul(GEN U, GEN x, GEN y)588{589GEN Ux = gel(U,1), Uy = gel(U,2);590if (lg(Ux) == 1) return ZM_ZC_mul(Uy,y);591if (lg(Uy) == 1) return ZM_ZC_mul(Ux,x);592return ZC_add(ZM_ZC_mul(Ux,x), ZM_ZC_mul(Uy,y));593}594595GEN596bnrisprincipalmod(GEN bnr, GEN x, GEN MOD, long flag)597{598pari_sp av = avma;599GEN E, G, clgp, bnf, nf, bid, ex, cycray, alpha, El;600int trivialbid;601602checkbnr(bnr);603El = bnr_get_El(bnr);604cycray = bnr_get_cyc(bnr);605if (MOD && flag) pari_err_FLAG("bnrisprincipalmod [MOD!=NULL and flag!=0]");606if (lg(cycray) == 1 && !(flag & nf_GEN)) return cgetg(1,t_COL);607if (MOD) cycray = ZV_snf_gcd(cycray, MOD);608609bnf = bnr_get_bnf(bnr); nf = bnf_get_nf(bnf);610bid = bnr_get_bid(bnr);611trivialbid = lg(bid_get_cyc(bid)) == 1;612if (trivialbid) ex = isprincipal(bnf, x);613else614{615GEN v = bnfisprincipal0(bnf, x, nf_FORCE|nf_GENMAT);616GEN e = gel(v,1), b = gel(v,2);617long i, j = lg(e);618for (i = 1; i < j; i++) /* modify b as if bnf.gen were El*bnf.gen */619if (typ(gel(El,i)) != t_INT && signe(gel(e,i))) /* <==> != 1 */620b = famat_mulpow_shallow(b, gel(El,i), negi(gel(e,i)));621if (!MOD && !(flag & nf_GEN)) MOD = gel(cycray,1);622ex = ZM2_ZC2_mul(bnr_get_U(bnr), e, ideallogmod(nf, b, bid, MOD));623}624ex = vecmodii(ex, cycray);625if (!(flag & (nf_GEN|nf_GENMAT))) return gerepileupto(av, ex);626627/* compute generator */628E = ZC_neg(ex);629clgp = bnr_get_clgp(bnr);630if (lg(clgp) == 4)631G = abgrp_get_gen(clgp);632else633{634G = get_Gen(bnf, bid, El);635E = ZM_ZC_mul(bnr_get_Ui(bnr), E);636}637alpha = isprincipalfact(bnf, x, G, E, nf_GENMAT|nf_GEN_IF_PRINCIPAL|nf_FORCE);638if (alpha == gen_0) pari_err_BUG("isprincipalray");639if (!trivialbid)640{641GEN v = gel(bnr,6), u2 = gel(v,1), u1 = gel(v,2), du2 = gel(v,3);642GEN y = ZM_ZC_mul(u2, ideallog(nf, alpha, bid));643if (!is_pm1(du2)) y = ZC_Z_divexact(y,du2);644y = ZC_reducemodmatrix(y, u1);645if (!ZV_equal0(y))646{647GEN U = shallowcopy(bnf_build_units(bnf));648settyp(U, t_COL);649alpha = famat_div_shallow(alpha, mkmat2(U,y));650}651}652alpha = famat_reduce(alpha);653if (!(flag & nf_GENMAT)) alpha = nffactorback(nf, alpha, NULL);654return gerepilecopy(av, mkvec2(ex,alpha));655}656657GEN658bnrisprincipal(GEN bnr, GEN x, long flag)659{ return bnrisprincipalmod(bnr, x, NULL, flag); }660661GEN662isprincipalray(GEN bnr, GEN x) { return bnrisprincipal(bnr,x,0); }663GEN664isprincipalraygen(GEN bnr, GEN x) { return bnrisprincipal(bnr,x,nf_GEN); }665666/* N! / N^N * (4/pi)^r2 * sqrt(|D|) */667GEN668minkowski_bound(GEN D, long N, long r2, long prec)669{670pari_sp av = avma;671GEN c = divri(mpfactr(N,prec), powuu(N,N));672if (r2) c = mulrr(c, powru(divur(4,mppi(prec)), r2));673c = mulrr(c, gsqrt(absi_shallow(D),prec));674return gerepileuptoleaf(av, c);675}676677/* N = [K:Q] > 1, D = disc(K) */678static GEN679zimmertbound(GEN D, long N, long R2)680{681pari_sp av = avma;682GEN w;683684if (N > 20) w = minkowski_bound(D, N, R2, DEFAULTPREC);685else686{687const double c[19][11] = {688{/*2*/ 0.6931, 0.45158},689{/*3*/ 1.71733859, 1.37420604},690{/*4*/ 2.91799837, 2.50091538, 2.11943331},691{/*5*/ 4.22701425, 3.75471588, 3.31196660},692{/*6*/ 5.61209925, 5.09730381, 4.60693851, 4.14303665},693{/*7*/ 7.05406203, 6.50550021, 5.97735406, 5.47145968},694{/*8*/ 8.54052636, 7.96438858, 7.40555445, 6.86558259, 6.34608077},695{/*9*/ 10.0630022, 9.46382812, 8.87952524, 8.31139202, 7.76081149},696{/*10*/11.6153797, 10.9966020, 10.3907654, 9.79895170, 9.22232770, 8.66213267},697{/*11*/13.1930961, 12.5573772, 11.9330458, 11.3210061, 10.7222412, 10.1378082},698{/*12*/14.7926394, 14.1420915, 13.5016616, 12.8721114, 12.2542699, 11.6490374,69911.0573775},700{/*13*/16.4112395, 15.7475710, 15.0929680, 14.4480777, 13.8136054, 13.1903162,70112.5790381},702{/*14*/18.0466672, 17.3712806, 16.7040780, 16.0456127, 15.3964878, 14.7573587,70314.1289364, 13.5119848},704{/*15*/19.6970961, 19.0111606, 18.3326615, 17.6620757, 16.9999233, 16.3467686,70515.7032228, 15.0699480},706{/*16*/21.3610081, 20.6655103, 19.9768082, 19.2953176, 18.6214885, 17.9558093,70717.2988108, 16.6510652, 16.0131906},708709{/*17*/23.0371259, 22.3329066, 21.6349299, 20.9435607, 20.2591899, 19.5822454,71018.9131878, 18.2525157, 17.6007672},711712{/*18*/24.7243611, 24.0121449, 23.3056902, 22.6053167, 21.9113705, 21.2242247,71320.5442836, 19.8719830, 19.2077941, 18.5522234},714715{/*19*/26.4217792, 25.7021950, 24.9879497, 24.2793271, 23.5766321, 22.8801952,71622.1903709, 21.5075437, 20.8321263, 20.1645647},717{/*20*/28.1285704, 27.4021674, 26.6807314, 25.9645140, 25.2537867, 24.5488420,71823.8499943, 23.1575823, 22.4719720, 21.7935548, 21.1227537}719};720w = mulrr(dbltor(exp(-c[N-2][R2])), gsqrt(absi_shallow(D),DEFAULTPREC));721}722return gerepileuptoint(av, ceil_safe(w));723}724725/* return \gamma_n^n if known, an upper bound otherwise */726GEN727Hermite_bound(long n, long prec)728{729GEN h,h1;730pari_sp av;731732switch(n)733{734case 1: return gen_1;735case 2: retmkfrac(utoipos(4), utoipos(3));736case 3: return gen_2;737case 4: return utoipos(4);738case 5: return utoipos(8);739case 6: retmkfrac(utoipos(64), utoipos(3));740case 7: return utoipos(64);741case 8: return utoipos(256);742case 24: return int2n(48);743}744av = avma;745h = powru(divur(2,mppi(prec)), n);746h1 = sqrr(ggamma(sstoQ(n+4,2),prec));747return gerepileuptoleaf(av, mulrr(h,h1));748}749750/* 1 if L (= nf != Q) primitive for sure, 0 if MAYBE imprimitive (may have a751* subfield K) */752static long753isprimitive(GEN nf)754{755long p, i, l, ep, N = nf_get_degree(nf);756GEN D, fa;757758p = ucoeff(factoru(N), 1,1); /* smallest prime | N */759if (p == N) return 1; /* prime degree */760761/* N = [L:Q] = product of primes >= p, same is true for [L:K]762* d_L = t d_K^[L:K] --> check that some q^p divides d_L */763D = nf_get_disc(nf);764fa = gel(absZ_factor_limit(D,0),2); /* list of v_q(d_L). Don't check large primes */765if (mod2(D)) i = 1;766else767{ /* q = 2 */768ep = itos(gel(fa,1));769if ((ep>>1) >= p) return 0; /* 2 | d_K ==> 4 | d_K */770i = 2;771}772l = lg(fa);773for ( ; i < l; i++)774{775ep = itos(gel(fa,i));776if (ep >= p) return 0;777}778return 1;779}780781static GEN782dft_bound(void)783{784if (DEBUGLEVEL>1) err_printf("Default bound for regulator: 0.2\n");785return dbltor(0.2);786}787788static GEN789regulatorbound(GEN bnf)790{791long N, R1, R2, R;792GEN nf, dK, p1, c1;793794nf = bnf_get_nf(bnf); N = nf_get_degree(nf);795if (!isprimitive(nf)) return dft_bound();796797dK = absi_shallow(nf_get_disc(nf));798nf_get_sign(nf, &R1, &R2); R = R1+R2-1;799c1 = (!R2 && N<12)? int2n(N & (~1UL)): powuu(N,N);800if (cmpii(dK,c1) <= 0) return dft_bound();801802p1 = sqrr(glog(gdiv(dK,c1),DEFAULTPREC));803p1 = divru(gmul2n(powru(divru(mulru(p1,3),N*(N*N-1)-6*R2),R),R2), N);804p1 = sqrtr(gdiv(p1, Hermite_bound(R, DEFAULTPREC)));805if (DEBUGLEVEL>1) err_printf("Mahler bound for regulator: %Ps\n",p1);806return gmax_shallow(p1, dbltor(0.2));807}808809static int810is_unit(GEN M, long r1, GEN x)811{812pari_sp av = avma;813GEN Nx = ground( embed_norm(RgM_zc_mul(M,x), r1) );814return gc_bool(av, is_pm1(Nx));815}816817/* FIXME: should use smallvectors */818static double819minimforunits(GEN nf, long BORNE, ulong w)820{821const long prec = MEDDEFAULTPREC;822long n, r1, i, j, k, *x, cnt = 0;823pari_sp av = avma;824GEN r, M;825double p, norme, normin;826double **q,*v,*y,*z;827double eps=0.000001, BOUND = BORNE * 1.00001;828829if (DEBUGLEVEL>=2)830{831err_printf("Searching minimum of T2-form on units:\n");832if (DEBUGLEVEL>2) err_printf(" BOUND = %ld\n",BORNE);833}834n = nf_get_degree(nf); r1 = nf_get_r1(nf);835minim_alloc(n+1, &q, &x, &y, &z, &v);836M = gprec_w(nf_get_M(nf), prec);837r = gaussred_from_QR(nf_get_G(nf), prec);838for (j=1; j<=n; j++)839{840v[j] = gtodouble(gcoeff(r,j,j));841for (i=1; i<j; i++) q[i][j] = gtodouble(gcoeff(r,i,j));842}843normin = (double)BORNE*(1-eps);844k=n; y[n]=z[n]=0;845x[n] = (long)(sqrt(BOUND/v[n]));846847for(;;x[1]--)848{849do850{851if (k>1)852{853long l = k-1;854z[l] = 0;855for (j=k; j<=n; j++) z[l] += q[l][j]*x[j];856p = (double)x[k] + z[k];857y[l] = y[k] + p*p*v[k];858x[l] = (long)floor(sqrt((BOUND-y[l])/v[l])-z[l]);859k = l;860}861for(;;)862{863p = (double)x[k] + z[k];864if (y[k] + p*p*v[k] <= BOUND) break;865k++; x[k]--;866}867}868while (k>1);869if (!x[1] && y[1]<=eps) break;870871if (DEBUGLEVEL>8) err_printf(".");872if (++cnt == 5000) return -1.; /* too expensive */873874p = (double)x[1] + z[1]; norme = y[1] + p*p*v[1];875if (is_unit(M, r1, x) && norme < normin)876{877/* exclude roots of unity */878if (norme < 2*n)879{880GEN t = nfpow_u(nf, zc_to_ZC(x), w);881if (typ(t) != t_COL || ZV_isscalar(t)) continue;882}883normin = norme*(1-eps);884if (DEBUGLEVEL>=2) err_printf("*");885}886}887if (DEBUGLEVEL>=2) err_printf("\n");888set_avma(av);889return normin;890}891892#undef NBMAX893static int894is_zero(GEN x, long bitprec) { return (gexpo(x) < -bitprec); }895896static int897is_complex(GEN x, long bitprec) { return !is_zero(imag_i(x), bitprec); }898899/* assume M_star t_REAL900* FIXME: what does this do ? To be rewritten */901static GEN902compute_M0(GEN M_star,long N)903{904long m1,m2,n1,n2,n3,lr,lr1,lr2,i,j,l,vx,vy,vz,vM;905GEN pol,p1,p2,p3,p4,p5,p6,p7,p8,p9,u,v,w,r,r1,r2,M0,M0_pro,S,P,M;906GEN f1,f2,f3,g1,g2,g3,pg1,pg2,pg3,pf1,pf2,pf3,X,Y,Z;907long bitprec = 24;908909if (N == 2) return gmul2n(sqrr(gacosh(gmul2n(M_star,-1),0)), -1);910vx = fetch_var(); X = pol_x(vx);911vy = fetch_var(); Y = pol_x(vy);912vz = fetch_var(); Z = pol_x(vz);913vM = fetch_var(); M = pol_x(vM);914915M0 = NULL; m1 = N/3;916for (n1=1; n1<=m1; n1++) /* 1 <= n1 <= n2 <= n3 < N */917{918m2 = (N-n1)>>1;919for (n2=n1; n2<=m2; n2++)920{921pari_sp av = avma; n3=N-n1-n2;922if (n1==n2 && n1==n3) /* n1 = n2 = n3 = m1 = N/3 */923{924p1 = divru(M_star, m1);925p4 = sqrtr_abs( mulrr(addsr(1,p1),subrs(p1,3)) );926p5 = subrs(p1,1);927u = gen_1;928v = gmul2n(addrr(p5,p4),-1);929w = gmul2n(subrr(p5,p4),-1);930M0_pro=gmul2n(mulur(m1,addrr(sqrr(logr_abs(v)),sqrr(logr_abs(w)))), -2);931if (DEBUGLEVEL>2)932err_printf("[ %ld, %ld, %ld ]: %.28Pg\n",n1,n2,n3,M0_pro);933if (!M0 || gcmp(M0_pro,M0) < 0) M0 = M0_pro;934}935else if (n1==n2 || n2==n3)936{ /* n3 > N/3 >= n1 */937long k = N - 2*n2;938p2 = deg1pol_shallow(stoi(-n2), M_star, vx); /* M* - n2 X */939p3 = gmul(powuu(k,k),940gpowgs(gsubgs(RgX_Rg_mul(p2, M_star), k*k), n2));941pol = gsub(p3, RgX_mul(monomial(powuu(n2,n2), n2, vx),942gpowgs(p2, N-n2)));943r = roots(pol, DEFAULTPREC); lr = lg(r);944for (i=1; i<lr; i++)945{946GEN n2S;947S = real_i(gel(r,i));948if (is_complex(gel(r,i), bitprec) || signe(S) <= 0) continue;949950n2S = mulur(n2,S);951p4 = subrr(M_star, n2S);952P = divrr(mulrr(n2S,p4), subrs(mulrr(M_star,p4),k*k));953p5 = subrr(sqrr(S), gmul2n(P,2));954if (gsigne(p5) < 0) continue;955956p6 = sqrtr(p5);957v = gmul2n(subrr(S,p6),-1);958if (gsigne(v) <= 0) continue;959960u = gmul2n(addrr(S,p6),-1);961w = gpow(P, sstoQ(-n2,k), 0);962p6 = mulur(n2, addrr(sqrr(logr_abs(u)), sqrr(logr_abs(v))));963M0_pro = gmul2n(addrr(p6, mulur(k, sqrr(logr_abs(w)))),-2);964if (DEBUGLEVEL>2)965err_printf("[ %ld, %ld, %ld ]: %.28Pg\n",n1,n2,n3,M0_pro);966if (!M0 || gcmp(M0_pro,M0) < 0) M0 = M0_pro;967}968}969else970{971f1 = gsub(gadd(gmulsg(n1,X),gadd(gmulsg(n2,Y),gmulsg(n3,Z))), M);972f2 = gmulsg(n1,gmul(Y,Z));973f2 = gadd(f2,gmulsg(n2,gmul(X,Z)));974f2 = gadd(f2,gmulsg(n3,gmul(X,Y)));975f2 = gsub(f2,gmul(M,gmul(X,gmul(Y,Z))));976f3 = gsub(gmul(gpowgs(X,n1),gmul(gpowgs(Y,n2),gpowgs(Z,n3))), gen_1);977/* f1 = n1 X + n2 Y + n3 Z - M */978/* f2 = n1 YZ + n2 XZ + n3 XY */979/* f3 = X^n1 Y^n2 Z^n3 - 1*/980g1=resultant(f1,f2); g1=primpart(g1);981g2=resultant(f1,f3); g2=primpart(g2);982g3=resultant(g1,g2); g3=primpart(g3);983pf1=gsubst(f1,vM,M_star); pg1=gsubst(g1,vM,M_star);984pf2=gsubst(f2,vM,M_star); pg2=gsubst(g2,vM,M_star);985pf3=gsubst(f3,vM,M_star); pg3=gsubst(g3,vM,M_star);986/* g3 = Res_Y,Z(f1,f2,f3) */987r = roots(pg3,DEFAULTPREC); lr = lg(r);988for (i=1; i<lr; i++)989{990w = real_i(gel(r,i));991if (is_complex(gel(r,i), bitprec) || signe(w) <= 0) continue;992p1=gsubst(pg1,vz,w);993p2=gsubst(pg2,vz,w);994p3=gsubst(pf1,vz,w);995p4=gsubst(pf2,vz,w);996p5=gsubst(pf3,vz,w);997r1 = roots(p1, DEFAULTPREC); lr1 = lg(r1);998for (j=1; j<lr1; j++)999{1000v = real_i(gel(r1,j));1001if (is_complex(gel(r1,j), bitprec) || signe(v) <= 01002|| !is_zero(gsubst(p2,vy,v), bitprec)) continue;10031004p7=gsubst(p3,vy,v);1005p8=gsubst(p4,vy,v);1006p9=gsubst(p5,vy,v);1007r2 = roots(p7, DEFAULTPREC); lr2 = lg(r2);1008for (l=1; l<lr2; l++)1009{1010u = real_i(gel(r2,l));1011if (is_complex(gel(r2,l), bitprec) || signe(u) <= 01012|| !is_zero(gsubst(p8,vx,u), bitprec)1013|| !is_zero(gsubst(p9,vx,u), bitprec)) continue;10141015M0_pro = mulur(n1, sqrr(logr_abs(u)));1016M0_pro = gadd(M0_pro, mulur(n2, sqrr(logr_abs(v))));1017M0_pro = gadd(M0_pro, mulur(n3, sqrr(logr_abs(w))));1018M0_pro = gmul2n(M0_pro,-2);1019if (DEBUGLEVEL>2)1020err_printf("[ %ld, %ld, %ld ]: %.28Pg\n",n1,n2,n3,M0_pro);1021if (!M0 || gcmp(M0_pro,M0) < 0) M0 = M0_pro;1022}1023}1024}1025}1026if (!M0) set_avma(av); else M0 = gerepilecopy(av, M0);1027}1028}1029for (i=1;i<=4;i++) (void)delete_var();1030return M0? M0: gen_0;1031}10321033static GEN1034lowerboundforregulator(GEN bnf, GEN units)1035{1036long i, N, R2, RU = lg(units)-1;1037GEN nf, M0, M, G, minunit;1038double bound;10391040if (!RU) return gen_1;1041nf = bnf_get_nf(bnf);1042N = nf_get_degree(nf);1043R2 = nf_get_r2(nf);10441045G = nf_get_G(nf);1046minunit = gnorml2(RgM_RgC_mul(G, gel(units,1))); /* T2(units[1]) */1047for (i=2; i<=RU; i++)1048{1049GEN t = gnorml2(RgM_RgC_mul(G, gel(units,i)));1050if (gcmp(t,minunit) < 0) minunit = t;1051}1052if (gexpo(minunit) > 30) return NULL;10531054bound = minimforunits(nf, itos(gceil(minunit)), bnf_get_tuN(bnf));1055if (bound < 0) return NULL;1056if (DEBUGLEVEL>1) err_printf("M* = %Ps\n", dbltor(bound));1057M0 = compute_M0(dbltor(bound), N);1058if (DEBUGLEVEL>1) err_printf("M0 = %.28Pg\n",M0);1059M = gmul2n(divru(gdiv(powrs(M0,RU),Hermite_bound(RU, DEFAULTPREC)),N),R2);1060if (cmprr(M, dbltor(0.04)) < 0) return NULL;1061M = sqrtr(M);1062if (DEBUGLEVEL>1)1063err_printf("(lower bound for regulator) M = %.28Pg\n",M);1064return M;1065}10661067/* upper bound for the index of bnf.fu in the full unit group */1068static GEN1069bound_unit_index(GEN bnf, GEN units)1070{1071pari_sp av = avma;1072GEN x = lowerboundforregulator(bnf, units);1073if (!x) { set_avma(av); x = regulatorbound(bnf); }1074return gerepileuptoint(av, ground(gdiv(bnf_get_reg(bnf), x)));1075}10761077/* Compute a square matrix of rank #beta attached to a family1078* (P_i), 1<=i<=#beta, of primes s.t. N(P_i) = 1 mod p, and1079* (P_i,beta[j]) = 1 for all i,j. nf = true nf */1080static void1081primecertify(GEN nf, GEN beta, ulong p, GEN bad)1082{1083long lb = lg(beta), rmax = lb - 1;1084GEN M, vQ, L;1085ulong q;1086forprime_t T;10871088if (p == 2)1089L = cgetg(1,t_VECSMALL);1090else1091L = mkvecsmall(p);1092(void)u_forprime_arith_init(&T, 1, ULONG_MAX, 1, p);1093M = cgetg(lb,t_MAT); setlg(M,1);1094while ((q = u_forprime_next(&T)))1095{1096GEN qq, gg, og;1097long lQ, i, j;1098ulong g, m;1099if (!umodiu(bad,q)) continue;11001101qq = utoipos(q);1102vQ = idealprimedec_limit_f(nf,qq,1);1103lQ = lg(vQ); if (lQ == 1) continue;11041105/* cf rootsof1_Fl */1106g = pgener_Fl_local(q, L);1107m = (q-1) / p;1108gg = utoipos( Fl_powu(g, m, q) ); /* order p in (Z/q)^* */1109og = mkmat2(mkcol(utoi(p)), mkcol(gen_1)); /* order of g */11101111if (DEBUGLEVEL>3) err_printf(" generator of (Zk/Q)^*: %lu\n", g);1112for (i = 1; i < lQ; i++)1113{1114GEN C = cgetg(lb, t_VECSMALL);1115GEN Q = gel(vQ,i); /* degree 1 */1116GEN modpr = zkmodprinit(nf, Q);1117long r;11181119for (j = 1; j < lb; j++)1120{1121GEN t = nf_to_Fp_coprime(nf, gel(beta,j), modpr);1122t = utoipos( Fl_powu(t[2], m, q) );1123C[j] = itou( Fp_log(t, gg, og, qq) ) % p;1124}1125r = lg(M);1126gel(M,r) = C; setlg(M, r+1);1127if (Flm_rank(M, p) != r) { setlg(M,r); continue; }11281129if (DEBUGLEVEL>2)1130{1131if (DEBUGLEVEL>3)1132{1133err_printf(" prime ideal Q: %Ps\n",Q);1134err_printf(" matrix log(b_j mod Q_i): %Ps\n", M);1135}1136err_printf(" new rank: %ld\n",r);1137}1138if (r == rmax) return;1139}1140}1141pari_err_BUG("primecertify");1142}11431144struct check_pr {1145long w; /* #mu(K) */1146GEN mu; /* generator of mu(K) */1147GEN fu;1148GEN cyc;1149GEN cycgen;1150GEN bad; /* p | bad <--> p | some element occurring in cycgen */1151};11521153static void1154check_prime(ulong p, GEN nf, struct check_pr *S)1155{1156pari_sp av = avma;1157long i,b, lc = lg(S->cyc), lf = lg(S->fu);1158GEN beta = cgetg(lf+lc, t_VEC);11591160if (DEBUGLEVEL>1) err_printf(" *** testing p = %lu\n",p);1161for (b=1; b<lc; b++)1162{1163if (umodiu(gel(S->cyc,b), p)) break; /* p \nmid cyc[b] */1164if (b==1 && DEBUGLEVEL>2) err_printf(" p divides h(K)\n");1165gel(beta,b) = gel(S->cycgen,b);1166}1167if (S->w % p == 0)1168{1169if (DEBUGLEVEL>2) err_printf(" p divides w(K)\n");1170gel(beta,b++) = S->mu;1171}1172for (i=1; i<lf; i++) gel(beta,b++) = gel(S->fu,i);1173setlg(beta, b); /* beta = [cycgen[i] if p|cyc[i], tu if p|w, fu] */1174if (DEBUGLEVEL>3) err_printf(" Beta list = %Ps\n",beta);1175primecertify(nf, beta, p, S->bad); set_avma(av);1176}11771178static void1179init_bad(struct check_pr *S, GEN nf, GEN gen)1180{1181long i, l = lg(gen);1182GEN bad = gen_1;11831184for (i=1; i < l; i++)1185bad = lcmii(bad, gcoeff(gel(gen,i),1,1));1186for (i = 1; i < l; i++)1187{1188GEN c = gel(S->cycgen,i);1189long j;1190if (typ(c) == t_MAT)1191{1192GEN g = gel(c,1);1193for (j = 1; j < lg(g); j++)1194{1195GEN h = idealhnf_shallow(nf, gel(g,j));1196bad = lcmii(bad, gcoeff(h,1,1));1197}1198}1199}1200S->bad = bad;1201}12021203long1204bnfcertify0(GEN bnf, long flag)1205{1206pari_sp av = avma;1207long N;1208GEN nf, cyc, B, U;1209ulong bound, p;1210struct check_pr S;1211forprime_t T;12121213bnf = checkbnf(bnf);1214nf = bnf_get_nf(bnf);1215N = nf_get_degree(nf); if (N==1) return 1;1216B = zimmertbound(nf_get_disc(nf), N, nf_get_r2(nf));1217if (is_bigint(B))1218pari_warn(warner,"Zimmert's bound is large (%Ps), certification will take a long time", B);1219if (!is_pm1(nf_get_index(nf)))1220{1221GEN D = nf_get_diff(nf), L;1222if (DEBUGLEVEL>1) err_printf("**** Testing Different = %Ps\n",D);1223L = bnfisprincipal0(bnf, D, nf_FORCE);1224if (DEBUGLEVEL>1) err_printf(" is %Ps\n", L);1225}1226if (DEBUGLEVEL)1227{1228err_printf("PHASE 1 [CLASS GROUP]: are all primes good ?\n");1229err_printf(" Testing primes <= %Ps\n", B);1230}1231bnftestprimes(bnf, B);1232if (flag) return 1;12331234U = bnf_build_units(bnf);1235cyc = bnf_get_cyc(bnf);1236S.w = bnf_get_tuN(bnf);1237S.mu = gel(U,1);1238S.fu = vecslice(U,2,lg(U)-1);1239S.cyc = cyc;1240S.cycgen = bnf_build_cycgen(bnf);1241init_bad(&S, nf, bnf_get_gen(bnf));12421243B = bound_unit_index(bnf, S.fu);1244if (DEBUGLEVEL)1245{1246err_printf("PHASE 2 [UNITS/RELATIONS]: are all primes good ?\n");1247err_printf(" Testing primes <= %Ps\n", B);1248}1249bound = itou_or_0(B);1250if (!bound) pari_err_OVERFLOW("bnfcertify [too many primes to check]");1251if (u_forprime_init(&T, 2, bound))1252while ( (p = u_forprime_next(&T)) ) check_prime(p, nf, &S);1253if (lg(cyc) > 1)1254{1255GEN f = Z_factor(cyc_get_expo(cyc)), P = gel(f,1);1256long i;1257if (DEBUGLEVEL>1) err_printf(" Primes dividing h(K)\n\n");1258for (i = lg(P)-1; i; i--)1259{1260p = itou(gel(P,i)); if (p <= bound) break;1261check_prime(p, nf, &S);1262}1263}1264return gc_long(av,1);1265}1266long1267bnfcertify(GEN bnf) { return bnfcertify0(bnf, 0); }12681269/*******************************************************************/1270/* */1271/* RAY CLASS FIELDS: CONDUCTORS AND DISCRIMINANTS */1272/* */1273/*******************************************************************/1274/* \chi(gen[i]) = zeta_D^chic[i])1275* denormalize: express chi(gen[i]) in terms of zeta_{cyc[i]} */1276GEN1277char_denormalize(GEN cyc, GEN D, GEN chic)1278{1279long i, l = lg(chic);1280GEN chi = cgetg(l, t_VEC);1281/* \chi(gen[i]) = e(chic[i] / D) = e(chi[i] / cyc[i])1282* hence chi[i] = chic[i]cyc[i]/ D mod cyc[i] */1283for (i = 1; i < l; ++i)1284{1285GEN di = gel(cyc, i), t = diviiexact(mulii(di, gel(chic,i)), D);1286gel(chi, i) = modii(t, di);1287}1288return chi;1289}1290static GEN1291bnrchar_i(GEN bnr, GEN g, GEN v)1292{1293long i, h, l = lg(g);1294GEN CH, D, U, U2, H, cyc, cycD, dv, dchi;1295checkbnr(bnr);1296switch(typ(g))1297{1298GEN G;1299case t_VEC:1300G = cgetg(l, t_MAT);1301for (i = 1; i < l; i++) gel(G,i) = isprincipalray(bnr, gel(g,i));1302g = G; break;1303case t_MAT:1304if (RgM_is_ZM(g)) break;1305default:1306pari_err_TYPE("bnrchar",g);1307}1308cyc = bnr_get_cyc(bnr);1309H = ZM_hnfall_i(shallowconcat(g,diagonal_shallow(cyc)), v? &U: NULL, 1);1310dv = NULL;1311if (v)1312{1313GEN w = Q_remove_denom(v, &dv);1314if (typ(v)!=t_VEC || lg(v)!=l || !RgV_is_ZV(w)) pari_err_TYPE("bnrchar",v);1315if (!dv) v = NULL;1316else1317{1318U = rowslice(U, 1, l-1);1319w = FpV_red(ZV_ZM_mul(w, U), dv);1320for (i = 1; i < l; i++)1321if (signe(gel(w,i))) pari_err_TYPE("bnrchar [inconsistent values]",v);1322v = vecslice(w,l,lg(w)-1);1323}1324}1325/* chi defined on subgroup H, chi(H[i]) = e(v[i] / dv)1326* unless v = NULL: chi|H = 1*/1327h = itos( ZM_det_triangular(H) ); /* #(clgp/H) = number of chars */1328if (h == 1) /* unique character, H = Id */1329{1330if (v)1331v = char_denormalize(cyc,dv,v);1332else1333v = zerovec(lg(cyc)-1); /* trivial char */1334return mkvec(v);1335}13361337/* chi defined on a subgroup of index h > 1; U H V = D diagonal,1338* Z^#H / (H) = Z^#H / (D) ~ \oplus (Z/diZ) */1339D = ZM_snfall_i(H, &U, NULL, 1);1340cycD = cyc_normalize(D); gel(cycD,1) = gen_1; /* cycD[i] = d1/di */1341dchi = gel(D,1);1342U2 = ZM_diag_mul(cycD, U);1343if (v)1344{1345GEN Ui = ZM_inv(U, NULL);1346GEN Z = hnf_solve(H, ZM_mul_diag(Ui, D));1347v = ZV_ZM_mul(ZV_ZM_mul(v, Z), U2);1348dchi = mulii(dchi, dv);1349U2 = ZM_Z_mul(U2, dv);1350}1351CH = cyc2elts(D);1352for (i = 1; i <= h; i++)1353{1354GEN c = zv_ZM_mul(gel(CH,i), U2);1355if (v) c = ZC_add(c, v);1356gel(CH,i) = char_denormalize(cyc, dchi, c);1357}1358return CH;1359}1360GEN1361bnrchar(GEN bnr, GEN g, GEN v)1362{1363pari_sp av = avma;1364return gerepilecopy(av, bnrchar_i(bnr,g,v));1365}13661367/* Let bnr1, bnr2 be such that mod(bnr2) | mod(bnr1), compute surjective map1368* p: Cl(bnr1) ->> Cl(bnr2).1369* Write (bnr gens) for the concatenation of the bnf [corrected by El] and bid1370* generators; and bnr.gen for the SNF generators. Then1371* bnr.gen = (bnf.gen*bnr.El | bid.gen) bnr.Ui1372* (bnf.gen*bnr.El | bid.gen) = bnr.gen * bnr.U */1373GEN1374bnrsurjection(GEN bnr1, GEN bnr2)1375{1376GEN bnf = bnr_get_bnf(bnr2), nf = bnf_get_nf(bnf);1377GEN M, U = bnr_get_U(bnr2), bid2 = bnr_get_bid(bnr2);1378GEN gen1 = bid_get_gen(bnr_get_bid(bnr1));1379GEN cyc2 = bnr_get_cyc(bnr2), e2 = cyc_get_expo(cyc2);1380long i, l = lg(bnf_get_cyc(bnf)), lb = lg(gen1);1381/* p(bnr1.gen) = p(bnr1 gens) * bnr1.Ui1382* = (bnr2 gens) * P * bnr1.Ui1383* = bnr2.gen * (bnr2.U * P * bnr1.Ui) */13841385/* p(bid1.gen) on bid2.gen */1386M = cgetg(lb, t_MAT);1387for (i = 1; i < lb; i++) gel(M,i) = ideallogmod(nf, gel(gen1,i), bid2, e2);1388/* [U[1], U[2]] * [Id, 0; N, M] = [U[1] + U[2]*N, U[2]*M] */1389M = ZM_mul(gel(U,2), M);1390if (l > 1)1391{ /* non trivial class group */1392/* p(bnf.gen * bnr1.El) in terms of bnf.gen * bnr2.El and bid2.gen */1393GEN El2 = bnr_get_El(bnr2), El1 = bnr_get_El(bnr1);1394long ngen2 = lg(bid_get_gen(bid2))-1;1395if (!ngen2)1396M = gel(U,1);1397else1398{1399GEN U1 = gel(U,1), U2 = gel(U,2), T = cgetg(l, t_MAT);1400/* T = U1 + U2 log(El2/El1) */1401for (i = 1; i < l; i++)1402{ /* bnf gen in bnr1 is bnf.gen * El1 = bnf gen in bnr 2 * El1/El2 */1403GEN c = gel(U1,i);1404if (typ(gel(El1,i)) != t_INT) /* else El1[i] = 1 => El2[i] = 1 */1405{1406GEN z = nfdiv(nf,gel(El1,i),gel(El2,i));1407c = ZC_add(c, ZM_ZC_mul(U2, ideallogmod(nf, z, bid2, e2)));1408}1409gel(T,i) = c;1410}1411M = shallowconcat(T, M);1412}1413}1414/* could reduce the matrix mod cyc2 */1415return mkvec3(ZM_mul(M, bnr_get_Ui(bnr1)), bnr_get_cyc(bnr1), cyc2);1416}14171418/* nchi a normalized character, S a surjective map ; return S(nchi)1419* still normalized wrt the original cyclic structure (S[2]) */1420static GEN1421ag_nchar_image(GEN S, GEN nchi)1422{1423GEN U, M = gel(S,1), Mc = diagonal_shallow(gel(S,3));1424long l = lg(M);14251426(void)ZM_hnfall_i(shallowconcat(M, Mc), &U, 1); /* identity */1427U = matslice(U,1,l-1, l,lg(U)-1);1428return char_simplify(gel(nchi,1), ZV_ZM_mul(gel(nchi,2), U));1429}1430static GEN1431ag_char_image(GEN S, GEN chi)1432{1433GEN nchi = char_normalize(chi, cyc_normalize(gel(S,2)));1434GEN DC = ag_nchar_image(S, nchi);1435return char_denormalize(gel(S,3), gel(DC,1), gel(DC,2));1436}14371438GEN1439bnrmap(GEN A, GEN B)1440{1441pari_sp av = avma;1442GEN KA, KB, M, c, C;1443if ((KA = checkbnf_i(A)))1444{1445checkbnr(A); checkbnr(B); KB = bnr_get_bnf(B);1446if (!gidentical(KA, KB))1447pari_err_TYPE("bnrmap [different fields]", mkvec2(KA,KB));1448return gerepilecopy(av, bnrsurjection(A,B));1449}1450if (lg(A) != 4 || typ(A) != t_VEC) pari_err_TYPE("bnrmap [not a map]", A);1451M = gel(A,1); c = gel(A,2); C = gel(A,3);1452if (typ(M) != t_MAT || !RgM_is_ZM(M) || typ(c) != t_VEC ||1453typ(C) != t_VEC || lg(c) != lg(M) || (lg(M) > 1 && lgcols(M) != lg(C)))1454pari_err_TYPE("bnrmap [not a map]", A);1455switch(typ(B))1456{1457case t_INT: /* subgroup */1458B = scalarmat_shallow(B, lg(C)-1);1459B = ZM_hnfmodid(B, C); break;1460case t_MAT: /* subgroup */1461if (!RgM_is_ZM(B)) pari_err_TYPE("bnrmap [not a subgroup]", B);1462B = ZM_hnfmodid(B, c); B = ag_subgroup_image(A, B); break;1463case t_VEC: /* character */1464if (!char_check(c, B))1465pari_err_TYPE("bnrmap [not a character mod mA]", B);1466B = ag_char_image(A, B); break;1467case t_COL: /* discrete log mod mA */1468if (lg(B) != lg(c) || !RgV_is_ZV(B))1469pari_err_TYPE("bnrmap [not a discrete log]", B);1470B = vecmodii(ZM_ZC_mul(M, B), C);1471return gerepileupto(av, B);1472}1473return gerepilecopy(av, B);1474}14751476/* Given normalized chi on bnr.clgp of conductor bnrc.mod,1477* compute primitive character chic on bnrc.clgp equivalent to chi,1478* still normalized wrt. bnr:1479* chic(genc[i]) = zeta_C^chic[i]), C = cyc_normalize(bnr.cyc)[1] */1480GEN1481bnrchar_primitive(GEN bnr, GEN nchi, GEN bnrc)1482{ return ag_nchar_image(bnrsurjection(bnr, bnrc), nchi); }14831484/* s: <gen> = Cl_f -> Cl_f2 -> 0, H subgroup of Cl_f (generators given as1485* HNF on [gen]). Return subgroup s(H) in Cl_f2 */1486static GEN1487imageofgroup(GEN bnr, GEN bnr2, GEN H)1488{1489if (!H) return diagonal_shallow(bnr_get_cyc(bnr2));1490return ag_subgroup_image(bnrsurjection(bnr, bnr2), H);1491}1492GEN1493bnrchar_primitive_raw(GEN bnr, GEN bnrc, GEN chi)1494{1495GEN S = bnrsurjection(bnr, bnrc);1496return ag_char_image(S, chi);1497}14981499/* convert A,B,C to [bnr, H] */1500GEN1501ABC_to_bnr(GEN A, GEN B, GEN C, GEN *H, int gen)1502{1503if (typ(A) == t_VEC)1504switch(lg(A))1505{1506case 7: /* bnr */1507*H = B; return A;1508case 11: /* bnf */1509if (!B) pari_err_TYPE("ABC_to_bnr [bnf+missing conductor]",A);1510*H = C; return Buchray(A,B, gen? nf_INIT | nf_GEN: nf_INIT);1511}1512pari_err_TYPE("ABC_to_bnr",A);1513*H = NULL; return NULL; /* LCOV_EXCL_LINE */1514}15151516/* OBSOLETE */1517GEN1518bnrconductor0(GEN A, GEN B, GEN C, long flag)1519{1520pari_sp av = avma;1521GEN H, bnr = ABC_to_bnr(A,B,C,&H, 0);1522return gerepilecopy(av, bnrconductor(bnr, H, flag));1523}15241525long1526bnrisconductor0(GEN A,GEN B,GEN C)1527{1528GEN H, bnr = ABC_to_bnr(A,B,C,&H, 0);1529return bnrisconductor(bnr, H);1530}15311532static GEN1533ideallog_to_bnr_i(GEN Ubid, GEN cyc, GEN z)1534{ return (lg(Ubid)==1)? zerocol(lg(cyc)-1): vecmodii(ZM_ZC_mul(Ubid,z), cyc); }1535/* return bnrisprincipal(bnr, (x)), assuming z = ideallog(x); allow a1536* t_MAT for z, understood as a collection of ideallog(x_i) */1537static GEN1538ideallog_to_bnr(GEN bnr, GEN z)1539{1540GEN U = gel(bnr_get_U(bnr), 2); /* bid part */1541GEN y, cyc = bnr_get_cyc(bnr);1542long i, l;1543if (typ(z) == t_COL) return ideallog_to_bnr_i(U, cyc, z);1544y = cgetg_copy(z, &l);1545for (i = 1; i < l; i++) gel(y,i) = ideallog_to_bnr_i(U, cyc, gel(z,i));1546return y;1547}1548static GEN1549bnr_log_gen_pr(GEN bnr, zlog_S *S, long e, long index)1550{ return ideallog_to_bnr(bnr, log_gen_pr(S, index, bnr_get_nf(bnr), e)); }1551static GEN1552bnr_log_gen_arch(GEN bnr, zlog_S *S, long index)1553{ return ideallog_to_bnr(bnr, log_gen_arch(S, index)); }15541555/* A \subset H ? Allow H = NULL = trivial subgroup */1556static int1557contains(GEN H, GEN A)1558{ return H? (hnf_solve(H, A) != NULL): gequal0(A); }15591560/* finite part of the conductor of H is S.P^e2*/1561static GEN1562cond0_e(GEN bnr, GEN H, zlog_S *S)1563{1564long j, k, l = lg(S->k), iscond0 = S->no2;1565GEN e = S->k, e2 = cgetg(l, t_COL);1566for (k = 1; k < l; k++)1567{1568for (j = itos(gel(e,k)); j > 0; j--)1569{1570if (!contains(H, bnr_log_gen_pr(bnr, S, j, k))) break;1571iscond0 = 0;1572}1573gel(e2,k) = utoi(j);1574}1575return iscond0? NULL: e2;1576}1577/* infinite part of the conductor of H in archp form */1578static GEN1579condoo_archp(GEN bnr, GEN H, zlog_S *S)1580{1581GEN archp = S->archp, archp2 = leafcopy(archp);1582long j, k, l = lg(archp);1583for (k = j = 1; k < l; k++)1584{1585if (!contains(H, bnr_log_gen_arch(bnr, S, k)))1586{1587archp2[j++] = archp[k];1588continue;1589}1590}1591if (j == l) return S->archp;1592setlg(archp2, j); return archp2;1593}1594/* MOD useless in this function */1595static GEN1596bnrconductor_factored_i(GEN bnr, GEN H, long raw)1597{1598GEN nf, bid, ideal, arch, archp, e, fa, cond = NULL;1599zlog_S S;16001601checkbnr(bnr);1602bid = bnr_get_bid(bnr); init_zlog(&S, bid);1603nf = bnr_get_nf(bnr);1604H = bnr_subgroup_check(bnr, H, NULL);1605e = cond0_e(bnr, H, &S); /* in terms of S.P */1606archp = condoo_archp(bnr, H, &S);1607ideal = e? factorbackprime(nf, S.P, e): bid_get_ideal(bid);1608if (archp == S.archp)1609{1610if (!e) cond = bnr_get_mod(bnr);1611arch = bid_get_arch(bid);1612}1613else1614arch = indices_to_vec01(archp, nf_get_r1(nf));1615if (!cond) cond = mkvec2(ideal, arch);1616if (raw) return cond;1617fa = e? famat_remove_trivial(mkmat2(S.P, e)): bid_get_fact(bid);1618return mkvec2(cond, fa);1619}1620GEN1621bnrconductor_factored(GEN bnr, GEN H)1622{ return bnrconductor_factored_i(bnr, H, 0); }1623GEN1624bnrconductor_raw(GEN bnr, GEN H)1625{ return bnrconductor_factored_i(bnr, H, 1); }16261627/* (see bnrdisc_i). Given a bnr, and a subgroup1628* H0 (possibly given as a character chi, in which case H0 = ker chi) of the1629* ray class group, compute the conductor of H if flag=0. If flag > 0, compute1630* also the corresponding H' and output1631* if flag = 1: [[ideal,arch],[hm,cyc,gen],H']1632* if flag = 2: [[ideal,arch],newbnr,H'] */1633GEN1634bnrconductormod(GEN bnr, GEN H0, GEN MOD)1635{1636GEN nf, bid, arch, archp, bnrc, e, H, cond = NULL;1637int ischi;1638zlog_S S;16391640checkbnr(bnr);1641bid = bnr_get_bid(bnr); init_zlog(&S, bid);1642nf = bnr_get_nf(bnr);1643H = bnr_subgroup_check(bnr, H0, NULL);1644e = cond0_e(bnr, H, &S);1645archp = condoo_archp(bnr, H, &S);1646if (archp == S.archp)1647{1648if (!e) cond = bnr_get_mod(bnr);1649arch = gel(bnr_get_mod(bnr), 2);1650}1651else1652arch = indices_to_vec01(archp, nf_get_r1(nf));16531654/* character or subgroup ? */1655ischi = H0 && typ(H0) == t_VEC;1656if (cond)1657{ /* same conductor */1658bnrc = bnr;1659if (ischi)1660H = H0;1661else if (!H)1662H = diagonal_shallow(bnr_get_cyc(bnr));1663}1664else1665{1666long fl = lg(bnr_get_clgp(bnr)) == 4? nf_INIT | nf_GEN: nf_INIT;1667GEN fa = famat_remove_trivial(mkmat2(S.P, e? e: S.k)), bid;1668bid = Idealstarmod(nf, mkvec2(fa, arch), nf_INIT | nf_GEN, MOD);1669bnrc = Buchraymod_i(bnr, bid, fl, MOD);1670cond = bnr_get_mod(bnrc);1671if (ischi)1672H = bnrchar_primitive_raw(bnr, bnrc, H0);1673else1674H = imageofgroup(bnr, bnrc, H);1675}1676return mkvec3(cond, bnrc, H);1677}1678/* OBSOLETE */1679GEN1680bnrconductor_i(GEN bnr, GEN H, long flag)1681{1682GEN v;1683if (flag == 0) return bnrconductor_raw(bnr, H);1684v = bnrconductormod(bnr, H, NULL);1685if (flag == 1) gel(v,2) = bnr_get_clgp(gel(v,2));1686return v;1687}1688/* OBSOLETE */1689GEN1690bnrconductor(GEN bnr, GEN H, long flag)1691{1692pari_sp av = avma;1693if (flag > 2 || flag < 0) pari_err_FLAG("bnrconductor");1694return gerepilecopy(av, bnrconductor_i(bnr, H, flag));1695}16961697long1698bnrisconductor(GEN bnr, GEN H0)1699{1700pari_sp av = avma;1701long j, k, l;1702GEN archp, e, H;1703zlog_S S;17041705checkbnr(bnr);1706init_zlog(&S, bnr_get_bid(bnr));1707if (!S.no2) return 0;1708H = bnr_subgroup_check(bnr, H0, NULL);17091710archp = S.archp;1711e = S.k; l = lg(e);1712for (k = 1; k < l; k++)1713{1714j = itos(gel(e,k));1715if (contains(H, bnr_log_gen_pr(bnr, &S, j, k))) return gc_long(av,0);1716}1717l = lg(archp);1718for (k = 1; k < l; k++)1719if (contains(H, bnr_log_gen_arch(bnr, &S, k))) return gc_long(av,0);1720return gc_long(av,1);1721}17221723/* return the norm group corresponding to the relative extension given by1724* polrel over bnr.bnf, assuming it is abelian and the modulus of bnr is a1725* multiple of the conductor */1726static GEN1727rnfnormgroup_i(GEN bnr, GEN polrel)1728{1729long i, j, degrel, degnf, k;1730GEN bnf, index, discnf, nf, G, detG, fa, gdegrel;1731GEN fac, col, cnd;1732forprime_t S;1733ulong p;17341735checkbnr(bnr); bnf = bnr_get_bnf(bnr);1736nf = bnf_get_nf(bnf);1737cnd = gel(bnr_get_mod(bnr), 1);1738polrel = RgX_nffix("rnfnormgroup", nf_get_pol(nf),polrel,1);1739if (!gequal1(leading_coeff(polrel)))1740pari_err_IMPL("rnfnormgroup for nonmonic polynomials");17411742degrel = degpol(polrel);1743if (umodiu(bnr_get_no(bnr), degrel)) return NULL;1744/* degrel-th powers are in norm group */1745gdegrel = utoipos(degrel);1746G = ZV_snf_gcd(bnr_get_cyc(bnr), gdegrel);1747detG = ZV_prod(G);1748k = abscmpiu(detG,degrel);1749if (k < 0) return NULL;1750if (!k) return diagonal(G);17511752G = diagonal_shallow(G);1753discnf = nf_get_disc(nf);1754index = nf_get_index(nf);1755degnf = nf_get_degree(nf);1756u_forprime_init(&S, 2, ULONG_MAX);1757while ( (p = u_forprime_next(&S)) )1758{1759long oldf, nfa;1760/* If all pr are unramified and have the same residue degree, p =prod pr1761* and including last pr^f or p^f is the same, but the last isprincipal1762* is much easier! oldf is used to track this */17631764if (!umodiu(index, p)) continue; /* can't be treated efficiently */17651766/* primes of degree 1 are enough, and simpler */1767fa = idealprimedec_limit_f(nf, utoipos(p), 1);1768nfa = lg(fa)-1;1769if (!nfa) continue;1770/* all primes above p included ? */1771oldf = (nfa == degnf)? -1: 0;1772for (i=1; i<=nfa; i++)1773{1774GEN pr = gel(fa,i), pp, T, polr, modpr;1775long f, nfac;1776/* if pr (probably) ramified, we have to use all (unramified) P | pr */1777if (idealval(nf,cnd,pr)) { oldf = 0; continue; }1778modpr = zk_to_Fq_init(nf, &pr, &T, &pp); /* T = NULL, pp ignored */1779polr = nfX_to_FqX(polrel, nf, modpr); /* in Fp[X] */1780polr = ZX_to_Flx(polr, p);1781if (!Flx_is_squarefree(polr, p)) { oldf = 0; continue; }17821783fac = gel(Flx_factor(polr, p), 1);1784f = degpol(gel(fac,1));1785if (f == degrel) continue; /* degrel-th powers already included */1786nfac = lg(fac)-1;1787/* check decomposition of pr has Galois type */1788for (j=2; j<=nfac; j++)1789if (degpol(gel(fac,j)) != f) return NULL;1790if (oldf < 0) oldf = f; else if (oldf != f) oldf = 0;17911792/* last prime & all pr^f, pr | p, included. Include p^f instead */1793if (oldf && i == nfa && degrel == nfa*f && !umodiu(discnf, p))1794pr = utoipos(p);17951796/* pr^f = N P, P | pr, hence is in norm group */1797col = bnrisprincipalmod(bnr,pr,gdegrel,0);1798if (f > 1) col = ZC_z_mul(col, f);1799G = ZM_hnf(shallowconcat(G, col));1800detG = ZM_det_triangular(G);1801k = abscmpiu(detG,degrel);1802if (k < 0) return NULL;1803if (!k) { cgiv(detG); return G; }1804}1805}1806return NULL;1807}1808GEN1809rnfnormgroup(GEN bnr, GEN polrel)1810{1811pari_sp av = avma;1812GEN G = rnfnormgroup_i(bnr, polrel);1813if (!G) { set_avma(av); return cgetg(1,t_MAT); }1814return gerepileupto(av, G);1815}18161817GEN1818nf_deg1_prime(GEN nf)1819{1820GEN z, T = nf_get_pol(nf), D = nf_get_disc(nf), f = nf_get_index(nf);1821long degnf = degpol(T);1822forprime_t S;1823pari_sp av;1824ulong p;1825u_forprime_init(&S, degnf, ULONG_MAX);1826av = avma;1827while ( (p = u_forprime_next(&S)) )1828{1829ulong r;1830if (!umodiu(D, p) || !umodiu(f, p)) continue;1831r = Flx_oneroot(ZX_to_Flx(T,p), p);1832if (r != p)1833{1834z = utoi(Fl_neg(r, p));1835z = deg1pol_shallow(gen_1, z, varn(T));1836return idealprimedec_kummer(nf, z, 1, utoipos(p));1837}1838set_avma(av);1839}1840return NULL;1841}18421843static long1844rnfisabelian_i(GEN nf, GEN pol)1845{1846GEN modpr, pr, T, Tnf, pp, ro, nfL, C, a, sig, eq;1847long i, j, l, v;1848ulong p, k, ka;18491850if (typ(nf) == t_POL)1851Tnf = nf;1852else {1853nf = checknf(nf);1854Tnf = nf_get_pol(nf);1855}1856v = varn(Tnf);1857if (degpol(Tnf) != 1 && typ(pol) == t_POL && RgX_is_QX(pol)1858&& rnfisabelian_i(pol_x(v), pol)) return 1;1859pol = RgX_nffix("rnfisabelian",Tnf,pol,1);1860eq = nf_rnfeq(nf,pol); /* init L := K[x]/(pol), nf attached to K */1861C = gel(eq,1); setvarn(C, v); /* L = Q[t]/(C) */1862a = gel(eq,2); setvarn(a, v); /* root of K.pol in L */1863nfL = C;1864ro = nfroots_if_split(&nfL, QXY_QXQ_evalx(pol, a, C));1865if (!ro) return 0;1866l = lg(ro)-1;1867/* small groups are abelian, as are groups of prime order */1868if (l < 6 || uisprime(l)) return 1;18691870pr = nf_deg1_prime(nfL);1871modpr = nf_to_Fq_init(nfL, &pr, &T, &pp);1872p = itou(pp);1873k = umodiu(gel(eq,3), p);1874ka = (k * itou(nf_to_Fq(nfL, a, modpr))) % p;1875sig= cgetg(l+1, t_VECSMALL);1876/* image of c = ro[1] + k a [distinguished root of C] by the l automorphisms1877* sig[i]: ro[1] -> ro[i] */1878for (i = 1; i <= l; i++)1879sig[i] = Fl_add(ka, itou(nf_to_Fq(nfL, gel(ro,i), modpr)), p);1880ro = Q_primpart(ro);1881for (i=2; i<=l; i++) { /* start at 2, since sig[1] = identity */1882gel(ro,i) = ZX_to_Flx(gel(ro,i), p);1883for (j=2; j<i; j++)1884if (Flx_eval(gel(ro,j), sig[i], p)1885!= Flx_eval(gel(ro,i), sig[j], p)) return 0;1886}1887return 1;1888}1889long1890rnfisabelian(GEN nf, GEN pol)1891{ pari_sp av = avma; return gc_long(av, rnfisabelian_i(nf, pol)); }18921893/* Given bnf and T defining an abelian relative extension, compute the1894* corresponding conductor and congruence subgroup. Return1895* [cond,bnr(cond),H] where cond=[ideal,arch] is the conductor. */1896GEN1897rnfconductor0(GEN bnf, GEN T, long flag)1898{1899pari_sp av = avma;1900GEN D, nf, module, bnr, H, lim, Tr, MOD;19011902if (flag < 0 || flag > 2) pari_err_FLAG("rnfconductor");1903bnf = checkbnf(bnf); nf = bnf_get_nf(bnf);1904Tr = rnfdisc_get_T(nf, T, &lim);1905T = nfX_to_monic(nf, Tr, NULL);1906if (!lim)1907D = rnfdisc_factored(nf, T, NULL);1908else1909{1910GEN P, E, Ez;1911long i, l, degT = degpol(T);1912D = idealfactor_partial(nf, nfX_disc(nf, Q_primpart(Tr)), lim);1913P = gel(D,1); l = lg(P);1914E = gel(D,2); Ez = ZV_to_zv(E);1915if (l > 1 && vecsmall_max(Ez) > 1)1916{ /* cheaply update tame primes */1917for (i = 1; i < l; i++)1918{ /* v_pr(f) = 1 + \sum_{0 < i < l} g_i/g_01919<= 1 + max_{i>0} g_i/(g_i-1) \sum_{0 < i < l} g_i -11920<= 1 + (p/(p-1)) * v_P(e(L/K, pr)), P | pr | p */1921GEN pr = gel(P,i), p = pr_get_p(pr), e = gen_1;1922long q, v = z_pvalrem(degT, p, &q);1923if (v)1924{ /* e = e_tame * e_wild, e_wild | p^v */1925long ee, pp = itou(p);1926long t = ugcd(umodiu(subiu(pr_norm(pr),1), q), q); /* e_tame | t */1927/* upper bound for 1 + p/(p-1) * v * e(L/Q,p) */1928ee = 1 + (pp * v * pr_get_e(pr) * upowuu(pp,v) * t) / (pp-1);1929e = utoi(minss(ee, Ez[i]));1930}1931gel(E,i) = e;1932}1933}1934}1935module = mkvec2(D, identity_perm(nf_get_r1(nf)));1936MOD = flag? utoipos(degpol(T)): NULL;1937bnr = Buchraymod_i(bnf, module, nf_INIT|nf_GEN, MOD);1938H = rnfnormgroup_i(bnr,T); if (!H) return gc_const(av,gen_0);1939return gerepilecopy(av, flag == 2? bnrconductor_factored(bnr, H)1940: bnrconductormod(bnr, H, MOD));1941}1942GEN1943rnfconductor(GEN bnf, GEN T) { return rnfconductor0(bnf, T, 0); }19441945static GEN1946prV_norms(GEN v)1947{1948long i, l;1949GEN w = cgetg_copy(v, &l);1950for (i = 1; i < l; i++) gel(w,i) = pr_norm(gel(v,i));1951return w;1952}19531954/* Given a number field bnf=bnr[1], a ray class group structure bnr, and a1955* subgroup H (HNF form) of the ray class group, compute [n, r1, dk]1956* attached to H. If flag & rnf_COND, abort (return NULL) if module is not the1957* conductor. If flag & rnf_REL, return relative data, else absolute */1958static GEN1959bnrdisc_i(GEN bnr, GEN H, long flag)1960{1961const long flcond = flag & rnf_COND;1962GEN nf, clhray, E, ED, dk;1963long k, d, l, n, r1;1964zlog_S S;19651966checkbnr(bnr);1967init_zlog(&S, bnr_get_bid(bnr));1968nf = bnr_get_nf(bnr);1969H = bnr_subgroup_check(bnr, H, &clhray);1970d = itos(clhray);1971if (!H) H = diagonal_shallow(bnr_get_cyc(bnr));1972E = S.k; ED = cgetg_copy(E, &l);1973for (k = 1; k < l; k++)1974{1975long j, e = itos(gel(E,k)), eD = e*d;1976GEN H2 = H;1977for (j = e; j > 0; j--)1978{1979GEN z = bnr_log_gen_pr(bnr, &S, j, k);1980long d2;1981H2 = ZM_hnf(shallowconcat(H2, z));1982d2 = itos( ZM_det_triangular(H2) );1983if (flcond && j==e && d2 == d) return NULL;1984if (d2 == 1) { eD -= j; break; }1985eD -= d2;1986}1987gel(ED,k) = utoi(eD); /* v_{P[k]}(relative discriminant) */1988}1989l = lg(S.archp); r1 = nf_get_r1(nf);1990for (k = 1; k < l; k++)1991{1992if (!contains(H, bnr_log_gen_arch(bnr, &S, k))) { r1--; continue; }1993if (flcond) return NULL;1994}1995/* d = relative degree1996* r1 = number of unramified real places;1997* [P,ED] = factorization of relative discriminant */1998if (flag & rnf_REL)1999{2000n = d;2001dk = factorbackprime(nf, S.P, ED);2002}2003else2004{2005n = d * nf_get_degree(nf);2006r1= d * r1;2007dk = factorback2(prV_norms(S.P), ED);2008if (((n-r1)&3) == 2) dk = negi(dk); /* (2r2) mod 4 = 2: r2(relext) is odd */2009dk = mulii(dk, powiu(absi_shallow(nf_get_disc(nf)), d));2010}2011return mkvec3(utoipos(n), utoi(r1), dk);2012}2013GEN2014bnrdisc(GEN bnr, GEN H, long flag)2015{2016pari_sp av = avma;2017GEN D = bnrdisc_i(bnr, H, flag);2018return D? gerepilecopy(av, D): gc_const(av, gen_0);2019}2020GEN2021bnrdisc0(GEN A, GEN B, GEN C, long flag)2022{2023GEN H, bnr = ABC_to_bnr(A,B,C,&H, 0);2024return bnrdisc(bnr,H,flag);2025}20262027/* Given a number field bnf=bnr[1], a ray class group structure bnr and a2028* vector chi representing a character on the generators bnr[2][3], compute2029* the conductor of chi. */2030GEN2031bnrconductorofchar(GEN bnr, GEN chi)2032{2033pari_sp av = avma;2034return gerepilecopy(av, bnrconductor_raw(bnr, chi));2035}20362037/* \sum U[i]*y[i], U[i],y[i] ZM, we allow lg(y) > lg(U). */2038static GEN2039ZMV_mul(GEN U, GEN y)2040{2041long i, l = lg(U);2042GEN z = NULL;2043if (l == 1) return cgetg(1,t_MAT);2044for (i = 1; i < l; i++)2045{2046GEN u = ZM_mul(gel(U,i), gel(y,i));2047z = z? ZM_add(z, u): u;2048}2049return z;2050}20512052/* t = [bid,U], h = #Cl(K) */2053static GEN2054get_classno(GEN t, GEN h)2055{2056GEN bid = gel(t,1), m = gel(t,2), cyc = bid_get_cyc(bid), U = bid_get_U(bid);2057return mulii(h, ZM_det_triangular(ZM_hnfmodid(ZMV_mul(U,m), cyc)));2058}20592060static void2061chk_listBU(GEN L, const char *s) {2062if (typ(L) != t_VEC) pari_err_TYPE(s,L);2063if (lg(L) > 1) {2064GEN z = gel(L,1);2065if (typ(z) != t_VEC) pari_err_TYPE(s,z);2066if (lg(z) == 1) return;2067z = gel(z,1); /* [bid,U] */2068if (typ(z) != t_VEC || lg(z) != 3) pari_err_TYPE(s,z);2069checkbid(gel(z,1));2070}2071}20722073/* Given lists of [bid, unit ideallogs], return lists of ray class numbers */2074GEN2075bnrclassnolist(GEN bnf,GEN L)2076{2077pari_sp av = avma;2078long i, l = lg(L);2079GEN V, h;20802081chk_listBU(L, "bnrclassnolist");2082if (l == 1) return cgetg(1, t_VEC);2083bnf = checkbnf(bnf);2084h = bnf_get_no(bnf);2085V = cgetg(l,t_VEC);2086for (i = 1; i < l; i++)2087{2088GEN v, z = gel(L,i);2089long j, lz = lg(z);2090gel(V,i) = v = cgetg(lz,t_VEC);2091for (j=1; j<lz; j++) gel(v,j) = get_classno(gel(z,j), h);2092}2093return gerepilecopy(av, V);2094}20952096static GEN2097Lbnrclassno(GEN L, GEN fac)2098{2099long i, l = lg(L);2100for (i=1; i<l; i++)2101if (gequal(gmael(L,i,1),fac)) return gmael(L,i,2);2102pari_err_BUG("Lbnrclassno");2103return NULL; /* LCOV_EXCL_LINE */2104}21052106static GEN2107factordivexact(GEN fa1,GEN fa2)2108{2109long i, j, k, c, l;2110GEN P, E, P1, E1, P2, E2, p1;21112112P1 = gel(fa1,1); E1 = gel(fa1,2); l = lg(P1);2113P2 = gel(fa2,1); E2 = gel(fa2,2);2114P = cgetg(l,t_COL);2115E = cgetg(l,t_COL);2116for (c = i = 1; i < l; i++)2117{2118j = RgV_isin(P2,gel(P1,i));2119if (!j) { gel(P,c) = gel(P1,i); gel(E,c) = gel(E1,i); c++; }2120else2121{2122p1 = subii(gel(E1,i), gel(E2,j)); k = signe(p1);2123if (k < 0) pari_err_BUG("factordivexact [not exact]");2124if (k > 0) { gel(P,c) = gel(P1,i); gel(E,c) = p1; c++; }2125}2126}2127setlg(P, c);2128setlg(E, c); return mkmat2(P, E);2129}2130/* remove index k */2131static GEN2132factorsplice(GEN fa, long k)2133{2134GEN p = gel(fa,1), e = gel(fa,2), P, E;2135long i, l = lg(p) - 1;2136P = cgetg(l, typ(p));2137E = cgetg(l, typ(e));2138for (i=1; i<k; i++) { P[i] = p[i]; E[i] = e[i]; }2139p++; e++;2140for ( ; i<l; i++) { P[i] = p[i]; E[i] = e[i]; }2141return mkvec2(P,E);2142}2143static GEN2144factorpow(GEN fa, long n)2145{2146if (!n) return trivial_fact();2147return mkmat2(gel(fa,1), gmulsg(n, gel(fa,2)));2148}2149static GEN2150factormul(GEN fa1,GEN fa2)2151{2152GEN p, pnew, e, enew, v, P, y = famat_mul_shallow(fa1,fa2);2153long i, c, lx;21542155p = gel(y,1); v = indexsort(p); lx = lg(p);2156e = gel(y,2);2157pnew = vecpermute(p, v);2158enew = vecpermute(e, v);2159P = gen_0; c = 0;2160for (i=1; i<lx; i++)2161{2162if (gequal(gel(pnew,i),P))2163gel(e,c) = addii(gel(e,c),gel(enew,i));2164else2165{2166c++; P = gel(pnew,i);2167gel(p,c) = P;2168gel(e,c) = gel(enew,i);2169}2170}2171setlg(p, c+1);2172setlg(e, c+1); return y;2173}21742175static long2176get_nz(GEN bnf, GEN ideal, GEN arch, long clhray)2177{2178GEN arch2, mod;2179long nz = 0, l = lg(arch), k, clhss;2180if (typ(arch) == t_VECSMALL)2181arch2 = indices_to_vec01(arch,nf_get_r1(bnf_get_nf(bnf)));2182else2183arch2 = leafcopy(arch);2184mod = mkvec2(ideal, arch2);2185for (k = 1; k < l; k++)2186{ /* FIXME: this is wasteful. Use the same algorithm as bnrconductor */2187if (signe(gel(arch2,k)))2188{2189gel(arch2,k) = gen_0; clhss = itos(bnrclassno(bnf,mod));2190gel(arch2,k) = gen_1;2191if (clhss == clhray) return -1;2192}2193else nz++;2194}2195return nz;2196}21972198static GEN2199get_NR1D(long Nf, long clhray, long degk, long nz, GEN fadkabs, GEN idealrel)2200{2201long n, R1;2202GEN dlk;2203if (nz < 0) return mkvec3(gen_0,gen_0,gen_0); /*EMPTY*/2204n = clhray * degk;2205R1 = clhray * nz;2206dlk = factordivexact(factorpow(Z_factor(utoipos(Nf)),clhray), idealrel);2207/* r2 odd, set dlk = -dlk */2208if (((n-R1)&3)==2) dlk = factormul(to_famat_shallow(gen_m1,gen_1), dlk);2209return mkvec3(utoipos(n),2210stoi(R1),2211factormul(dlk,factorpow(fadkabs,clhray)));2212}22132214/* t = [bid,U], h = #Cl(K) */2215static GEN2216get_discdata(GEN t, GEN h)2217{2218GEN bid = gel(t,1), fa = bid_get_fact(bid);2219GEN P = gel(fa,1), E = vec_to_vecsmall(gel(fa,2));2220return mkvec3(mkvec2(P, E), (GEN)itou(get_classno(t, h)), bid_get_mod(bid));2221}2222typedef struct _disc_data {2223long degk;2224GEN bnf, fadk, idealrelinit, V;2225} disc_data;22262227static GEN2228get_discray(disc_data *D, GEN V, GEN z, long N)2229{2230GEN idealrel = D->idealrelinit;2231GEN mod = gel(z,3), Fa = gel(z,1);2232GEN P = gel(Fa,1), E = gel(Fa,2);2233long k, nz, clhray = z[2], lP = lg(P);2234for (k=1; k<lP; k++)2235{2236GEN pr = gel(P,k), p = pr_get_p(pr);2237long e, ep = E[k], f = pr_get_f(pr);2238long S = 0, norm = N, Npr = upowuu(p[2],f), clhss;2239for (e=1; e<=ep; e++)2240{2241GEN fad;2242if (e < ep) { E[k] = ep-e; fad = Fa; }2243else fad = factorsplice(Fa, k);2244norm /= Npr;2245clhss = (long)Lbnrclassno(gel(V,norm), fad);2246if (e==1 && clhss==clhray) { E[k] = ep; return cgetg(1, t_VEC); }2247if (clhss == 1) { S += ep-e+1; break; }2248S += clhss;2249}2250E[k] = ep;2251idealrel = factormul(idealrel, to_famat_shallow(p, utoi(f * S)));2252}2253nz = get_nz(D->bnf, gel(mod,1), gel(mod,2), clhray);2254return get_NR1D(N, clhray, D->degk, nz, D->fadk, idealrel);2255}22562257/* Given a list of bids and attached unit log matrices, return the2258* list of discrayabs. Only keep moduli which are conductors. */2259GEN2260discrayabslist(GEN bnf, GEN L)2261{2262pari_sp av = avma;2263long i, l = lg(L);2264GEN nf, V, D, h;2265disc_data ID;22662267chk_listBU(L, "discrayabslist");2268if (l == 1) return cgetg(1, t_VEC);2269ID.bnf = bnf = checkbnf(bnf);2270nf = bnf_get_nf(bnf);2271h = bnf_get_no(bnf);2272ID.degk = nf_get_degree(nf);2273ID.fadk = absZ_factor(nf_get_disc(nf));2274ID.idealrelinit = trivial_fact();2275V = cgetg(l, t_VEC);2276D = cgetg(l, t_VEC);2277for (i = 1; i < l; i++)2278{2279GEN z = gel(L,i), v, d;2280long j, lz = lg(z);2281gel(V,i) = v = cgetg(lz,t_VEC);2282gel(D,i) = d = cgetg(lz,t_VEC);2283for (j=1; j<lz; j++) {2284gel(d,j) = get_discdata(gel(z,j), h);2285gel(v,j) = get_discray(&ID, D, gel(d,j), i);2286}2287}2288return gerepilecopy(av, V);2289}22902291/* a zsimp is [fa, cyc, v]2292* fa: vecsmall factorisation,2293* cyc: ZV (concatenation of (Z_K/pr^k)^* SNFs), the generators2294* are positive at all real places [defined implicitly by weak approximation]2295* v: ZC (log of units on (Z_K/pr^k)^* components) */2296static GEN2297zsimp(void)2298{2299GEN empty = cgetg(1, t_VECSMALL);2300return mkvec3(mkvec2(empty,empty), cgetg(1,t_VEC), cgetg(1,t_MAT));2301}23022303/* fa a vecsmall factorization, append p^e */2304static GEN2305fasmall_append(GEN fa, long p, long e)2306{2307GEN P = gel(fa,1), E = gel(fa,2);2308retmkvec2(vecsmall_append(P,p), vecsmall_append(E,e));2309}23102311static GEN2312sprk_get_cyc(GEN s) { return gel(s,1); }23132314/* sprk = sprkinit(pr,k), b zsimp with modulus coprime to pr */2315static GEN2316zsimpjoin(GEN b, GEN sprk, GEN U_pr, long prcode, long e)2317{2318GEN fa, cyc = sprk_get_cyc(sprk);2319if (lg(gel(b,2)) == 1) /* trivial group */2320fa = mkvec2(mkvecsmall(prcode),mkvecsmall(e));2321else2322{2323fa = fasmall_append(gel(b,1), prcode, e);2324cyc = shallowconcat(gel(b,2), cyc); /* no SNF ! */2325U_pr = vconcat(gel(b,3),U_pr);2326}2327return mkvec3(fa, cyc, U_pr);2328}2329/* B a zsimp, sgnU = [cyc[f_oo], sgn_{f_oo}(units)] */2330static GEN2331bnrclassno_1(GEN B, ulong h, GEN sgnU)2332{2333long lx = lg(B), j;2334GEN L = cgetg(lx,t_VEC);2335for (j=1; j<lx; j++)2336{2337pari_sp av = avma;2338GEN b = gel(B,j), cyc = gel(b,2), qm = gel(b,3);2339ulong z;2340cyc = shallowconcat(cyc, gel(sgnU,1));2341qm = vconcat(qm, gel(sgnU,2));2342z = itou( mului(h, ZM_det_triangular(ZM_hnfmodid(qm, cyc))) );2343set_avma(av);2344gel(L,j) = mkvec2(gel(b,1), mkvecsmall(z));2345}2346return L;2347}23482349static void2350vecselect_p(GEN A, GEN B, GEN p, long init, long lB)2351{2352long i; setlg(B, lB);2353for (i=init; i<lB; i++) B[i] = A[p[i]];2354}2355/* B := p . A = row selection according to permutation p. Treat only lower2356* right corner init x init */2357static void2358rowselect_p(GEN A, GEN B, GEN p, long init)2359{2360long i, lB = lg(A), lp = lg(p);2361for (i=1; i<init; i++) setlg(B[i],lp);2362for ( ; i<lB; i++) vecselect_p(gel(A,i),gel(B,i),p,init,lp);2363}2364static ulong2365hdet(ulong h, GEN m)2366{2367pari_sp av = avma;2368GEN z = mului(h, ZM_det_triangular(ZM_hnf(m)));2369return gc_ulong(av, itou(z));2370}2371static GEN2372bnrclassno_all(GEN B, ulong h, GEN sgnU)2373{2374long lx, k, kk, j, r1, jj, nba, nbarch;2375GEN _2, L, m, H, mm, rowsel;23762377if (typ(sgnU) == t_VEC) return bnrclassno_1(B,h,sgnU);2378lx = lg(B); if (lx == 1) return B;23792380r1 = nbrows(sgnU); _2 = const_vec(r1, gen_2);2381L = cgetg(lx,t_VEC); nbarch = 1L<<r1;2382for (j=1; j<lx; j++)2383{2384pari_sp av = avma;2385GEN b = gel(B,j), cyc = gel(b,2), qm = gel(b,3);2386long nc = lg(cyc)-1;2387/* [ qm cyc 0 ]2388* [ sgnU 0 2 ] */2389m = ZM_hnfmodid(vconcat(qm, sgnU), shallowconcat(cyc,_2));2390mm = RgM_shallowcopy(m);2391rowsel = cgetg(nc+r1+1,t_VECSMALL);2392H = cgetg(nbarch+1,t_VECSMALL);2393for (k = 0; k < nbarch; k++)2394{2395nba = nc+1;2396for (kk=k,jj=1; jj<=r1; jj++,kk>>=1)2397if (kk&1) rowsel[nba++] = nc + jj;2398setlg(rowsel, nba);2399rowselect_p(m, mm, rowsel, nc+1);2400H[k+1] = hdet(h, mm);2401}2402H = gerepileuptoleaf(av, H);2403gel(L,j) = mkvec2(gel(b,1), H);2404}2405return L;2406}24072408static int2409is_module(GEN v)2410{2411if (lg(v) != 3 || (typ(v) != t_MAT && typ(v) != t_VEC)) return 0;2412return typ(gel(v,1)) == t_VECSMALL && typ(gel(v,2)) == t_VECSMALL;2413}2414GEN2415decodemodule(GEN nf, GEN fa)2416{2417long n, nn, k;2418pari_sp av = avma;2419GEN G, E, id, pr;24202421nf = checknf(nf);2422if (!is_module(fa)) pari_err_TYPE("decodemodule [not a factorization]", fa);2423n = nf_get_degree(nf); nn = n*n; id = NULL;2424G = gel(fa,1);2425E = gel(fa,2);2426for (k=1; k<lg(G); k++)2427{2428long code = G[k], p = code / nn, j = (code%n)+1;2429GEN P = idealprimedec(nf, utoipos(p)), e = stoi(E[k]);2430if (lg(P) <= j) pari_err_BUG("decodemodule [incorrect hash code]");2431pr = gel(P,j);2432id = id? idealmulpowprime(nf,id, pr,e)2433: idealpow(nf, pr,e);2434}2435if (!id) { set_avma(av); return matid(n); }2436return gerepileupto(av,id);2437}24382439/* List of ray class fields. Do all from scratch, bound < 2^30. No subgroups.2440*2441* Output: a vector V, V[k] contains the ideals of norm k. Given such an ideal2442* m, the component is as follows:2443*2444* + if arch = NULL, run through all possible archimedean parts; archs are2445* ordered using inverse lexicographic order, [0,..,0], [1,0,..,0], [0,1,..,0],2446* Component is [m,V] where V is a vector with 2^r1 entries, giving for each2447* arch the triple [N,R1,D], with N, R1, D as in discrayabs; D is in factored2448* form.2449*2450* + otherwise [m,N,R1,D] */2451GEN2452discrayabslistarch(GEN bnf, GEN arch, ulong bound)2453{2454int allarch = (arch==NULL), flbou = 0;2455long degk, j, k, l, nba, nbarch, r1, c, sqbou;2456pari_sp av0 = avma, av, av1;2457GEN nf, p, Z, fa, Disc, U, sgnU, EMPTY, empty, archp;2458GEN res, Ray, discall, idealrel, idealrelinit, fadkabs, BOUND;2459ulong i, h;2460forprime_t S;24612462if (bound == 0)2463pari_err_DOMAIN("discrayabslistarch","bound","==",gen_0,utoi(bound));2464res = discall = NULL; /* -Wall */24652466bnf = checkbnf(bnf);2467nf = bnf_get_nf(bnf);2468r1 = nf_get_r1(nf);2469degk = nf_get_degree(nf);2470fadkabs = absZ_factor(nf_get_disc(nf));2471h = itou(bnf_get_no(bnf));24722473if (allarch)2474{2475if (r1>15) pari_err_IMPL("r1>15 in discrayabslistarch");2476arch = const_vec(r1, gen_1);2477}2478else if (lg(arch)-1 != r1)2479pari_err_TYPE("Idealstar [incorrect archimedean component]",arch);2480U = log_prk_units_init(bnf);2481archp = vec01_to_indices(arch);2482nba = lg(archp)-1;2483sgnU = zm_to_ZM( nfsign_units(bnf, archp, 1) );2484if (!allarch) sgnU = mkvec2(const_vec(nba,gen_2), sgnU);24852486empty = cgetg(1,t_VEC);2487/* what follows was rewritten from Ideallist */2488BOUND = utoipos(bound);2489p = cgetipos(3);2490u_forprime_init(&S, 2, bound);2491av = avma;2492sqbou = (long)sqrt((double)bound) + 1;2493Z = const_vec(bound, empty);2494gel(Z,1) = mkvec(zsimp());2495if (DEBUGLEVEL>1) err_printf("Starting zidealstarunits computations\n");2496/* The goal is to compute Ray (lists of bnrclassno). Z contains "zsimps",2497* simplified bid, from which bnrclassno is easy to compute.2498* Once p > sqbou, delete Z[i] for i > sqbou and compute directly Ray */2499Ray = Z;2500while ((p[2] = u_forprime_next(&S)))2501{2502if (!flbou && p[2] > sqbou)2503{2504flbou = 1;2505if (DEBUGLEVEL>1) err_printf("\nStarting bnrclassno computations\n");2506Z = gerepilecopy(av,Z);2507Ray = cgetg(bound+1, t_VEC);2508for (i=1; i<=bound; i++) gel(Ray,i) = bnrclassno_all(gel(Z,i),h,sgnU);2509Z = vecslice(Z, 1, sqbou);2510}2511fa = idealprimedec_limit_norm(nf,p,BOUND);2512for (j=1; j<lg(fa); j++)2513{2514GEN pr = gel(fa,j);2515long prcode, f = pr_get_f(pr);2516ulong q, Q = upowuu(p[2], f);25172518/* p, f-1, j-1 as a single integer in "base degk" (f,j <= degk)*/2519prcode = (p[2]*degk + f-1)*degk + j-1;2520q = Q;2521/* FIXME: if Q = 2, should start at l = 2 */2522for (l = 1;; l++) /* Q <= bound */2523{2524ulong iQ;2525GEN sprk = log_prk_init(nf, pr, l, NULL);2526GEN U_pr = log_prk_units(nf, U, sprk);2527for (iQ = Q, i = 1; iQ <= bound; iQ += Q, i++)2528{2529GEN pz, p2, p1 = gel(Z,i);2530long lz = lg(p1);2531if (lz == 1) continue;25322533p2 = cgetg(lz,t_VEC); c = 0;2534for (k=1; k<lz; k++)2535{2536GEN z = gel(p1,k), v = gmael(z,1,1); /* primes in zsimp's fact. */2537long lv = lg(v);2538/* If z has a power of pr in its modulus, skip it */2539if (i != 1 && lv > 1 && v[lv-1] == prcode) break;2540gel(p2,++c) = zsimpjoin(z,sprk,U_pr,prcode,l);2541}2542setlg(p2, c+1);2543pz = gel(Ray,iQ);2544if (flbou) p2 = bnrclassno_all(p2,h,sgnU);2545if (lg(pz) > 1) p2 = shallowconcat(pz,p2);2546gel(Ray,iQ) = p2;2547}2548Q = itou_or_0( muluu(Q, q) );2549if (!Q || Q > bound) break;2550}2551}2552if (gc_needed(av,1))2553{2554if(DEBUGMEM>1) pari_warn(warnmem,"[1]: discrayabslistarch");2555gerepileall(av, flbou? 2: 1, &Z, &Ray);2556}2557}2558if (!flbou) /* occurs iff bound = 1,2,4 */2559{2560if (DEBUGLEVEL>1) err_printf("\nStarting bnrclassno computations\n");2561Ray = cgetg(bound+1, t_VEC);2562for (i=1; i<=bound; i++) gel(Ray,i) = bnrclassno_all(gel(Z,i),h,sgnU);2563}2564Ray = gerepilecopy(av, Ray);25652566if (DEBUGLEVEL>1) err_printf("Starting discrayabs computations\n");2567if (allarch) nbarch = 1L<<r1;2568else2569{2570nbarch = 1;2571discall = cgetg(2,t_VEC);2572}2573EMPTY = mkvec3(gen_0,gen_0,gen_0);2574idealrelinit = trivial_fact();2575av1 = avma;2576Disc = const_vec(bound, empty);2577for (i=1; i<=bound; i++)2578{2579GEN sousdisc, sous = gel(Ray,i);2580long ls = lg(sous);2581gel(Disc,i) = sousdisc = cgetg(ls,t_VEC);2582for (j=1; j<ls; j++)2583{2584GEN b = gel(sous,j), clhrayall = gel(b,2), Fa = gel(b,1);2585GEN P = gel(Fa,1), E = gel(Fa,2);2586long lP = lg(P), karch;25872588if (allarch) discall = cgetg(nbarch+1,t_VEC);2589for (karch=0; karch<nbarch; karch++)2590{2591long nz, clhray = clhrayall[karch+1];2592if (allarch)2593{2594long ka, k2;2595nba = 0;2596for (ka=karch,k=1; k<=r1; k++,ka>>=1)2597if (ka & 1) nba++;2598for (k2=1,k=1; k<=r1; k++,k2<<=1)2599if (karch&k2 && clhrayall[karch-k2+1] == clhray)2600{ res = EMPTY; goto STORE; }2601}2602idealrel = idealrelinit;2603for (k=1; k<lP; k++) /* cf get_discray */2604{2605long e, ep = E[k], pf = P[k] / degk, f = (pf%degk) + 1, S = 0;2606ulong normi = i, Npr;2607p = utoipos(pf / degk);2608Npr = upowuu(p[2],f);2609for (e=1; e<=ep; e++)2610{2611long clhss;2612GEN fad;2613if (e < ep) { E[k] = ep-e; fad = Fa; }2614else fad = factorsplice(Fa, k);2615normi /= Npr;2616clhss = Lbnrclassno(gel(Ray,normi),fad)[karch+1];2617if (e==1 && clhss==clhray) { E[k] = ep; res = EMPTY; goto STORE; }2618if (clhss == 1) { S += ep-e+1; break; }2619S += clhss;2620}2621E[k] = ep;2622idealrel = factormul(idealrel, to_famat_shallow(p, utoi(f * S)));2623}2624if (!allarch && nba)2625nz = get_nz(bnf, decodemodule(nf,Fa), arch, clhray);2626else2627nz = r1 - nba;2628res = get_NR1D(i, clhray, degk, nz, fadkabs, idealrel);2629STORE: gel(discall,karch+1) = res;2630}2631res = allarch? mkvec2(Fa, discall)2632: mkvec4(Fa, gel(res,1), gel(res,2), gel(res,3));2633gel(sousdisc,j) = res;2634if (gc_needed(av1,1))2635{2636long jj;2637if(DEBUGMEM>1) pari_warn(warnmem,"[2]: discrayabslistarch");2638for (jj=j+1; jj<ls; jj++) gel(sousdisc,jj) = gen_0; /* dummy */2639Disc = gerepilecopy(av1, Disc);2640sousdisc = gel(Disc,i);2641}2642}2643}2644return gerepilecopy(av0, Disc);2645}26462647int2648subgroup_conductor_ok(GEN H, GEN L)2649{ /* test conductor */2650long i, l = lg(L);2651for (i = 1; i < l; i++)2652if ( hnf_solve(H, gel(L,i)) ) return 0;2653return 1;2654}2655static GEN2656conductor_elts(GEN bnr)2657{2658long le, la, i, k;2659GEN e, L;2660zlog_S S;26612662if (!bnrisconductor(bnr, NULL)) return NULL;2663init_zlog(&S, bnr_get_bid(bnr));2664e = S.k; le = lg(e); la = lg(S.archp);2665L = cgetg(le + la - 1, t_VEC);2666i = 1;2667for (k = 1; k < le; k++)2668gel(L,i++) = bnr_log_gen_pr(bnr, &S, itos(gel(e,k)), k);2669for (k = 1; k < la; k++)2670gel(L,i++) = bnr_log_gen_arch(bnr, &S, k);2671return L;2672}26732674/* Let C a congruence group in bnr, compute its subgroups whose index is2675* described by bound (see subgrouplist) as subgroups of Clk(bnr).2676* Restrict to subgroups having the same conductor as bnr */2677GEN2678subgrouplist_cond_sub(GEN bnr, GEN C, GEN bound)2679{2680pari_sp av = avma;2681long l, i, j;2682GEN D, Mr, U, T, subgrp, L, cyc = bnr_get_cyc(bnr);26832684L = conductor_elts(bnr); if (!L) return cgetg(1,t_VEC);2685Mr = diagonal_shallow(cyc);2686D = ZM_snfall_i(hnf_solve(C, Mr), &U, NULL, 1);2687T = ZM_mul(C, RgM_inv(U));2688subgrp = subgrouplist(D, bound);2689l = lg(subgrp);2690for (i = j = 1; i < l; i++)2691{2692GEN H = ZM_hnfmodid(ZM_mul(T, gel(subgrp,i)), cyc);2693if (subgroup_conductor_ok(H, L)) gel(subgrp, j++) = H;2694}2695setlg(subgrp, j);2696return gerepilecopy(av, subgrp);2697}26982699static GEN2700subgroupcond(GEN bnr, GEN indexbound)2701{2702pari_sp av = avma;2703GEN L = conductor_elts(bnr);27042705if (!L) return cgetg(1, t_VEC);2706L = subgroupcondlist(bnr_get_cyc(bnr), indexbound, L);2707if (indexbound && typ(indexbound) != t_VEC)2708{ /* sort by increasing index if not single value */2709long i, l = lg(L);2710GEN D = cgetg(l,t_VEC);2711for (i=1; i<l; i++) gel(D,i) = ZM_det_triangular(gel(L,i));2712L = vecreverse( vecpermute(L, indexsort(D)) );2713}2714return gerepilecopy(av, L);2715}27162717GEN2718subgrouplist0(GEN cyc, GEN indexbound, long all)2719{2720if (!all && checkbnr_i(cyc)) return subgroupcond(cyc,indexbound);2721if (typ(cyc) != t_VEC || !RgV_is_ZV(cyc)) cyc = member_cyc(cyc);2722return subgrouplist(cyc,indexbound);2723}27242725GEN2726bnrdisclist0(GEN bnf, GEN L, GEN arch)2727{2728if (typ(L)!=t_INT) return discrayabslist(bnf,L);2729return discrayabslistarch(bnf,arch,itos(L));2730}27312732/****************************************************************************/2733/* Galois action on a BNR */2734/****************************************************************************/2735GEN2736bnrautmatrix(GEN bnr, GEN aut)2737{2738pari_sp av = avma;2739GEN bnf = bnr_get_bnf(bnr), nf = bnf_get_nf(bnf), bid = bnr_get_bid(bnr);2740GEN M, Gen = get_Gen(bnf, bid, bnr_get_El(bnr)), cyc = bnr_get_cyc(bnr);2741long i, l = lg(Gen);27422743M = cgetg(l, t_MAT); aut = nfgaloismatrix(nf, aut);2744/* Gen = clg.gen*U, clg.gen = Gen*Ui */2745for (i = 1; i < l; i++)2746gel(M,i) = isprincipalray(bnr, nfgaloismatrixapply(nf, aut, gel(Gen,i)));2747M = ZM_mul(M, bnr_get_Ui(bnr));2748l = lg(M);2749for (i = 1; i < l; i++) gel(M,i) = vecmodii(gel(M,i), cyc);2750return gerepilecopy(av, M);2751}27522753GEN2754bnrgaloismatrix(GEN bnr, GEN aut)2755{2756checkbnr(bnr);2757switch (typ(aut))2758{2759case t_POL:2760case t_COL:2761return bnrautmatrix(bnr, aut);2762case t_VEC:2763{2764pari_sp av = avma;2765long i, l = lg(aut);2766GEN v;2767if (l == 9)2768{2769GEN g = gal_get_gen(aut);2770if (typ(g) == t_VEC) { aut = galoispermtopol(aut, g); l = lg(aut); }2771}2772v = cgetg(l, t_VEC);2773for(i = 1; i < l; i++) gel(v,i) = bnrautmatrix(bnr, gel(aut,i));2774return gerepileupto(av, v);2775}2776default:2777pari_err_TYPE("bnrgaloismatrix", aut);2778return NULL; /*LCOV_EXCL_LINE*/2779}2780}27812782GEN2783bnrgaloisapply(GEN bnr, GEN mat, GEN x)2784{2785pari_sp av=avma;2786GEN cyc;2787checkbnr(bnr);2788cyc = bnr_get_cyc(bnr);2789if (typ(mat)!=t_MAT || !RgM_is_ZM(mat))2790pari_err_TYPE("bnrgaloisapply",mat);2791if (typ(x)!=t_MAT || !RgM_is_ZM(x))2792pari_err_TYPE("bnrgaloisapply",x);2793return gerepileupto(av, ZM_hnfmodid(ZM_mul(mat, x), cyc));2794}27952796static GEN2797check_bnrgal(GEN bnr, GEN M)2798{2799checkbnr(bnr);2800if (typ(M)==t_MAT)2801return mkvec(M);2802else if (typ(M)==t_VEC && lg(M)==9 && typ(gal_get_gen(M))==t_VEC)2803{2804pari_sp av = avma;2805GEN V = galoispermtopol(M, gal_get_gen(M));2806return gerepileupto(av, bnrgaloismatrix(bnr, V));2807}2808else if (!is_vec_t(typ(M)))2809pari_err_TYPE("bnrisgalois",M);2810return M;2811}28122813long2814bnrisgalois(GEN bnr, GEN M, GEN H)2815{2816pari_sp av = avma;2817long i, l;2818if (typ(H)!=t_MAT || !RgM_is_ZM(H))2819pari_err_TYPE("bnrisgalois",H);2820M = check_bnrgal(bnr, M); l = lg(M);2821for (i=1; i<l; i++)2822{2823long res = ZM_equal(bnrgaloisapply(bnr,gel(M,i), H), H);2824if (!res) return gc_long(av,0);2825}2826return gc_long(av,1);2827}282828292830