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 - modules - stark.c (source / functions) Hit Total Coverage
Test: PARI/GP v2.18.1 lcov report (development 30702-bddb8d6928) Lines: 1676 1810 92.6 %
Date: 2026-02-23 02:23:56 Functions: 129 131 98.5 %
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; either version 2 of the License, or (at your option) any later
       8             : version. It is distributed in the hope that it will be useful, but WITHOUT
       9             : ANY WARRANTY WHATSOEVER.
      10             : 
      11             : Check the License for details. You should have received a copy of it, along
      12             : with the package; see the file 'COPYING'. If not, write to the Free Software
      13             : Foundation, Inc., 51 Franklin Street, Fifth Floor, Boston, MA 02110-1301 USA. */
      14             : 
      15             : /*******************************************************************/
      16             : /*                                                                 */
      17             : /*        COMPUTATION OF STARK UNITS OF TOTALLY REAL FIELDS        */
      18             : /*                                                                 */
      19             : /*******************************************************************/
      20             : #include "pari.h"
      21             : #include "paripriv.h"
      22             : 
      23             : #define DEBUGLEVEL DEBUGLEVEL_stark
      24             : 
      25             : /* ComputeCoeff */
      26             : typedef struct {
      27             :   GEN L0, L1, L11, L2; /* VECSMALL of p */
      28             :   GEN L1ray, L11ray; /* precomputed isprincipalray(pr), pr | p */
      29             :   GEN rayZ; /* precomputed isprincipalray(i), i < condZ */
      30             :   long condZ; /* generates cond(bnr) \cap Z, assumed small */
      31             : } LISTray;
      32             : 
      33             : /* Char evaluation */
      34             : typedef struct {
      35             :   long ord;
      36             :   GEN *val, chi;
      37             : } CHI_t;
      38             : 
      39             : /* RecCoeff */
      40             : typedef struct {
      41             :   GEN M, beta, B, U, nB;
      42             :   long v, G, N;
      43             : } RC_data;
      44             : 
      45             : /********************************************************************/
      46             : /*                    Miscellaneous functions                       */
      47             : /********************************************************************/
      48             : static GEN
      49       13563 : chi_get_c(GEN chi) { return gmael(chi,1,2); }
      50             : static long
      51       39111 : chi_get_deg(GEN chi) { return itou(gmael(chi,1,1)); }
      52             : 
      53             : /* Compute the image of logelt by character chi, zeta_ord(chi)^n; return n */
      54             : static ulong
      55       10545 : CharEval_n(GEN chi, GEN logelt)
      56             : {
      57       10545 :   GEN gn = ZV_dotproduct(chi_get_c(chi), logelt);
      58       10545 :   return umodiu(gn, chi_get_deg(chi));
      59             : }
      60             : /* Compute the image of logelt by character chi, as a complex number */
      61             : static GEN
      62       10505 : CharEval(GEN chi, GEN logelt)
      63             : {
      64       10505 :   ulong n = CharEval_n(chi, logelt), d = chi_get_deg(chi);
      65       10505 :   long nn = Fl_center(n,d,d>>1);
      66       10505 :   GEN x = gel(chi,2);
      67       10505 :   x = gpowgs(x, labs(nn));
      68       10505 :   if (nn < 0) x = conj_i(x);
      69       10505 :   return x;
      70             : }
      71             : 
      72             : /* return n such that C(elt) = z^n */
      73             : static ulong
      74      461785 : CHI_eval_n(CHI_t *C, GEN logelt)
      75             : {
      76      461785 :   GEN n = ZV_dotproduct(C->chi, logelt);
      77      461785 :   return umodiu(n, C->ord);
      78             : }
      79             : /* return C(elt) */
      80             : static GEN
      81      460517 : CHI_eval(CHI_t *C, GEN logelt)
      82             : {
      83      460517 :   return C->val[CHI_eval_n(C, logelt)];
      84             : }
      85             : 
      86             : static void
      87        3018 : init_CHI(CHI_t *c, GEN CHI, GEN z)
      88             : {
      89        3018 :   long i, d = chi_get_deg(CHI);
      90        3018 :   GEN *v = (GEN*)new_chunk(d);
      91        3018 :   v[0] = gen_1;
      92        3018 :   if (d != 1)
      93             :   {
      94        3018 :     v[1] = z;
      95       25572 :     for (i=2; i<d; i++) v[i] = gmul(v[i-1], z);
      96             :   }
      97        3018 :   c->chi = chi_get_c(CHI);
      98        3018 :   c->ord = d;
      99        3018 :   c->val = v;
     100        3018 : }
     101             : /* as t_POLMOD */
     102             : static void
     103        1854 : init_CHI_alg(CHI_t *c, GEN CHI) {
     104        1854 :   long d = chi_get_deg(CHI);
     105             :   GEN z;
     106        1854 :   switch(d)
     107             :   {
     108           0 :     case 1: z = gen_1; break;
     109         698 :     case 2: z = gen_m1; break;
     110        1156 :     default: z = mkpolmod(pol_x(0), polcyclo(d,0));
     111             :   }
     112        1854 :   init_CHI(c,CHI, z);
     113        1854 : }
     114             : /* as t_COMPLEX */
     115             : static void
     116        1164 : init_CHI_C(CHI_t *c, GEN CHI) {
     117        1164 :   init_CHI(c,CHI, gel(CHI,2));
     118        1164 : }
     119             : 
     120             : typedef struct {
     121             :   long r; /* rank = lg(gen) */
     122             :   GEN j; /* current elt is gen[1]^j[1] ... gen[r]^j[r] */
     123             :   GEN cyc; /* t_VECSMALL of elementary divisors */
     124             : } GROUP_t;
     125             : 
     126             : static int
     127      578838 : NextElt(GROUP_t *G)
     128             : {
     129      578838 :   long i = 1;
     130      578838 :   if (G->r == 0) return 0; /* no more elt */
     131      640149 :   while (++G->j[i] == G->cyc[i]) /* from 0 to cyc[i]-1 */
     132             :   {
     133       61901 :     G->j[i] = 0;
     134       61901 :     if (++i > G->r) return 0; /* no more elt */
     135             :   }
     136      578248 :   return i; /* we have multiplied by gen[i] */
     137             : }
     138             : 
     139             : /* enumerate all group elements; trivial elt comes last */
     140             : GEN
     141       38176 : cyc2elts(GEN cyc)
     142             : {
     143             :   long i, n;
     144             :   GEN z;
     145             :   GROUP_t G;
     146             : 
     147       38176 :   G.cyc = typ(cyc)==t_VECSMALL? cyc: gtovecsmall(cyc);
     148       38176 :   n = zv_prod(G.cyc);
     149       38176 :   G.r = lg(cyc)-1;
     150       38176 :   G.j = zero_zv(G.r);
     151             : 
     152       38176 :   z = cgetg(n+1, t_VEC);
     153       38176 :   gel(z,n) = leafcopy(G.j); /* trivial elt comes last */
     154      577872 :   for  (i = 1; i < n; i++)
     155             :   {
     156      539696 :     (void)NextElt(&G);
     157      539696 :     gel(z,i) = leafcopy(G.j);
     158             :   }
     159       38176 :   return z;
     160             : }
     161             : 
     162             : /* nchi: a character given by a vector [d, (c_i)], e.g. from char_normalize
     163             :  * such that chi(x) = e((c . log(x)) / d) where log(x) on bnr.gen */
     164             : static GEN
     165        2682 : get_Char(GEN nchi, long prec)
     166        2682 : { return mkvec2(nchi, rootsof1_cx(gel(nchi,1), prec)); }
     167             : 
     168             : /* prime divisors of conductor */
     169             : static GEN
     170         350 : divcond(GEN bnr) {GEN bid = bnr_get_bid(bnr); return gel(bid_get_fact(bid),1);}
     171             : 
     172             : /* vector of prime ideals dividing bnr but not bnrc */
     173             : static GEN
     174         130 : get_prdiff(GEN D, GEN Dc)
     175             : {
     176         130 :   long n, i, l  = lg(D);
     177         130 :   GEN diff = cgetg(l, t_COL);
     178         358 :   for (n = i = 1; i < l; i++)
     179         228 :     if (!tablesearch(Dc, gel(D,i), &cmp_prime_ideal)) gel(diff,n++) = gel(D,i);
     180         130 :   setlg(diff, n); return diff;
     181             : }
     182             : 
     183             : #define ch_prec(x) realprec(gel(x,1))
     184             : #define ch_C(x)    gel(x,1)
     185             : #define ch_bnr(x)  gel(x,2)
     186             : #define ch_3(x)    gel(x,3)
     187             : #define ch_q(x)    gel(x,3)[1]
     188             : #define ch_CHI(x)  gel(x,4)
     189             : #define ch_diff(x) gel(x,5)
     190             : #define ch_CHI0(x) gel(x,6)
     191             : #define ch_small(x) gel(x,7)
     192             : #define ch_comp(x) gel(x,7)[1]
     193             : #define ch_phideg(x) gel(x,7)[2]
     194             : static long
     195         959 : ch_deg(GEN dtcr) { return chi_get_deg(ch_CHI(dtcr)); }
     196             : 
     197             : /********************************************************************/
     198             : /*                    1rst part: find the field K                   */
     199             : /********************************************************************/
     200             : static GEN AllStark(GEN data, long flag, long prec);
     201             : 
     202             : /* Columns of C [HNF] give the generators of a subgroup of the finite abelian
     203             :  * group A [ in terms of implicit generators ], compute data to work in A/C:
     204             :  * 1) order
     205             :  * 2) structure
     206             :  * 3) the matrix A ->> A/C
     207             :  * 4) the subgroup C */
     208             : static GEN
     209         950 : InitQuotient(GEN C)
     210             : {
     211         950 :   GEN U, D = ZM_snfall_i(C, &U, NULL, 1), h = ZV_prod(D);
     212         950 :   return mkvec5(h, D, U, C, cyc_normalize(D));
     213             : }
     214             : 
     215             : /* lift chi character on A/C [Qt from InitQuotient] to character on A [cyc]*/
     216             : static GEN
     217        2633 : LiftChar(GEN Qt, GEN cyc, GEN chi)
     218             : {
     219        2633 :   GEN ncyc = gel(Qt,5), U = gel(Qt,3), nchi = char_normalize(chi, ncyc);
     220        2633 :   GEN c = ZV_ZM_mul(gel(nchi,2), U), d = gel(nchi,1);
     221        2633 :   return char_denormalize(cyc, d, c);
     222             : }
     223             : 
     224             : /* Let C be a subgroup, system of representatives of the quotient */
     225             : static GEN
     226         235 : subgroup_classes(GEN C)
     227             : {
     228         235 :   GEN U, D = ZM_snfall_i(C, &U, NULL, 1), e = cyc2elts(D);
     229         235 :   long i, l = lg(e);
     230             : 
     231         235 :   if (ZM_isidentity(U))
     232        1450 :     for (i = 1; i < l; i++) (void)vecsmall_to_vec_inplace(gel(e,i));
     233             :   else
     234             :   {
     235          10 :     GEN Ui = ZM_inv(U,NULL);
     236          60 :     for (i = 1; i < l; i++) gel(e,i) = ZM_zc_mul(Ui, gel(e,i));
     237             :   }
     238         235 :   return e;
     239             : }
     240             : 
     241             : /* Let s: A -> B given by [P,cycA,cycB] A and B, compute the kernel of s. */
     242             : GEN
     243         325 : abmap_kernel(GEN S)
     244             : {
     245         325 :   GEN U, P = gel(S,1), cycA = gel(S,2), DB = diagonal_shallow(gel(S,3));
     246         325 :   long nA = lg(cycA)-1, rk;
     247             : 
     248         325 :   rk = nA + lg(DB) - lg(ZM_hnfall_i(shallowconcat(P, DB), &U, 1));
     249         325 :   return ZM_hnfmodid(matslice(U, 1,nA, 1,rk), cycA);
     250             : }
     251             : /* let H be a subgroup of A; return s(H) */
     252             : GEN
     253       11934 : abmap_subgroup_image(GEN S, GEN H)
     254       11934 : { return ZM_hnfmodid(ZM_mul(gel(S,1), H),  gel(S,3)); }
     255             : 
     256             : /* Let m and n be two moduli such that n|m and let C be a congruence
     257             :    group modulo n, compute the corresponding congruence group modulo m
     258             :    ie the kernel of the map S: Clk(m) ->> Clk(n)/C */
     259             : static GEN
     260         325 : ComputeKernel(GEN S, GEN dtQ)
     261             : {
     262         325 :   pari_sp av = avma;
     263         325 :   GEN P = ZM_mul(gel(dtQ,3), gel(S,1));
     264         325 :   return gc_upto(av, abmap_kernel(mkvec3(P, gel(S,2), gel(dtQ,2))));
     265             : }
     266             : 
     267             : static long
     268         835 : cyc_is_cyclic(GEN cyc) { return lg(cyc) <= 2 || equali1(gel(cyc,2)); }
     269             : 
     270             : /* Let H be a subgroup of cl(bnr)/sugbroup, return 1 if
     271             :    cl(bnr)/subgoup/H is cyclic and the signature of the
     272             :    corresponding field is equal to sig and no finite prime
     273             :    dividing cond(bnr) is totally split in this field. Return 0
     274             :    otherwise. */
     275             : static long
     276         370 : IsGoodSubgroup(GEN H, GEN bnr, GEN map)
     277             : {
     278         370 :   pari_sp av = avma;
     279             :   GEN S, mod, modH, p1, U, P, PH, bnrH, iH, qH;
     280             :   long j;
     281             : 
     282         370 :   p1 = InitQuotient(H);
     283             :   /* quotient is non cyclic */
     284         370 :   if (!cyc_is_cyclic(gel(p1,2))) return gc_long(av,0);
     285             : 
     286         180 :   (void)ZM_hnfall_i(shallowconcat(map,H), &U, 0);
     287         180 :   setlg(U, lg(H));
     288         660 :   for (j = 1; j < lg(U); j++) setlg(gel(U,j), lg(H));
     289         180 :   p1 = ZM_hnfmodid(U, bnr_get_cyc(bnr)); /* H as a subgroup of bnr */
     290         180 :   modH = bnrconductor_raw(bnr, p1);
     291         180 :   mod  = bnr_get_mod(bnr);
     292             : 
     293             :   /* is the signature correct? */
     294         180 :   if (!gequal(gel(modH,2), gel(mod,2))) return gc_long(av, 0);
     295             : 
     296             :   /* finite part are the same: OK */
     297         130 :   if (gequal(gel(modH,1), gel(mod,1))) return gc_long(av, 1);
     298             : 
     299             :   /* need to check the splitting of primes dividing mod but not modH */
     300          45 :   bnrH = Buchray(bnr, modH, nf_INIT);
     301          45 :   P = divcond(bnr);
     302          45 :   PH = divcond(bnrH);
     303          45 :   S = bnrsurjection(bnr, bnrH);
     304             :   /* H as a subgroup of bnrH */
     305          45 :   iH = abmap_subgroup_image(S, p1);
     306          45 :   qH = InitQuotient(iH);
     307         145 :   for (j = 1; j < lg(P); j++)
     308             :   {
     309         115 :     GEN pr = gel(P, j), e;
     310             :     /* if pr divides modH, it is ramified, so it's good */
     311         115 :     if (tablesearch(PH, pr, cmp_prime_ideal)) continue;
     312             :     /* inertia degree of pr in bnr(modH)/H is charorder(cycH, e) */
     313          40 :     e = ZM_ZC_mul(gel(qH,3), isprincipalray(bnrH, pr));
     314          40 :     e = ZV_ZV_mod(e, gel(qH,2));
     315          40 :     if (ZV_equal0(e)) return gc_long(av,0); /* f = 1 */
     316             :   }
     317          30 :   return gc_long(av,1);
     318             : }
     319             : 
     320             : /* compute list of nontrivial characters trivial on H, modulo complex
     321             :  * conjugation. If flag is set, impose a nontrivial conductor at infinity */
     322             : static GEN
     323         295 : AllChars(GEN bnr, GEN dtQ, long flag)
     324             : {
     325         295 :   GEN v, vchi, cyc = bnr_get_cyc(bnr);
     326         295 :   long n, i, hD = itos(gel(dtQ,1));
     327             :   hashtable *S;
     328             : 
     329         295 :   v = cgetg(hD+1, t_VEC); /* nonconjugate chars */
     330         295 :   vchi = cyc2elts(gel(dtQ,2));
     331         295 :   S = hash_create(hD, (ulong(*)(void*))&hash_GEN,
     332             :                   (int(*)(void*,void*))&ZV_equal, 1);
     333        2928 :   for (i = n = 1; i < hD; i++) /* remove i = hD: trivial char */
     334             :   { /* lift a character of D in Clk(m) */
     335        2633 :     GEN F, lchi = LiftChar(dtQ, cyc, zv_to_ZV(gel(vchi,i))), cchi = NULL;
     336             : 
     337        2633 :     if (hash_search(S, lchi)) continue;
     338        2009 :     F = bnrconductor_raw(bnr, lchi);
     339        2009 :     if (flag && gequal0(gel(F,2))) continue; /* f_oo(chi) trivial ? */
     340             : 
     341         969 :     if (abscmpiu(charorder(cyc,lchi), 2) > 0)
     342             :     { /* nonreal chi: add its conjugate character to S */
     343         624 :       cchi = charconj(cyc, lchi);
     344         624 :       hash_insert(S, cchi, (void*)1);
     345             :     }
     346         969 :     gel(v, n++) = cchi? mkvec3(lchi, F, cchi): mkvec2(lchi, F);
     347             :   }
     348         295 :   setlg(v, n); return v;
     349             : }
     350             : 
     351             : static GEN InitChar(GEN bnr, GEN CR, long flag, long prec);
     352             : static void CharNewPrec(GEN data, long prec);
     353             : 
     354             : /* Given a conductor and a subgroups, return the corresponding complexity and
     355             :  * precision required using quickpol. Fill data[5] with dataCR */
     356             : static long
     357         235 : CplxModulus(GEN data, long *newprec)
     358             : {
     359         235 :   long dprec = DEFAULTPREC;
     360         235 :   pari_sp av = avma;
     361             :   for (;;)
     362           0 :   {
     363         235 :     GEN cpl, pol = AllStark(data, 1, dprec);
     364         235 :     cpl = RgX_fpnorml2(pol, LOWDEFAULTPREC);
     365         235 :     dprec = maxss(dprec, nbits2extraprec(gexpo(pol))) + EXTRAPREC64;
     366         235 :     if (!gequal0(cpl)) { *newprec = dprec; return gexpo(cpl); }
     367           0 :     set_avma(av);
     368           0 :     if (DEBUGLEVEL>1) pari_warn(warnprec, "CplxModulus", dprec);
     369           0 :     CharNewPrec(data, dprec);
     370             :   }
     371             : }
     372             : 
     373             : /* return A \cap B in abelian group defined by cyc. NULL = whole group */
     374             : static GEN
     375         405 : subgp_intersect(GEN cyc, GEN A, GEN B)
     376             : {
     377             :   GEN H, U;
     378             :   long k, lH;
     379         405 :   if (!A) return B;
     380         160 :   if (!B) return A;
     381         160 :   H = ZM_hnfall_i(shallowconcat(A,B), &U, 1);
     382         160 :   setlg(U, lg(A)); lH = lg(H);
     383         570 :   for (k = 1; k < lg(U); k++) setlg(gel(U,k), lH);
     384         160 :   return ZM_hnfmodid(ZM_mul(A,U), cyc);
     385             : }
     386             : 
     387             : static void CharNewPrec(GEN dataCR, long prec);
     388             : /* Let (f,C) be a conductor without infinite part and a congruence group mod f.
     389             :  * Compute (m,D) such that D is a congruence group of conductor m, f | m,
     390             :  * divisible by all the infinite places but one, D is a subgroup of index 2 of
     391             :  * Cm = Ker: Clk(m) -> Clk(f)/C. Consider further the subgroups H of Clk(m)/D
     392             :  * with cyclic quotient Clk(m)/H such that no place dividing m is totally split
     393             :  * in the extension KH corresponding to H: we want their intersection to be
     394             :  * trivial. These H correspond to (the kernels of Galois orbits of) characters
     395             :  * chi of Clk(m)/D such that chi(log_gen_arch(m_oo)) != 1 and for all pr | m
     396             :  * we either have
     397             :  * - chi(log_gen_pr(pr,1)) != 1 [pr | cond(chi) => ramified in KH]
     398             :  * - or [pr \nmid cond(chi)] chi lifted to Clk(m/pr^oo) is not trivial at pr.
     399             :  * We want the map from Clk(m)/D given by the vector of such characters to have
     400             :  * trivial kernel. Return bnr(m), D, Ck(m)/D and Clk(m)/Cm */
     401             : static GEN
     402         230 : FindModulus(GEN bnr, GEN dtQ, long *newprec)
     403             : {
     404         230 :   const long LIMNORM = 400;
     405         230 :   long n, i, maxnorm, minnorm, N, pr, rb, iscyc, olde = LONG_MAX;
     406         230 :   pari_sp av = avma;
     407         230 :   GEN bnf, nf, f, varch, m, rep = NULL;
     408             : 
     409         230 :   bnf = bnr_get_bnf(bnr);
     410         230 :   nf  = bnf_get_nf(bnf);
     411         230 :   N   = nf_get_degree(nf);
     412         230 :   f   = gel(bnr_get_mod(bnr), 1);
     413             : 
     414             :   /* if cpl < rb, it is not necessary to try another modulus */
     415         230 :   rb = expi( powii(mulii(nf_get_disc(nf), ZM_det_triangular(f)),
     416             :                    gmul2n(bnr_get_no(bnr), 3)) );
     417             : 
     418             :   /* Initialization of the possible infinite part */
     419         230 :   varch = cgetg(N+1,t_VEC);
     420         790 :   for (i = 1; i <= N; i++)
     421             :   {
     422         560 :     GEN a = const_vec(N,gen_1);
     423         560 :     gel(a, N+1-i) = gen_0;
     424         560 :     gel(varch, i) = a;
     425             :   }
     426         230 :   m = cgetg(3, t_VEC);
     427             : 
     428             :   /* Go from minnorm up to maxnorm; if necessary, increase these values.
     429             :    * If the extension is cyclic then a suitable conductor exists and we go on
     430             :    * until we find it. Else, stop at norm LIMNORM. */
     431         230 :   minnorm = 1;
     432         230 :   maxnorm = 50;
     433         230 :   iscyc = cyc_is_cyclic(gel(dtQ,2));
     434             : 
     435         230 :   if (DEBUGLEVEL>1)
     436           0 :     err_printf("Looking for a modulus of norm: ");
     437             : 
     438             :   for(;;)
     439           0 :   {
     440         230 :     GEN listid = ideallist0(nf, maxnorm, 4+8); /* ideals of norm <= maxnorm */
     441         230 :     pari_sp av1 = avma;
     442        1045 :     for (n = minnorm; n <= maxnorm; n++, set_avma(av1))
     443             :     {
     444        1045 :       GEN idnormn = gel(listid,n);
     445        1045 :       long nbidnn  = lg(idnormn) - 1;
     446        1045 :       if (DEBUGLEVEL>1) err_printf(" %ld", n);
     447        1640 :       for (i = 1; i <= nbidnn; i++)
     448             :       { /* finite part of the conductor */
     449             :         long s;
     450             : 
     451         825 :         gel(m,1) = idealmul(nf, f, gel(idnormn,i));
     452        2295 :         for (s = 1; s <= N; s++)
     453             :         { /* infinite part */
     454             :           GEN candD, Cm, bnrm;
     455             :           long lD, c;
     456             : 
     457        1700 :           gel(m,2) = gel(varch,s);
     458             :           /* compute Clk(m), check if m is a conductor */
     459        1700 :           bnrm = Buchray(bnf, m, nf_INIT);
     460        1700 :           if (!bnrisconductor(bnrm, NULL)) continue;
     461             : 
     462             :           /* compute Im(C) in Clk(m)... */
     463         320 :           Cm = ComputeKernel(bnrsurjection(bnrm, bnr), dtQ);
     464             :           /* ... and its subgroups of index 2 with conductor m */
     465         320 :           candD = subgrouplist_cond_sub(bnrm, Cm, mkvec(gen_2));
     466         320 :           lD = lg(candD);
     467         325 :           for (c = 1; c < lD; c++)
     468             :           {
     469         235 :             GEN data, CR, D = gel(candD,c), QD = InitQuotient(D);
     470         235 :             GEN ord = gel(QD,1), cyc = gel(QD,2), map = gel(QD,3);
     471             :             long e;
     472             : 
     473         235 :             if (!cyc_is_cyclic(cyc)) /* cyclic => suitable, else test */
     474             :             {
     475          55 :               GEN lH = subgrouplist(cyc, NULL), IK = NULL;
     476          55 :               long j, ok = 0;
     477         410 :               for (j = 1; j < lg(lH); j++)
     478             :               {
     479         405 :                 GEN H = gel(lH, j), IH = subgp_intersect(cyc, IK, H);
     480             :                 /* if H > IK, no need to test H */
     481         405 :                 if (IK && gidentical(IH, IK)) continue;
     482         370 :                 if (IsGoodSubgroup(H, bnrm, map))
     483             :                 {
     484         115 :                   IK = IH; /* intersection of all good subgroups */
     485         115 :                   if (equalii(ord, ZM_det_triangular(IK))) { ok = 1; break; }
     486             :                 }
     487             :               }
     488          55 :               if (!ok) continue;
     489             :             }
     490         230 :             CR = InitChar(bnrm, AllChars(bnrm, QD, 1), 0, DEFAULTPREC);
     491         230 :             data = mkvec4(bnrm, D, subgroup_classes(Cm), CR);
     492         230 :             if (DEBUGLEVEL>1)
     493           0 :               err_printf("\nTrying modulus = %Ps and subgroup = %Ps\n",
     494             :                          bnr_get_mod(bnrm), D);
     495         230 :             e = CplxModulus(data, &pr);
     496         230 :             if (DEBUGLEVEL>1) err_printf("cpl = 2^%ld\n", e);
     497         230 :             if (e < olde)
     498             :             {
     499         230 :               guncloneNULL(rep); rep = gclone(data);
     500         230 :               *newprec = pr; olde = e;
     501             :             }
     502         230 :             if (olde < rb) goto END; /* OK */
     503           0 :             if (DEBUGLEVEL>1) err_printf("Trying to find another modulus...");
     504             :           }
     505             :         }
     506         595 :         if (rep) goto END; /* OK */
     507             :       }
     508             :     }
     509             :     /* if necessary compute more ideals */
     510           0 :     minnorm = maxnorm;
     511           0 :     maxnorm <<= 1;
     512           0 :     if (!iscyc && maxnorm > LIMNORM) return NULL;
     513             :   }
     514         230 : END:
     515         230 :   if (DEBUGLEVEL>1)
     516           0 :     err_printf("No, we're done!\nModulus = %Ps and subgroup = %Ps\n",
     517           0 :                bnr_get_mod(gel(rep,1)), gel(rep,2));
     518         230 :   CharNewPrec(rep, *newprec); return gc_GEN(av, rep);
     519             : }
     520             : 
     521             : /********************************************************************/
     522             : /*                      2nd part: compute W(X)                      */
     523             : /********************************************************************/
     524             : 
     525             : /* find ilambda s.t. Diff*f*ilambda integral and coprime to f
     526             :    and ilambda >> 0 at foo, fa = factorization of f */
     527             : static GEN
     528         590 : get_ilambda(GEN nf, GEN fa, GEN foo)
     529             : {
     530         590 :   GEN x, w, E2, P = gel(fa,1), E = gel(fa,2), D = nf_get_diff(nf);
     531         590 :   long i, l = lg(P);
     532         590 :   if (l == 1) return gen_1;
     533         492 :   w = cgetg(l, t_VEC);
     534         492 :   E2 = cgetg(l, t_COL);
     535        1082 :   for (i = 1; i < l; i++)
     536             :   {
     537         590 :     GEN pr = gel(P,i), t = pr_get_tau(pr);
     538         590 :     long e = itou(gel(E,i)), v = idealval(nf, D, pr);
     539         590 :     if (v) { D = idealdivpowprime(nf, D, pr, utoipos(v)); e += v; }
     540         590 :     gel(E2,i) = stoi(e+1);
     541         590 :     if (typ(t) == t_MAT) t = gel(t,1);
     542         590 :     gel(w,i) = gdiv(nfpow(nf, t, stoi(e)), powiu(pr_get_p(pr),e));
     543             :   }
     544         492 :   x = mkmat2(P, E2);
     545         492 :   return idealchinese(nf, mkvec2(x, foo), w);
     546             : }
     547             : /* compute the list of W(chi) such that Ld(s,chi) = W(chi) Ld(1 - s, chi*),
     548             :  * for all chi in LCHI. All chi have the same conductor (= cond(bnr)). */
     549             : static GEN
     550         744 : ArtinNumber(GEN bnr, GEN LCHI, long prec)
     551             : {
     552         744 :   long ic, i, j, nz, nChar = lg(LCHI)-1;
     553             :   pari_sp av2;
     554             :   GEN sqrtnc, cond, condZ, cond0, cond1, nf, T, cyc, vN, vB, diff, vt, idh;
     555             :   GEN zid, gen, z, nchi, indW, W, classe, s0, s, den, ilambda, sarch;
     556             :   CHI_t **lC;
     557             :   GROUP_t G;
     558             : 
     559         744 :   lC = (CHI_t**)new_chunk(nChar + 1);
     560         744 :   indW = cgetg(nChar + 1, t_VECSMALL);
     561         744 :   W = cgetg(nChar + 1, t_VEC);
     562        2544 :   for (ic = 0, i = 1; i <= nChar; i++)
     563             :   {
     564        1800 :     GEN CHI = gel(LCHI,i);
     565        1800 :     if (chi_get_deg(CHI) <= 2) { gel(W,i) = gen_1; continue; }
     566        1164 :     ic++; indW[ic] = i;
     567        1164 :     lC[ic] = (CHI_t*)new_chunk(sizeof(CHI_t));
     568        1164 :     init_CHI_C(lC[ic], CHI);
     569             :   }
     570         744 :   if (!ic) return W;
     571         590 :   nChar = ic;
     572             : 
     573         590 :   nf    = bnr_get_nf(bnr);
     574         590 :   diff  = nf_get_diff(nf);
     575         590 :   T     = nf_get_Tr(nf);
     576         590 :   cond  = bnr_get_mod(bnr);
     577         590 :   cond0 = gel(cond,1); condZ = gcoeff(cond0,1,1);
     578         590 :   cond1 = gel(cond,2);
     579             : 
     580         590 :   sqrtnc = gsqrt(idealnorm(nf, cond0), prec);
     581         590 :   ilambda = get_ilambda(nf, bid_get_fact(bnr_get_bid(bnr)), cond1);
     582         590 :   idh = idealmul(nf, ilambda, idealmul(nf, diff, cond0)); /* integral */
     583         590 :   ilambda = Q_remove_denom(ilambda, &den);
     584         590 :   z = den? rootsof1_cx(den, prec): NULL;
     585             : 
     586             :   /* compute a system of generators of (Ok/cond)^*, we'll make them
     587             :    * cond1-positive in the main loop */
     588         590 :   zid = Idealstar(nf, cond0, nf_GEN);
     589         590 :   cyc = abgrp_get_cyc(zid);
     590         590 :   gen = abgrp_get_gen(zid);
     591         590 :   nz = lg(gen) - 1;
     592         590 :   sarch = nfarchstar(nf, cond0, vec01_to_indices(cond1));
     593             : 
     594         590 :   nchi = cgetg(nChar+1, t_VEC);
     595        1754 :   for (ic = 1; ic <= nChar; ic++) gel(nchi,ic) = cgetg(nz + 1, t_VECSMALL);
     596        1210 :   for (i = 1; i <= nz; i++)
     597             :   {
     598         620 :     if (is_bigint(gel(cyc,i)))
     599           0 :       pari_err_OVERFLOW("ArtinNumber [conductor too large]");
     600         620 :     gel(gen,i) = set_sign_mod_divisor(nf, NULL, gel(gen,i), sarch);
     601         620 :     classe = isprincipalray(bnr, gel(gen,i));
     602        1888 :     for (ic = 1; ic <= nChar; ic++) {
     603        1268 :       GEN n = gel(nchi,ic);
     604        1268 :       n[i] = CHI_eval_n(lC[ic], classe);
     605             :     }
     606             :   }
     607             : 
     608             :   /* Sum chi(beta) * exp(2i * Pi * Tr(beta * ilambda) where beta
     609             :      runs through the classes of (Ok/cond0)^* and beta cond1-positive */
     610         590 :   vt = gel(T,1); /* ( Tr(w_i) )_i */
     611         590 :   if (typ(ilambda) == t_COL)
     612         416 :     vt = ZV_ZM_mul(vt, zk_multable(nf, ilambda));
     613             :   else
     614         174 :     vt = ZC_Z_mul(vt, ilambda);
     615             :   /*vt = den . (Tr(w_i * ilambda))_i */
     616         590 :   G.cyc = gtovecsmall(cyc);
     617         590 :   G.r = nz;
     618         590 :   G.j = zero_zv(nz);
     619         590 :   vN = zero_Flm_copy(nz, nChar);
     620             : 
     621         590 :   av2 = avma;
     622         590 :   vB = const_vec(nz, gen_1);
     623         590 :   s0 = z? powgi(z, modii(gel(vt,1), den)): gen_1; /* for beta = 1 */
     624         590 :   s = const_vec(nChar, s0);
     625             : 
     626       39142 :   while ( (i = NextElt(&G)) )
     627             :   {
     628       38552 :     GEN b = gel(vB,i);
     629       38552 :     b = nfmuli(nf, b, gel(gen,i));
     630       38552 :     b = typ(b) == t_COL? FpC_red(b, condZ): modii(b, condZ);
     631       77368 :     for (j=1; j<=i; j++) gel(vB,j) = b;
     632             : 
     633      162388 :     for (ic = 1; ic <= nChar; ic++)
     634             :     {
     635      123836 :       GEN v = gel(vN,ic), n = gel(nchi,ic);
     636      123836 :       v[i] = Fl_add(v[i], n[i], lC[ic]->ord);
     637      124260 :       for (j=1; j<i; j++) v[j] = v[i];
     638             :     }
     639             : 
     640       38552 :     gel(vB,i) = b = set_sign_mod_divisor(nf, NULL, b, sarch);
     641       38552 :     if (!z)
     642           0 :       s0 = gen_1;
     643             :     else
     644             :     {
     645       38552 :       b = typ(b) == t_COL? ZV_dotproduct(vt, b): mulii(gel(vt,1),b);
     646       38552 :       s0 = powgi(z, modii(b,den));
     647             :     }
     648      162388 :     for (ic = 1; ic <= nChar; ic++)
     649             :     {
     650      123836 :       GEN v = gel(vN,ic), val = lC[ic]->val[ v[i] ];
     651      123836 :       gel(s,ic) = gadd(gel(s,ic), gmul(val, s0));
     652             :     }
     653             : 
     654       38552 :     if (gc_needed(av2, 1))
     655             :     {
     656          39 :       if (DEBUGMEM > 1) pari_warn(warnmem,"ArtinNumber");
     657          39 :       (void)gc_all(av2, 2, &s, &vB);
     658             :     }
     659             :   }
     660             : 
     661         590 :   classe = isprincipalray(bnr, idh);
     662         590 :   z = powIs(- (lg(gel(sarch,1))-1));
     663             : 
     664        1754 :   for (ic = 1; ic <= nChar; ic++)
     665             :   {
     666        1164 :     s0 = gmul(gel(s,ic), CHI_eval(lC[ic], classe));
     667        1164 :     gel(W, indW[ic]) = gmul(gdiv(s0, sqrtnc), z);
     668             :   }
     669         590 :   return W;
     670             : }
     671             : 
     672             : static GEN
     673         524 : AllArtinNumbers(GEN CR, long prec)
     674             : {
     675         524 :   pari_sp av = avma;
     676         524 :   GEN vChar = gel(CR,1), dataCR = gel(CR,2);
     677         524 :   long j, k, cl = lg(dataCR) - 1, J = lg(vChar)-1;
     678         524 :   GEN W = cgetg(cl+1,t_VEC), WbyCond, LCHI;
     679             : 
     680        1202 :   for (j = 1; j <= J; j++)
     681             :   {
     682         678 :     GEN LChar = gel(vChar,j), ldata = vecpermute(dataCR, LChar);
     683         678 :     GEN dtcr = gel(ldata,1), bnr = ch_bnr(dtcr);
     684         678 :     long l = lg(LChar);
     685             : 
     686         678 :     if (DEBUGLEVEL>1)
     687           0 :       err_printf("* Root Number: cond. no %ld/%ld (%ld chars)\n", j, J, l-1);
     688         678 :     LCHI = cgetg(l, t_VEC);
     689        2412 :     for (k = 1; k < l; k++) gel(LCHI,k) = ch_CHI0(gel(ldata,k));
     690         678 :     WbyCond = ArtinNumber(bnr, LCHI, prec);
     691        2412 :     for (k = 1; k < l; k++) gel(W,LChar[k]) = gel(WbyCond,k);
     692             :   }
     693         524 :   return gc_GEN(av, W);
     694             : }
     695             : 
     696             : /* compute the constant W of the functional equation of
     697             :    Lambda(chi). If flag = 1 then chi is assumed to be primitive */
     698             : GEN
     699          66 : bnrrootnumber(GEN bnr, GEN chi, long flag, long prec)
     700             : {
     701          66 :   pari_sp av = avma;
     702             :   GEN cyc, W;
     703             : 
     704          66 :   if (flag < 0 || flag > 1) pari_err_FLAG("bnrrootnumber");
     705          66 :   checkbnr(bnr);
     706          66 :   if (flag)
     707             :   {
     708           0 :     cyc = bnr_get_cyc(bnr);
     709           0 :     if (!char_check(cyc,chi)) pari_err_TYPE("bnrrootnumber [character]", chi);
     710             :   }
     711             :   else
     712             :   {
     713          66 :     bnr_char_sanitize(&bnr, &chi);
     714          66 :     cyc = bnr_get_cyc(bnr);
     715             :   }
     716          66 :   chi = char_normalize(chi, cyc_normalize(cyc));
     717          66 :   chi = get_Char(chi, prec);
     718          66 :   W = ArtinNumber(bnr, mkvec(chi), prec);
     719          66 :   return gc_GEN(av, gel(W,1));
     720             : }
     721             : 
     722             : /********************************************************************/
     723             : /*               3rd part: initialize the characters                */
     724             : /********************************************************************/
     725             : 
     726             : /* Let chi be a character, A(chi) corresponding to the primes dividing diff
     727             :  * at s = flag. If s = 0, returns [r, A] where r is the order of vanishing
     728             :  * at s = 0 corresponding to diff. */
     729             : static GEN
     730        1520 : AChi(GEN dtcr, long *r, long flag, long prec)
     731             : {
     732        1520 :   GEN A, diff = ch_diff(dtcr), bnrc = ch_bnr(dtcr), chi  = ch_CHI0(dtcr);
     733        1520 :   long i, l = lg(diff);
     734             : 
     735        1520 :   A = gen_1; *r = 0;
     736        1595 :   for (i = 1; i < l; i++)
     737             :   {
     738          75 :     GEN B, pr = gel(diff,i), z = CharEval(chi, isprincipalray(bnrc, pr));
     739          75 :     if (flag)
     740           0 :       B = gsubsg(1, gdiv(z, pr_norm(pr)));
     741          75 :     else if (gequal1(z))
     742             :     {
     743          15 :       B = glog(pr_norm(pr), prec);
     744          15 :       (*r)++;
     745             :     }
     746             :     else
     747          60 :       B = gsubsg(1, z);
     748          75 :     A = gmul(A, B);
     749             :   }
     750        1520 :   return A;
     751             : }
     752             : /* simplified version of Achi: return 1 if L(0,chi) = 0 */
     753             : static int
     754         765 : L_vanishes_at_0(GEN D)
     755             : {
     756         765 :   GEN diff = ch_diff(D), bnrc = ch_bnr(D), chi  = ch_CHI0(D);
     757         765 :   long i, l = lg(diff);
     758         795 :   for (i = 1; i < l; i++)
     759             :   {
     760          40 :     GEN pr = gel(diff,i);
     761          40 :     if (!CharEval_n(chi, isprincipalray(bnrc, pr))) return 1;
     762             :   }
     763         755 :   return 0;
     764             : }
     765             : 
     766             : static GEN
     767         408 : _data3(GEN arch, long r2)
     768             : {
     769         408 :   GEN z = cgetg(4, t_VECSMALL);
     770         408 :   long i, r1 = lg(arch) - 1, q = 0;
     771        1291 :   for (i = 1; i <= r1; i++) if (signe(gel(arch,i))) q++;
     772         408 :   z[1] = q;
     773         408 :   z[2] = r1 - q;
     774         408 :   z[3] = r2; return z;
     775             : }
     776             : static void
     777        1200 : ch_get3(GEN dtcr, long *a, long *b, long *c)
     778        1200 : { GEN v = ch_3(dtcr); *a = v[1]; *b = v[2]; *c = v[3]; }
     779             : /* 2^(2*r2) pi^N */
     780             : static GEN
     781         519 : get_P(GEN nf, long prec)
     782             : {
     783         519 :   long r2 = nf_get_r2(nf), N = nf_get_degree(nf);
     784         519 :   GEN P = powru(mppi(prec), N); shiftr_inplace(P, 2*r2); return P;
     785             : }
     786             : /* sort chars according to conductor */
     787             : static GEN
     788         289 : sortChars(GEN ch)
     789             : {
     790         289 :   long j, l = lg(ch);
     791         289 :   GEN F = cgetg(l, t_VEC);
     792        1258 :   for (j = 1; j < l; j++) gel(F, j) = gmael(ch,j,2);
     793         289 :   return vec_equiv(F);
     794             : }
     795             : 
     796             : /* Given a list [chi, F = cond(chi)] of characters over Cl(bnr), return
     797             :  * [vChar, dataCR], where vChar contains the equivalence classes of
     798             :  * characters with the same conductor, and dataCR contains for each character:
     799             :  * - bnr(F)
     800             :  * - the constant C(F) [t_REAL]
     801             :  * - [q, r1 - q, r2, rc] where
     802             :  *      q = number of real places in F
     803             :  *      rc = max{r1 + r2 - q + 1, r2 + q}
     804             :  * - diff(chi) primes dividing m but not F
     805             :  * - chi in bnr(m)
     806             :  * - chi in bnr(F).
     807             :  * If all is unset, only compute characters s.t. L(chi,0) != 0 */
     808             : static GEN
     809         289 : InitChar(GEN bnr, GEN ch, long all, long prec)
     810             : {
     811         289 :   GEN bnf = checkbnf(bnr), nf = bnf_get_nf(bnf), mod = bnr_get_mod(bnr);
     812         289 :   GEN P, dataCR, ncyc, dk = nf_get_disc(nf), vChar = sortChars(ch);
     813         289 :   long n, l, r2 = nf_get_r2(nf), prec2 = precdbl(prec) + EXTRAPREC64;
     814         289 :   long lv = lg(vChar);
     815             : 
     816         289 :   P = get_P(nf, prec2);
     817         289 :   ncyc = cyc_normalize(bnr_get_cyc(bnr));
     818             : 
     819         289 :   dataCR = cgetg_copy(ch, &l);
     820         697 :   for (n = 1; n < lv; n++)
     821             :   {
     822         408 :     GEN D, bnrc, v = gel(vChar, n); /* indices of chars of given cond */
     823         408 :     long a, i = v[1], lc = lg(v);
     824         408 :     GEN F = gmael(ch,i,2), map = NULL;
     825         408 :     gel(dataCR, i) = D = cgetg(8, t_VEC);
     826         408 :     ch_C(D) = sqrtr_abs(divir(mulii(dk, ZM_det_triangular(gel(F,1))), P));
     827         408 :     ch_3(D) = _data3(gel(F,2), r2);
     828         408 :     if (gequal(F, mod))
     829             :     {
     830         278 :       ch_bnr(D) = bnrc = bnr;
     831         278 :       ch_diff(D) = cgetg(1, t_VEC);
     832             :     }
     833             :     else
     834             :     {
     835         130 :       ch_bnr(D) = bnrc = Buchray(bnf, F, nf_INIT);
     836         130 :       ch_diff(D) = get_prdiff(divcond(bnr), divcond(bnrc));
     837         130 :       map = bnrsurjection(bnr, bnrc);
     838             :     }
     839        1377 :     for (a = 1; a < lc; a++)
     840             :     {
     841         969 :       long i = v[a];
     842         969 :       GEN chi = gmael(ch,i,1);
     843             : 
     844         969 :       if (a > 1) gel(dataCR, i) = D = leafcopy(D);
     845         969 :       chi = char_normalize(chi,ncyc);
     846         969 :       ch_CHI(D) = get_Char(chi, prec2);
     847         969 :       if (!map)
     848         792 :         ch_CHI0(D) = ch_CHI(D);
     849             :       else
     850             :       {
     851         177 :         chi = abmap_nchar_image(map, chi);
     852         177 :         ch_CHI0(D) = get_Char(chi, prec2);
     853             :       }
     854             :       /* set last */
     855         969 :       ch_small(D) = mkvecsmall2(all || !L_vanishes_at_0(D),
     856         969 :                                 eulerphiu(itou(gel(chi,1))));
     857             :     }
     858             :   }
     859         289 :   return mkvec2(vChar, dataCR);
     860             : }
     861             : 
     862             : /* recompute dataCR with the new precision, in place. Don't modify bnr
     863             :  * components in place  */
     864             : static void
     865         230 : CharNewPrec(GEN data, long prec)
     866             : {
     867         230 :   long j, l, prec2 = precdbl(prec) + EXTRAPREC64;
     868         230 :   GEN P, nf, dataCR = gmael(data,4,2), D = gel(dataCR,1);
     869             : 
     870         230 :   if (ch_prec(D) >= prec2) return;
     871         230 :   nf = bnr_get_nf(ch_bnr(D));
     872         230 :   if (nf_get_prec(nf) < prec) nf = nfnewprec_shallow(nf, prec); /* not prec2 */
     873         230 :   P = get_P(nf, prec2); l = lg(dataCR);
     874         965 :   for (j = 1; j < l; j++)
     875             :   {
     876             :     GEN f0, bnr;
     877         735 :     D = gel(dataCR,j); f0 = gel(bnr_get_mod(ch_bnr(D)), 1);
     878         735 :     ch_C(D) = sqrtr_abs(divir(mulii(nf_get_disc(nf),ZM_det_triangular(f0)), P));
     879         735 :     ch_bnr(D) = bnr = shallowcopy(ch_bnr(D));
     880         735 :     gel(bnr,1) = shallowcopy(gel(bnr,1)); gmael(bnr,1,7) = nf;
     881         735 :     ch_CHI(D) = get_Char(gel(ch_CHI(D),1), prec2);
     882         735 :     ch_CHI0(D)= get_Char(gel(ch_CHI0(D),1), prec2);
     883             :   }
     884             : }
     885             : 
     886             : /********************************************************************/
     887             : /*             4th part: compute the coefficients an(chi)           */
     888             : /*                                                                  */
     889             : /* matan entries are arrays of ints containing the coefficients of  */
     890             : /* an(chi) as a polmod modulo polcyclo(order(chi))                     */
     891             : /********************************************************************/
     892             : 
     893             : static void
     894      711717 : _0toCoeff(int *rep, long deg)
     895             : {
     896             :   long i;
     897     3255540 :   for (i=0; i<deg; i++) rep[i] = 0;
     898      711717 : }
     899             : 
     900             : /* transform a polmod into Coeff */
     901             : static void
     902      285514 : Polmod2Coeff(int *rep, GEN polmod, long deg)
     903             : {
     904             :   long i;
     905      285514 :   if (typ(polmod) == t_POLMOD)
     906             :   {
     907      210511 :     GEN pol = gel(polmod,2);
     908      210511 :     long d = degpol(pol);
     909             : 
     910      210511 :     pol += 2;
     911      906833 :     for (i=0; i<=d; i++) rep[i] = itos(gel(pol,i));
     912      492765 :     for (   ; i<deg; i++) rep[i] = 0;
     913             :   }
     914             :   else
     915             :   {
     916       75003 :     rep[0] = itos(polmod);
     917       80865 :     for (i=1; i<deg; i++) rep[i] = 0;
     918             :   }
     919      285514 : }
     920             : 
     921             : /* initialize a deg * n matrix of ints */
     922             : static int**
     923        3051 : InitMatAn(long n, long deg, long flag)
     924             : {
     925             :   long i, j;
     926        3051 :   int *a, **A = (int**)pari_malloc((n+1)*sizeof(int*));
     927        3051 :   A[0] = NULL;
     928     4410085 :   for (i = 1; i <= n; i++)
     929             :   {
     930     4407034 :     a = (int*)pari_malloc(deg*sizeof(int));
     931     4407034 :     A[i] = a; a[0] = (i == 1 || flag);
     932    14284015 :     for (j = 1; j < deg; j++) a[j] = 0;
     933             :   }
     934        3051 :   return A;
     935             : }
     936             : 
     937             : static void
     938        4775 : FreeMat(int **A, long n)
     939             : {
     940             :   long i;
     941     4421214 :   for (i = 0; i <= n; i++)
     942     4416439 :     if (A[i]) pari_free((void*)A[i]);
     943        4775 :   pari_free((void*)A);
     944        4775 : }
     945             : 
     946             : /* initialize Coeff reduction */
     947             : static int**
     948        1724 : InitReduction(long d, long deg)
     949             : {
     950             :   long j;
     951        1724 :   pari_sp av = avma;
     952             :   int **A;
     953             :   GEN polmod, pol;
     954             : 
     955        1724 :   A   = (int**)pari_malloc(deg*sizeof(int*));
     956        1724 :   pol = polcyclo(d, 0);
     957        8078 :   for (j = 0; j < deg; j++)
     958             :   {
     959        6354 :     A[j] = (int*)pari_malloc(deg*sizeof(int));
     960        6354 :     polmod = gmodulo(pol_xn(deg+j, 0), pol);
     961        6354 :     Polmod2Coeff(A[j], polmod, deg);
     962             :   }
     963             : 
     964        1724 :   set_avma(av); return A;
     965             : }
     966             : 
     967             : #if 0
     968             : void
     969             : pan(int **an, long n, long deg)
     970             : {
     971             :   long i,j;
     972             :   for (i = 1; i <= n; i++)
     973             :   {
     974             :     err_printf("n = %ld: ",i);
     975             :     for (j = 0; j < deg; j++) err_printf("%d ",an[i][j]);
     976             :     err_printf("\n");
     977             :   }
     978             : }
     979             : #endif
     980             : 
     981             : /* returns 0 if c is zero, 1 otherwise. */
     982             : static int
     983     5482544 : IsZero(int* c, long deg)
     984             : {
     985             :   long i;
     986    18697316 :   for (i = 0; i < deg; i++)
     987    14620005 :     if (c[i]) return 0;
     988     4077311 :   return 1;
     989             : }
     990             : 
     991             : /* set c0 <-- c0 * c1 */
     992             : static void
     993     1073993 : MulCoeff(int *c0, int* c1, int** reduc, long deg)
     994             : {
     995             :   long i,j;
     996             :   int c, *T;
     997             : 
     998     1073993 :   if (IsZero(c0,deg)) return;
     999             : 
    1000      576051 :   T = (int*)new_chunk(2*deg);
    1001     6077589 :   for (i = 0; i < 2*deg; i++)
    1002             :   {
    1003     5501538 :     c = 0;
    1004    59924229 :     for (j = 0; j <= i; j++)
    1005    54422691 :       if (j < deg && j > i - deg) c += c0[j] * c1[i-j];
    1006     5501538 :     T[i] = c;
    1007             :   }
    1008     3326820 :   for (i = 0; i < deg; i++)
    1009             :   {
    1010     2750769 :     c = T[i];
    1011    28586730 :     for (j = 0; j < deg; j++) c += reduc[j][i] * T[deg+j];
    1012     2750769 :     c0[i] = c;
    1013             :   }
    1014             : }
    1015             : 
    1016             : /* c0 <- c0 + c1 * c2 */
    1017             : static void
    1018     4408551 : AddMulCoeff(int *c0, int *c1, int* c2, int** reduc, long deg)
    1019             : {
    1020             :   long i, j;
    1021             :   pari_sp av;
    1022             :   int c, *t;
    1023             : 
    1024     4408551 :   if (IsZero(c2,deg)) return;
    1025      829182 :   if (!c1) /* c1 == 1 */
    1026             :   {
    1027      509212 :     for (i = 0; i < deg; i++) c0[i] += c2[i];
    1028      209391 :     return;
    1029             :   }
    1030      619791 :   av = avma;
    1031      619791 :   t = (int*)new_chunk(2*deg); /* = c1 * c2, not reduced */
    1032     4799547 :   for (i = 0; i < 2*deg; i++)
    1033             :   {
    1034     4179756 :     c = 0;
    1035    30633450 :     for (j = 0; j <= i; j++)
    1036    26453694 :       if (j < deg && j > i - deg) c += c1[j] * c2[i-j];
    1037     4179756 :     t[i] = c;
    1038             :   }
    1039     2709669 :   for (i = 0; i < deg; i++)
    1040             :   {
    1041     2089878 :     c = t[i];
    1042    14271786 :     for (j = 0; j < deg; j++) c += reduc[j][i] * t[deg+j];
    1043     2089878 :     c0[i] += c;
    1044             :   }
    1045      619791 :   set_avma(av);
    1046             : }
    1047             : 
    1048             : /* evaluate the Coeff. No Garbage collector */
    1049             : static GEN
    1050     2618219 : EvalCoeff(GEN z, int* c, long deg)
    1051             : {
    1052             :   long i,j;
    1053     2618219 :   GEN r, e = NULL;
    1054             : 
    1055             : #if 0
    1056             :   /* standard Horner */
    1057             :   e = stoi(c[deg - 1]);
    1058             :   for (i = deg - 2; i >= 0; i--)
    1059             :     e = gadd(stoi(c[i]), gmul(z, e));
    1060             : #else
    1061             :   /* specific attention to sparse polynomials */
    1062     3760968 :   for (i = deg-1; i >=0; i=j-1)
    1063             :   {
    1064     8858441 :     for (j=i; c[j] == 0; j--)
    1065     7715692 :       if (j==0)
    1066             :       {
    1067     2145160 :         if (!e) return NULL;
    1068      222319 :         if (i!=j) z = gpowgs(z,i-j+1);
    1069      222319 :         return gmul(e,z);
    1070             :       }
    1071     1142749 :     if (e)
    1072             :     {
    1073      447371 :       r = (i==j)? z: gpowgs(z,i-j+1);
    1074      447371 :       e = gadd(gmul(e,r), stoi(c[j]));
    1075             :     }
    1076             :     else
    1077      695378 :       e = stoi(c[j]);
    1078             :   }
    1079             : #endif
    1080      473059 :   return e;
    1081             : }
    1082             : 
    1083             : /* a2 <- copy the n x m array */
    1084             : static void
    1085      257113 : CopyCoeff(int** a, int** a2, long n, long m)
    1086             : {
    1087             :   long i,j;
    1088             : 
    1089     3772379 :   for (i = 1; i <= n; i++)
    1090             :   {
    1091     3515266 :     int *b = a[i], *b2 = a2[i];
    1092    13899040 :     for (j = 0; j < m; j++) b2[j] = b[j];
    1093             :   }
    1094      257113 : }
    1095             : 
    1096             : static void
    1097      257113 : an_AddMul(int **an,int **an2, long np, long n, long deg, GEN chi, int **reduc)
    1098             : {
    1099      257113 :   GEN chi2 = chi;
    1100             :   long q, qk, k;
    1101      257113 :   int *c, *c2 = (int*)new_chunk(deg);
    1102             : 
    1103      257113 :   CopyCoeff(an, an2, n/np, deg);
    1104      257113 :   for (q=np;;)
    1105             :   {
    1106      281709 :     if (gequal1(chi2)) c = NULL; else { Polmod2Coeff(c2, chi2, deg); c = c2; }
    1107     4690260 :     for(k = 1, qk = q; qk <= n; k++, qk += q)
    1108     4408551 :       AddMulCoeff(an[qk], c, an2[k], reduc, deg);
    1109      281709 :     if (! (q = umuluu_le(q,np, n)) ) break;
    1110             : 
    1111       24596 :     chi2 = gmul(chi2, chi);
    1112             :   }
    1113      257113 : }
    1114             : 
    1115             : /* correct the coefficients an(chi) according with diff(chi) in place */
    1116             : static void
    1117        1724 : CorrectCoeff(GEN dtcr, int** an, int** reduc, long n, long deg)
    1118             : {
    1119        1724 :   pari_sp av = avma;
    1120             :   long lg, j;
    1121             :   pari_sp av1;
    1122             :   int **an2;
    1123             :   GEN bnrc, diff;
    1124             :   CHI_t C;
    1125             : 
    1126        1724 :   diff = ch_diff(dtcr); lg = lg(diff) - 1;
    1127        1724 :   if (!lg) return;
    1128             : 
    1129         130 :   if (DEBUGLEVEL>2) err_printf("diff(CHI) = %Ps", diff);
    1130         130 :   bnrc = ch_bnr(dtcr);
    1131         130 :   init_CHI_alg(&C, ch_CHI0(dtcr));
    1132             : 
    1133         130 :   an2 = InitMatAn(n, deg, 0);
    1134         130 :   av1 = avma;
    1135         283 :   for (j = 1; j <= lg; j++)
    1136             :   {
    1137         153 :     GEN pr = gel(diff,j);
    1138         153 :     long Np = upr_norm(pr);
    1139         153 :     GEN chi  = CHI_eval(&C, isprincipalray(bnrc, pr));
    1140         153 :     an_AddMul(an,an2,Np,n,deg,chi,reduc);
    1141         153 :     set_avma(av1);
    1142             :   }
    1143         130 :   FreeMat(an2, n); set_avma(av);
    1144             : }
    1145             : 
    1146             : /* compute the coefficients an in the general case */
    1147             : static int**
    1148        1197 : ComputeCoeff(GEN dtcr, LISTray *R, long n, long deg)
    1149             : {
    1150        1197 :   pari_sp av = avma, av2;
    1151             :   long i, l;
    1152             :   int **an, **reduc, **an2;
    1153             :   GEN L;
    1154             :   CHI_t C;
    1155             : 
    1156        1197 :   init_CHI_alg(&C, ch_CHI(dtcr));
    1157        1197 :   an  = InitMatAn(n, deg, 0);
    1158        1197 :   an2 = InitMatAn(n, deg, 0);
    1159        1197 :   reduc  = InitReduction(C.ord, deg);
    1160        1197 :   av2 = avma;
    1161             : 
    1162        1197 :   L = R->L1; l = lg(L);
    1163      258157 :   for (i=1; i<l; i++, set_avma(av2))
    1164             :   {
    1165      256960 :     long np = L[i];
    1166      256960 :     GEN chi  = CHI_eval(&C, gel(R->L1ray,i));
    1167      256960 :     an_AddMul(an,an2,np,n,deg,chi,reduc);
    1168             :   }
    1169        1197 :   FreeMat(an2, n);
    1170             : 
    1171        1197 :   CorrectCoeff(dtcr, an, reduc, n, deg);
    1172        1197 :   FreeMat(reduc, deg-1);
    1173        1197 :   set_avma(av); return an;
    1174             : }
    1175             : 
    1176             : /********************************************************************/
    1177             : /*              5th part: compute L-functions at s=1                */
    1178             : /********************************************************************/
    1179             : static void
    1180         318 : deg11(LISTray *R, long p, GEN bnr, GEN pr) {
    1181         318 :   vecsmalltrunc_append(R->L1, p);
    1182         318 :   vectrunc_append(R->L1ray, isprincipalray(bnr, pr));
    1183         318 : }
    1184             : static void
    1185       17553 : deg12(LISTray *R, long p, GEN bnr, GEN pr) {
    1186       17553 :   vecsmalltrunc_append(R->L11, p);
    1187       17553 :   vectrunc_append(R->L11ray, isprincipalray(bnr, pr));
    1188       17553 : }
    1189             : static void
    1190          25 : deg0(LISTray *R, long p) { vecsmalltrunc_append(R->L0, p); }
    1191             : static void
    1192       18969 : deg2(LISTray *R, long p) { vecsmalltrunc_append(R->L2, p); }
    1193             : 
    1194             : static void
    1195         147 : InitPrimesQuad(GEN bnr, ulong N0, LISTray *R)
    1196             : {
    1197         147 :   pari_sp av = avma;
    1198         147 :   GEN bnf = bnr_get_bnf(bnr), F = gel(bnr_get_mod(bnr), 1);
    1199         147 :   GEN v, N, nf = bnf_get_nf(bnf), dk = nf_get_disc(nf);
    1200         147 :   long l = 1 + primepi_upper_bound(N0);
    1201         147 :   ulong i, p, FZ = itou(gcoeff(F,1,1)), FZ2 = itou(gcoeff(F,2,2));
    1202             :   forprime_t T;
    1203             : 
    1204         147 :   FZ2 = ugcd(FZ, FZ2); /* content(F) */
    1205         147 :   R->L0 = vecsmalltrunc_init(l);
    1206         147 :   R->L2 = vecsmalltrunc_init(l); R->condZ = FZ;
    1207         147 :   R->L1 = vecsmalltrunc_init(l); R->L1ray = vectrunc_init(l);
    1208         147 :   R->L11= vecsmalltrunc_init(l); R->L11ray= vectrunc_init(l);
    1209         147 :   N = utoipos(2);
    1210         147 :   u_forprime_init(&T, 2, N0);
    1211       37012 :   while ( (p = u_forprime_next(&T)) )
    1212             :   {
    1213       36865 :     N[2] = p;
    1214       36865 :     switch (kroiu(dk, p))
    1215             :     {
    1216       18979 :     case -1: /* inert */
    1217       18979 :       if (FZ % p == 0) deg0(R,p); else deg2(R,p);
    1218       18979 :       break;
    1219       17688 :     case 1: /* split */
    1220       17688 :       if (FZ2 % p == 0) deg0(R,p);
    1221             :       else
    1222             :       {
    1223       17688 :         GEN Lpr = idealprimedec(nf, N);
    1224       17688 :         if (FZ % p) deg12(R, p, bnr, gel(Lpr,1));
    1225             :         else
    1226             :         {
    1227         135 :           long t = ZC_prdvd(gel(F,2), gel(Lpr,1))? 2: 1;
    1228         135 :           deg11(R, p, bnr, gel(Lpr,t));
    1229             :         }
    1230             :       }
    1231       17688 :       break;
    1232         198 :     default: /* ramified */
    1233         198 :       if (FZ % p == 0)
    1234          15 :         deg0(R,p);
    1235             :       else
    1236         183 :         deg11(R, p, bnr, idealprimedec_galois(nf, N));
    1237         198 :       break;
    1238             :     }
    1239             :   }
    1240             :   /* precompute isprincipalray(x), x in Z */
    1241         147 :   v = coprimes_zv(FZ);
    1242         147 :   R->rayZ = cgetg(FZ, t_VEC);
    1243        2812 :   for (i = 1; i < FZ; i++)
    1244             :   {
    1245        2665 :     N[2] = i;
    1246        2665 :     gel(R->rayZ,i) = v[i]? isprincipalray(bnr, N): gen_0;
    1247             :   }
    1248         147 :   (void)gc_all(av, 7, &(R->L0), &(R->L2), &(R->rayZ),
    1249             :               &(R->L1), &(R->L1ray), &(R->L11), &(R->L11ray));
    1250         147 : }
    1251             : 
    1252             : static void
    1253         377 : InitPrimes(GEN bnr, ulong N0, LISTray *R)
    1254             : {
    1255         377 :   GEN bnf = bnr_get_bnf(bnr), cond = gel(bnr_get_mod(bnr), 1);
    1256         377 :   long p, l, condZ, N = lg(cond)-1;
    1257         377 :   GEN DL, prime, BOUND, nf = bnf_get_nf(bnf);
    1258             :   forprime_t T;
    1259             : 
    1260         377 :   R->condZ = condZ = itos(gcoeff(cond,1,1));
    1261         377 :   l = primepi_upper_bound(N0) * N;
    1262         377 :   DL = cgetg(N+1, t_VEC);
    1263         377 :   R->L1 = vecsmalltrunc_init(l);
    1264         377 :   R->L1ray = vectrunc_init(l);
    1265         377 :   u_forprime_init(&T, 2, N0);
    1266         377 :   prime = utoipos(2);
    1267         377 :   BOUND = utoi(N0);
    1268       88295 :   while ( (p = u_forprime_next(&T)) )
    1269             :   {
    1270       87918 :     pari_sp av = avma;
    1271             :     long j, k, lP;
    1272             :     GEN P;
    1273       87918 :     prime[2] = p;
    1274       87918 :     if (DEBUGLEVEL>1 && (p & 2047) == 1) err_printf("%ld ", p);
    1275       87918 :     P = idealprimedec_limit_norm(nf, prime, BOUND); lP = lg(P);
    1276      175431 :     for (j = 1; j < lP; j++)
    1277             :     {
    1278       87513 :       GEN pr  = gel(P,j), dl = NULL;
    1279       87513 :       if (condZ % p || !idealval(nf, cond, pr))
    1280             :       {
    1281       87166 :         dl = gclone( isprincipalray(bnr, pr) );
    1282       87166 :         vecsmalltrunc_append(R->L1, upowuu(p, pr_get_f(pr)));
    1283             :       }
    1284       87513 :       gel(DL,j) = dl;
    1285             :     }
    1286       87918 :     set_avma(av);
    1287      175431 :     for (k = 1; k < j; k++)
    1288             :     {
    1289       87513 :       if (!DL[k]) continue;
    1290       87166 :       vectrunc_append(R->L1ray, ZC_copy(gel(DL,k)));
    1291       87166 :       gunclone(gel(DL,k));
    1292             :     }
    1293             :   }
    1294         377 : }
    1295             : 
    1296             : static GEN /* cf polcoef */
    1297      333624 : _sercoeff(GEN x, long n)
    1298             : {
    1299      333624 :   long i = n - valser(x);
    1300      333624 :   return (i < 0)? gen_0: gel(x,i+2);
    1301             : }
    1302             : 
    1303             : static void
    1304      333624 : affect_coeff(GEN q, long n, GEN y, long t)
    1305             : {
    1306      333624 :   GEN x = _sercoeff(q,-n);
    1307      333624 :   if (x == gen_0) gel(y,n) = NULL;
    1308      165713 :   else { affgr(x, gel(y,n)); shiftr_inplace(gel(y,n), t); }
    1309      333624 : }
    1310             : /* (x-i)(x-(i+1)) */
    1311             : static GEN
    1312       85300 : d2(long i) { return deg2pol_shallow(gen_1, utoineg(2*i+1), muluu(i,i+1), 0); }
    1313             : /* x-i */
    1314             : static GEN
    1315      255948 : d1(long i) { return deg1pol_shallow(gen_1, stoi(-i), 0); }
    1316             : 
    1317             : typedef struct {
    1318             :   GEN c1, aij, bij, cS, cT, powracpi;
    1319             :   long i0, a,b,c, r, rc1, rc2;
    1320             : } ST_t;
    1321             : 
    1322             : /* compute the principal part at the integers s = 0, -1, -2, ..., -i0
    1323             :  * of Gamma((s+1)/2)^a Gamma(s/2)^b Gamma(s)^c / (s - z) with z = 0 and 1 */
    1324             : static void
    1325         241 : ppgamma(ST_t *T, long prec)
    1326             : {
    1327             :   GEN G, G1, G2, A, E, O, x, sqpi, aij, bij;
    1328         241 :   long c = T->c, r = T->r, i0 = T->i0, i, j, s, t, dx;
    1329             :   pari_sp av;
    1330             : 
    1331         241 :   T->aij = aij = cgetg(i0+1, t_VEC);
    1332         241 :   T->bij = bij = cgetg(i0+1, t_VEC);
    1333       85541 :   for (i = 1; i <= i0; i++)
    1334             :   {
    1335             :     GEN p1, p2;
    1336       85300 :     gel(aij,i) = p1 = cgetg(r+1, t_VEC);
    1337       85300 :     gel(bij,i) = p2 = cgetg(r+1, t_VEC);
    1338      252112 :     for (j=1; j<=r; j++) { gel(p1,j) = cgetr(prec); gel(p2,j) = cgetr(prec); }
    1339             :   }
    1340         241 :   av = avma; x = pol_x(0);
    1341         241 :   sqpi = sqrtr_abs(mppi(prec)); /* Gamma(1/2) */
    1342             : 
    1343         241 :   G1 = gexp(integser(psi1series(r-1, 0, prec)), prec); /* Gamma(1 + x) */
    1344         241 :   G = shallowcopy(G1); setvalser(G,-1); /* Gamma(x) */
    1345             : 
    1346             :   /* expansion of log(Gamma(u) / Gamma(1/2)) at u = 1/2 */
    1347         241 :   G2 = cgetg(r+2, t_SER);
    1348         241 :   G2[1] = evalsigne(1) | _evalvalser(1) | evalvarn(0);
    1349         241 :   gel(G2,2) = gneg(gadd(gmul2n(mplog2(prec), 1), mpeuler(prec)));
    1350         506 :   for (i = 1; i < r; i++) gel(G2,i+2) = mulri(gel(G1,i+2), int2um1(i));
    1351         241 :   G2 = gmul(sqpi, gexp(G2, prec)); /* Gamma(1/2 + x) */
    1352             : 
    1353             :  /* We simplify to get one of the following two expressions
    1354             :   * if (b > a) : sqrt(Pi)^a 2^{a-au} Gamma(u)^{a+c} Gamma(  u/2  )^{|b-a|}
    1355             :   * if (b <= a): sqrt(Pi)^b 2^{b-bu} Gamma(u)^{b+c} Gamma((u+1)/2)^{|b-a|} */
    1356         241 :   if (T->b > T->a)
    1357             :   {
    1358          48 :     t = T->a; s = T->b; dx = 1;
    1359          48 :     E = ser_unscale(G, ghalf);
    1360          48 :     O = gmul2n(gdiv(ser_unscale(G2, ghalf), d1(1)), 1); /* Gamma((x-1)/2) */
    1361             :   }
    1362             :   else
    1363             :   {
    1364         193 :     t = T->b; s = T->a; dx = 0;
    1365         193 :     E = ser_unscale(G2, ghalf);
    1366         193 :     O = ser_unscale(G, ghalf);
    1367             :   }
    1368             :   /* (sqrt(Pi) 2^{1-x})^t Gamma(x)^{t+c} */
    1369         241 :   A = gmul(gmul(powru(gmul2n(sqpi,1), t), gpowgs(G, t+c)),
    1370             :            gpow(gen_2, RgX_to_ser(gmulgs(x,-t), r+2), prec));
    1371             :   /* A * Gamma((x - dx + 1)/2)^{s-t} */
    1372         241 :   E = gmul(A, gpowgs(E, s-t));
    1373             :   /* A * Gamma((x - dx)/2)^{s-t} */
    1374         241 :   O = gdiv(gmul(A, gpowgs(O, s-t)), gpowgs(gsubgs(x, 1), t+c));
    1375       42891 :   for (i = 0; i < i0/2; i++)
    1376             :   {
    1377       42650 :     GEN p1, q1, A1 = gel(aij,2*i+1), B1 = gel(bij,2*i+1);
    1378       42650 :     GEN p2, q2, A2 = gel(aij,2*i+2), B2 = gel(bij,2*i+2);
    1379       42650 :     long t1 = i * (s+t), t2 = t1 + t;
    1380             : 
    1381       42650 :     p1 = gdiv(E, d1(2*i));
    1382       42650 :     q1 = gdiv(E, d1(2*i+1));
    1383       42650 :     p2 = gdiv(O, d1(2*i+1));
    1384       42650 :     q2 = gdiv(O, d1(2*i+2));
    1385      126056 :     for (j = 1; j <= r; j++)
    1386             :     {
    1387       83406 :       affect_coeff(p1, j, A1, t1); affect_coeff(q1, j, B1, t1);
    1388       83406 :       affect_coeff(p2, j, A2, t2); affect_coeff(q2, j, B2, t2);
    1389             :     }
    1390       42650 :     E = gdiv(E, gmul(gpowgs(d1(2*i+1+dx), s-t), gpowgs(d2(2*i+1), t+c)));
    1391       42650 :     O = gdiv(O, gmul(gpowgs(d1(2*i+2+dx), s-t), gpowgs(d2(2*i+2), t+c)));
    1392             :   }
    1393         241 :   set_avma(av);
    1394         241 : }
    1395             : 
    1396             : /* chi != 1. Return L(1, chi) if fl & 1, else [r, c] where L(s, chi) ~ c s^r
    1397             :  * at s = 0. */
    1398             : static GEN
    1399         959 : GetValue(GEN dtcr, GEN W, GEN S, GEN T, long fl, long prec)
    1400             : {
    1401             :   GEN cf, z;
    1402             :   long q, b, c, r;
    1403         959 :   int isreal = (ch_deg(dtcr) <= 2);
    1404             : 
    1405         959 :   ch_get3(dtcr, &q, &b, &c);
    1406         959 :   if (fl & 1)
    1407             :   { /* S(chi) + W(chi).T(chi)) / (C(chi) sqrt(Pi)^{r1 - q}) */
    1408         168 :     cf = gmul(ch_C(dtcr), powruhalf(mppi(prec), b));
    1409             : 
    1410         168 :     z = gadd(S, gmul(W, T));
    1411         168 :     if (isreal) z = real_i(z);
    1412         168 :     z = gdiv(z, cf);
    1413         168 :     if (fl & 2) z = gmul(z, AChi(dtcr, &r, 1, prec));
    1414             :   }
    1415             :   else
    1416             :   { /* (W(chi).S(conj(chi)) + T(chi)) / (sqrt(Pi)^q 2^{r1 - q}) */
    1417         791 :     cf = gmul2n(powruhalf(mppi(prec), q), b);
    1418             : 
    1419         791 :     z = gadd(gmul(W, conj_i(S)), conj_i(T));
    1420         791 :     if (isreal) z = real_i(z);
    1421         791 :     z = gdiv(z, cf); r = 0;
    1422         791 :     if (fl & 2) z = gmul(z, AChi(dtcr, &r, 0, prec));
    1423         791 :     z = mkvec2(utoi(b + c + r), z);
    1424             :   }
    1425         959 :   return z;
    1426             : }
    1427             : 
    1428             : /* return the order and the first nonzero term of L(s, chi0)
    1429             :    at s = 0. If flag != 0, adjust the value to get L_S(s, chi0). */
    1430             : static GEN
    1431          30 : GetValue1(GEN bnr, long flag, long prec)
    1432             : {
    1433          30 :   GEN bnf = bnr_get_bnf(bnr), nf = bnf_get_nf(bnf);
    1434          30 :   GEN h = bnf_get_no(bnf), R = bnf_get_reg(bnf);
    1435          30 :   GEN c = gdivgs(mpmul(h, R), -bnf_get_tuN(bnf));
    1436          30 :   long r = lg(nf_get_roots(nf)) - 2; /* r1 + r2 - 1 */;
    1437          30 :   if (flag)
    1438             :   {
    1439           0 :     GEN diff = divcond(bnr);
    1440           0 :     long i, l = lg(diff);
    1441           0 :     r += l - 1;
    1442           0 :     for (i = 1; i < l; i++) c = gmul(c, glog(pr_norm(gel(diff,i)), prec));
    1443             :   }
    1444          30 :   return mkvec2(utoi(r), c);
    1445             : }
    1446             : 
    1447             : /********************************************************************/
    1448             : /*                6th part: recover the coefficients                */
    1449             : /********************************************************************/
    1450             : static long
    1451        1842 : TestOne(GEN plg, RC_data *d)
    1452             : {
    1453        1842 :   long j, v = d->v;
    1454        1842 :   GEN z = gsub(d->beta, gel(plg,v));
    1455        1842 :   if (expo(z) >= d->G) return 0;
    1456        5040 :   for (j = 1; j < lg(plg); j++)
    1457        3512 :     if (j != v && mpcmp(d->B, mpabs_shallow(gel(plg,j))) < 0) return 0;
    1458        1528 :   return 1;
    1459             : }
    1460             : 
    1461             : static GEN
    1462         298 : chk_reccoeff_init(FP_chk_fun *chk, GEN r, GEN mat)
    1463             : {
    1464         298 :   RC_data *d = (RC_data*)chk->data;
    1465         298 :   (void)r; d->U = mat; return d->nB;
    1466             : }
    1467             : 
    1468             : static GEN
    1469         286 : chk_reccoeff(void *data, GEN x)
    1470             : {
    1471         286 :   RC_data *d = (RC_data*)data;
    1472         286 :   GEN v = gmul(d->U, x), z = gel(v,1);
    1473             : 
    1474         286 :   if (!gequal1(z)) return NULL;
    1475         286 :   *++v = evaltyp(t_COL) | _evallg( lg(d->M) ); /* pop 1st elt */
    1476         286 :   if (TestOne(gmul(d->M, v), d)) return v;
    1477           0 :   return NULL;
    1478             : }
    1479             : 
    1480             : /* Using Cohen's method */
    1481             : static GEN
    1482         298 : RecCoeff3(GEN nf, RC_data *d, long prec)
    1483             : {
    1484             :   GEN A, M, nB, cand, p1, B2, C2, tB, beta2, nf2, Bd;
    1485         298 :   GEN beta = d->beta, B = d->B;
    1486         298 :   long N = d->N, v = d->v, e, BIG;
    1487         298 :   long i, j, k, ct = 0, prec2;
    1488         298 :   FP_chk_fun chk = { &chk_reccoeff, &chk_reccoeff_init, NULL, NULL, 0 };
    1489         298 :   chk.data = (void*)d;
    1490             : 
    1491         298 :   d->G = minss(-10, -prec >> 4);
    1492         298 :   BIG = maxss(32, -2*d->G);
    1493         298 :   tB  = sqrtnr(real2n(BIG-N,DEFAULTPREC), N-1);
    1494         298 :   Bd  = grndtoi(gmin_shallow(B, tB), &e);
    1495         298 :   if (e > 0) return NULL; /* failure */
    1496         298 :   Bd = addiu(Bd, 1);
    1497         298 :   prec2 = nbits2prec( expi(Bd) + 192 );
    1498         298 :   prec2 = maxss(precdbl(prec), prec2);
    1499         298 :   B2 = sqri(Bd);
    1500         298 :   C2 = shifti(B2, BIG<<1);
    1501             : 
    1502         298 : LABrcf: ct++;
    1503         298 :   beta2 = gprec_w(beta, prec2);
    1504         298 :   nf2 = nfnewprec_shallow(nf, prec2);
    1505         298 :   d->M = M = nf_get_M(nf2);
    1506             : 
    1507         298 :   A = cgetg(N+2, t_MAT);
    1508        1227 :   for (i = 1; i <= N+1; i++) gel(A,i) = cgetg(N+2, t_COL);
    1509             : 
    1510         298 :   gcoeff(A, 1, 1) = gadd(gmul(C2, gsqr(beta2)), B2);
    1511         929 :   for (j = 2; j <= N+1; j++)
    1512             :   {
    1513         631 :     p1 = gmul(C2, gmul(gneg_i(beta2), gcoeff(M, v, j-1)));
    1514         631 :     gcoeff(A, 1, j) = gcoeff(A, j, 1) = p1;
    1515             :   }
    1516         929 :   for (i = 2; i <= N+1; i++)
    1517        1630 :     for (j = i; j <= N+1; j++)
    1518             :     {
    1519         999 :       p1 = gen_0;
    1520        3207 :       for (k = 1; k <= N; k++)
    1521             :       {
    1522        2208 :         GEN p2 = gmul(gcoeff(M, k, j-1), gcoeff(M, k, i-1));
    1523        2208 :         if (k == v) p2 = gmul(C2, p2);
    1524        2208 :         p1 = gadd(p1,p2);
    1525             :       }
    1526         999 :       gcoeff(A, i, j) = gcoeff(A, j, i) = p1;
    1527             :     }
    1528             : 
    1529         298 :   nB = mului(N+1, B2);
    1530         298 :   d->nB = nB;
    1531         298 :   cand = fincke_pohst(A, nB, -1, prec2, &chk);
    1532             : 
    1533         298 :   if (!cand)
    1534             :   {
    1535           0 :     if (ct > 3) return NULL;
    1536           0 :     prec2 = precdbl(prec2);
    1537           0 :     if (DEBUGLEVEL>1) pari_warn(warnprec,"RecCoeff", prec2);
    1538           0 :     goto LABrcf;
    1539             :   }
    1540             : 
    1541         298 :   cand = gel(cand,1);
    1542         298 :   if (lg(cand) == 2) return gel(cand,1);
    1543             : 
    1544         155 :   if (DEBUGLEVEL>1) err_printf("RecCoeff3: no solution found!\n");
    1545         155 :   return NULL;
    1546             : }
    1547             : 
    1548             : /* Using linear dependence relations */
    1549             : static GEN
    1550        1540 : RecCoeff2(GEN nf,  RC_data *d,  long prec)
    1551             : {
    1552             :   pari_sp av;
    1553        1540 :   GEN vec, M = nf_get_M(nf), beta = d->beta;
    1554        1540 :   long bit, min, max, lM = lg(M);
    1555             : 
    1556        1540 :   d->G = minss(-20, -prec >> 4);
    1557             : 
    1558        1540 :   vec  = vec_prepend(row(M, d->v), gneg(beta));
    1559        1540 :   min = (long)prec * 0.75;
    1560        1540 :   max = (long)prec * 0.98;
    1561        1540 :   av = avma;
    1562        1941 :   for (bit = max; bit >= min; bit-=32, set_avma(av))
    1563             :   {
    1564             :     long e;
    1565        1643 :     GEN v = lindep_bit(vec, bit), z = gel(v,1);
    1566        1643 :     if (!signe(z)) continue;
    1567        1556 :     *++v = evaltyp(t_COL) | _evallg(lM); /* pop 1st elt */
    1568        1556 :     v = grndtoi(gdiv(v, z), &e);
    1569        1556 :     if (e > 0) break;
    1570        1556 :     if (TestOne(RgM_RgC_mul(M, v), d)) return v;
    1571             :   }
    1572             :   /* failure */
    1573         298 :   return RecCoeff3(nf,d,prec);
    1574             : }
    1575             : 
    1576             : /* Attempts to find a polynomial with coefficients in nf such that
    1577             :    its coefficients are close to those of pol at the place v and
    1578             :    less than a certain bound at all the other places. This bound is obtained by assuming
    1579             :    that the roots at the other places are, in absolute values, less than 2 if flag == 0,
    1580             :    and less than 1 if flag = 1. */
    1581             : static GEN
    1582         390 : RecCoeff(GEN nf,  GEN pol,  long v, long flag, long prec)
    1583             : {
    1584         390 :   long j, md, cl = degpol(pol);
    1585         390 :   pari_sp av = avma;
    1586             :   RC_data d;
    1587             : 
    1588             :   /* if precision(pol) is too low, abort */
    1589        2645 :   for (j = 2; j <= cl+1; j++)
    1590             :   {
    1591        2255 :     GEN t = gel(pol, j);
    1592        2255 :     if (precision(t) - gexpo(t) < 34) return NULL;
    1593             :   }
    1594             : 
    1595         390 :   md = cl/2;
    1596         390 :   pol = leafcopy(pol);
    1597             : 
    1598         390 :   d.N = nf_get_degree(nf);
    1599         390 :   d.v = v;
    1600             : 
    1601        1775 :   for (j = 1; j <= cl; j++)
    1602             :   { /* start with the coefficients in the middle,
    1603             :        since they are the harder to recognize! */
    1604        1540 :     long cf = md + (j%2? j/2: -j/2);
    1605        1540 :     GEN t, bound = binomial(utoipos(cl), cf);
    1606             : 
    1607        1540 :     if (flag == 0) bound = shifti(bound, cl-cf);
    1608        1540 :     if (DEBUGLEVEL>1) err_printf("RecCoeff (cf = %ld, B = %Ps)\n", cf, bound);
    1609        1540 :     d.beta = real_i( gel(pol,cf+2) );
    1610        1540 :     d.B    = bound;
    1611        1540 :     if (! (t = RecCoeff2(nf, &d, prec)) ) return NULL;
    1612        1385 :     gel(pol, cf+2) = nf_to_scalar_or_polmod(nf, t);
    1613             :   }
    1614         235 :   gel(pol,cl+2) = gen_1;
    1615         235 :   return gc_GEN(av, pol);
    1616             : }
    1617             : 
    1618             : /* an[q * i] *= chi for all (i,p)=1 */
    1619             : static void
    1620       79231 : an_mul(int **an, long p, long q, long n, long deg, GEN chi, int **reduc)
    1621             : {
    1622             :   pari_sp av;
    1623             :   long c,i;
    1624             :   int *T;
    1625             : 
    1626       79231 :   if (gequal1(chi)) return;
    1627       72981 :   av = avma;
    1628       72981 :   T = (int*)new_chunk(deg); Polmod2Coeff(T,chi, deg);
    1629     1533859 :   for (c = 1, i = q; i <= n; i += q, c++)
    1630     1460878 :     if (c == p) c = 0; else MulCoeff(an[i], T, reduc, deg);
    1631       72981 :   set_avma(av);
    1632             : }
    1633             : /* an[q * i] = 0 for all (i,p)=1 */
    1634             : static void
    1635       71708 : an_set0_coprime(int **an, long p, long q, long n, long deg)
    1636             : {
    1637             :   long c,i;
    1638      863840 :   for (c = 1, i = q; i <= n; i += q, c++)
    1639      792132 :     if (c == p) c = 0; else _0toCoeff(an[i], deg);
    1640       71708 : }
    1641             : /* an[q * i] = 0 for all i */
    1642             : static void
    1643          75 : an_set0(int **an, long p, long n, long deg)
    1644             : {
    1645             :   long i;
    1646       42710 :   for (i = p; i <= n; i += p) _0toCoeff(an[i], deg);
    1647          75 : }
    1648             : 
    1649             : /* compute the coefficients an for the quadratic case */
    1650             : static int**
    1651         527 : computean(GEN dtcr, LISTray *R, long n, long deg)
    1652             : {
    1653         527 :   pari_sp av = avma, av2;
    1654             :   long i, p, q, condZ, l;
    1655             :   int **an, **reduc;
    1656             :   GEN L, chi, chi1;
    1657             :   CHI_t C;
    1658             : 
    1659         527 :   init_CHI_alg(&C, ch_CHI(dtcr));
    1660         527 :   condZ= R->condZ;
    1661             : 
    1662         527 :   an = InitMatAn(n, deg, 1);
    1663         527 :   reduc = InitReduction(C.ord, deg);
    1664         527 :   av2 = avma;
    1665             : 
    1666             :   /* all pr | p divide cond */
    1667         527 :   L = R->L0; l = lg(L);
    1668         602 :   for (i=1; i<l; i++) an_set0(an,L[i],n,deg);
    1669             : 
    1670             :   /* 1 prime of degree 2 */
    1671         527 :   L = R->L2; l = lg(L);
    1672       71196 :   for (i=1; i<l; i++, set_avma(av2))
    1673             :   {
    1674       70669 :     p = L[i];
    1675       70669 :     if (condZ == 1) chi = C.val[0]; /* 1 */
    1676       70537 :     else            chi = CHI_eval(&C, gel(R->rayZ, p%condZ));
    1677       70669 :     chi1 = chi;
    1678       70669 :     for (q=p;;)
    1679             :     {
    1680       71708 :       an_set0_coprime(an, p,q,n,deg); /* v_p(q) odd */
    1681       71708 :       if (! (q = umuluu_le(q,p, n)) ) break;
    1682             : 
    1683        3289 :       an_mul(an,p,q,n,deg,chi,reduc);
    1684        3289 :       if (! (q = umuluu_le(q,p, n)) ) break;
    1685        1039 :       chi = gmul(chi, chi1);
    1686             :     }
    1687             :   }
    1688             : 
    1689             :   /* 1 prime of degree 1 */
    1690         527 :   L = R->L1; l = lg(L);
    1691        1730 :   for (i=1; i<l; i++, set_avma(av2))
    1692             :   {
    1693        1203 :     p = L[i];
    1694        1203 :     chi = CHI_eval(&C, gel(R->L1ray,i));
    1695        1203 :     chi1 = chi;
    1696        1203 :     for(q=p;;)
    1697             :     {
    1698        5550 :       an_mul(an,p,q,n,deg,chi,reduc);
    1699        5550 :       if (! (q = umuluu_le(q,p, n)) ) break;
    1700        4347 :       chi = gmul(chi, chi1);
    1701             :     }
    1702             :   }
    1703             : 
    1704             :   /* 2 primes of degree 1 */
    1705         527 :   L = R->L11; l = lg(L);
    1706       65777 :   for (i=1; i<l; i++, set_avma(av2))
    1707             :   {
    1708             :     GEN ray1, ray2, chi11, chi12, chi2;
    1709             : 
    1710       65250 :     p = L[i]; ray1 = gel(R->L11ray,i); /* use pr1 pr2 = (p) */
    1711       65250 :     if (condZ == 1)
    1712          96 :       ray2 = ZC_neg(ray1);
    1713             :     else
    1714       65154 :       ray2 = ZC_sub(gel(R->rayZ, p%condZ),  ray1);
    1715       65250 :     chi11 = CHI_eval(&C, ray1);
    1716       65250 :     chi12 = CHI_eval(&C, ray2);
    1717             : 
    1718       65250 :     chi1 = gadd(chi11, chi12);
    1719       65250 :     chi2 = chi12;
    1720       65250 :     for(q=p;;)
    1721             :     {
    1722       70392 :       an_mul(an,p,q,n,deg,chi1,reduc);
    1723       70392 :       if (! (q = umuluu_le(q,p, n)) ) break;
    1724        5142 :       chi2 = gmul(chi2, chi12);
    1725        5142 :       chi1 = gadd(chi2, gmul(chi1, chi11));
    1726             :     }
    1727             :   }
    1728             : 
    1729         527 :   CorrectCoeff(dtcr, an, reduc, n, deg);
    1730         527 :   FreeMat(reduc, deg-1);
    1731         527 :   set_avma(av); return an;
    1732             : }
    1733             : 
    1734             : /* return the vector of A^i/i for i = 1...n */
    1735             : static GEN
    1736         167 : mpvecpowdiv(GEN A, long n)
    1737             : {
    1738         167 :   pari_sp av = avma;
    1739             :   long i;
    1740         167 :   GEN v = powersr(A, n);
    1741         167 :   GEN w = cgetg(n+1, t_VEC);
    1742         167 :   gel(w,1) = rcopy(gel(v,2));
    1743      261998 :   for (i=2; i<=n; i++) gel(w,i) = divru(gel(v,i+1), i);
    1744         167 :   return gc_upto(av, w);
    1745             : }
    1746             : 
    1747             : static void GetST0(GEN bnr, GEN *pS, GEN *pT, GEN CR, long prec);
    1748             : /* allocate memory for GetST answer */
    1749             : static void
    1750         319 : ST_alloc(GEN *pS, GEN *pT, long l, long prec)
    1751             : {
    1752             :   long j;
    1753         319 :   *pS = cgetg(l, t_VEC);
    1754         319 :   *pT = cgetg(l, t_VEC);
    1755        1456 :   for (j = 1; j < l; j++)
    1756             :   {
    1757        1137 :     gel(*pS,j) = cgetc(prec);
    1758        1137 :     gel(*pT,j) = cgetc(prec);
    1759             :   }
    1760         319 : }
    1761             : 
    1762             : /* compute S and T for the quadratic case. The following cases are:
    1763             :  * 1) bnr complex;
    1764             :  * 2) bnr real and no infinite place divide cond_chi (TODO);
    1765             :  * 3) bnr real and one infinite place divide cond_chi;
    1766             :  * 4) bnr real and both infinite places divide cond_chi (TODO) */
    1767             : static void
    1768         177 : QuadGetST(GEN bnr, GEN *pS, GEN *pT, GEN CR, long prec)
    1769             : {
    1770             :   pari_sp av, av1, av2;
    1771             :   long ncond, n, j, k, n0;
    1772         177 :   GEN vChar = gel(CR,1), dataCR = gel(CR,2), S, T, an, cs, N0, C;
    1773             :   LISTray LIST;
    1774             : 
    1775             :   /* initializations */
    1776         177 :   ST_alloc(pS, pT, lg(dataCR), prec); T = *pT; S = *pS;
    1777         177 :   av = avma;
    1778         177 :   ncond = lg(vChar)-1;
    1779         177 :   C    = cgetg(ncond+1, t_VEC);
    1780         177 :   N0   = cgetg(ncond+1, t_VECSMALL);
    1781         177 :   cs   = cgetg(ncond+1, t_VECSMALL);
    1782         177 :   n0 = 0;
    1783         344 :   for (j = 1; j <= ncond; j++)
    1784             :   {
    1785         197 :     GEN dtcr = gel(dataCR, mael(vChar,j,1)), c = ch_C(dtcr);
    1786             :     long r1, r2;
    1787             : 
    1788         197 :     gel(C,j) = c;
    1789         197 :     nf_get_sign(bnr_get_nf(ch_bnr(dtcr)), &r1, &r2);
    1790         197 :     if (r1 == 2) /* real quadratic */
    1791             :     {
    1792         185 :       cs[j] = 2 + ch_q(dtcr);
    1793         185 :       if (cs[j] == 2 || cs[j] == 4)
    1794             :       { /* NOT IMPLEMENTED YET */
    1795          30 :         GetST0(bnr, pS, pT, CR, prec);
    1796          30 :         return;
    1797             :       }
    1798             :       /* FIXME: is this value of N0 correct for the general case ? */
    1799         155 :       N0[j] = (long)prec * 0.35 * gtodouble(c);
    1800             :     }
    1801             :     else /* complex quadratic */
    1802             :     {
    1803          12 :       cs[j] = 1;
    1804          12 :       N0[j] = (long)prec * 0.7 * gtodouble(c);
    1805             :     }
    1806         167 :     if (n0 < N0[j]) n0 = N0[j];
    1807             :   }
    1808         147 :   if (DEBUGLEVEL>1) err_printf("N0 = %ld\n", n0);
    1809         147 :   InitPrimesQuad(bnr, n0, &LIST);
    1810             : 
    1811         147 :   av1 = avma;
    1812             :   /* loop over conductors */
    1813         314 :   for (j = 1; j <= ncond; j++)
    1814             :   {
    1815         167 :     GEN c0 = gel(C,j), c1 = divur(1, c0), c2 = divur(2, c0);
    1816         167 :     GEN ec1 = mpexp(c1), ec2 = mpexp(c2), LChar = gel(vChar,j);
    1817             :     GEN vf0, vf1, cf0, cf1;
    1818         167 :     const long nChar = lg(LChar)-1, NN = N0[j];
    1819             : 
    1820         167 :     if (DEBUGLEVEL>1)
    1821           0 :       err_printf("* conductor no %ld/%ld (N = %ld)\n\tInit: ", j,ncond,NN);
    1822         167 :     if (realprec(ec1) > prec) ec1 = rtor(ec1, prec);
    1823         167 :     if (realprec(ec2) > prec) ec2 = rtor(ec2, prec);
    1824         167 :     switch(cs[j])
    1825             :     {
    1826          12 :     case 1:
    1827          12 :       cf0 = gen_1;
    1828          12 :       cf1 = c0;
    1829          12 :       vf0 = mpveceint1(rtor(c1, prec), ec1, NN);
    1830          12 :       vf1 = mpvecpowdiv(invr(ec1), NN); break;
    1831             : 
    1832         155 :     case 3:
    1833         155 :       cf0 = sqrtr(mppi(prec));
    1834         155 :       cf1 = gmul2n(cf0, 1);
    1835         155 :       cf0 = gmul(cf0, c0);
    1836         155 :       vf0 = mpvecpowdiv(invr(ec2), NN);
    1837         155 :       vf1 = mpveceint1(rtor(c2, prec), ec2, NN); break;
    1838             : 
    1839           0 :     default:
    1840           0 :       cf0 = cf1 = NULL; /* FIXME: not implemented */
    1841           0 :       vf0 = vf1 = NULL;
    1842             :     }
    1843         699 :     for (k = 1; k <= nChar; k++)
    1844             :     {
    1845         532 :       long u = LChar[k], d, c;
    1846         532 :       GEN dtcr = gel(dataCR, u), z, s, t;
    1847             :       int **matan;
    1848             : 
    1849         532 :       if (!ch_comp(dtcr)) continue;
    1850         527 :       if (DEBUGLEVEL>1) err_printf("\tchar no: %ld (%ld/%ld)\n", u,k,nChar);
    1851         527 :       d = ch_phideg(dtcr);
    1852         527 :       z = gel(ch_CHI(dtcr), 2); s = t = gen_0; av2 = avma;
    1853         527 :       matan = computean(gel(dataCR,u), &LIST, NN, d);
    1854      892214 :       for (n = 1, c = 0; n <= NN; n++)
    1855      891687 :         if ((an = EvalCoeff(z, matan[n], d)))
    1856             :         {
    1857      234861 :           s = gadd(s, gmul(an, gel(vf0,n)));
    1858      234861 :           t = gadd(t, gmul(an, gel(vf1,n)));
    1859      234861 :           if (++c == 256) { (void)gc_all(av2,2, &s,&t); c = 0; }
    1860             :         }
    1861         527 :       affgc(gmul(cf0, s), gel(S,u));
    1862         527 :       affgc(gmul(cf1, conj_i(t)), gel(T,u));
    1863         527 :       FreeMat(matan,NN); set_avma(av2);
    1864             :     }
    1865         167 :     if (DEBUGLEVEL>1) err_printf("\n");
    1866         167 :     set_avma(av1);
    1867             :   }
    1868         147 :   set_avma(av);
    1869             : }
    1870             : 
    1871             : /* s += t*u. All 3 of them t_REAL, except we allow s or u = NULL (for 0) */
    1872             : static GEN
    1873    37342414 : _addmulrr(GEN s, GEN t, GEN u)
    1874             : {
    1875    37342414 :   if (u)
    1876             :   {
    1877    37146332 :     GEN v = mulrr(t, u);
    1878    37146332 :     return s? addrr(s, v): v;
    1879             :   }
    1880      196082 :   return s;
    1881             : }
    1882             : /* s += t. Both real, except we allow s or t = NULL (for exact 0) */
    1883             : static GEN
    1884    76280236 : _addrr(GEN s, GEN t)
    1885    76280236 : { return t? (s? addrr(s, t): t) : s; }
    1886             : 
    1887             : /* S & T for the general case. This is time-critical: optimize */
    1888             : static void
    1889      374677 : get_cS_cT(ST_t *T, long n)
    1890             : {
    1891             :   pari_sp av;
    1892             :   GEN csurn, nsurc, lncsurn, A, B, s, t, Z, aij, bij;
    1893             :   long i, j, r, i0;
    1894             : 
    1895      374677 :   if (gel(T->cS,n)) return;
    1896             : 
    1897      175712 :   av = avma;
    1898      175712 :   aij = T->aij; i0= T->i0;
    1899      175712 :   bij = T->bij; r = T->r;
    1900      175712 :   Z = cgetg(r+1, t_VEC);
    1901      175712 :   gel(Z,1) = NULL; /* unused */
    1902             : 
    1903      175712 :   csurn = divru(T->c1, n);
    1904      175712 :   nsurc = invr(csurn);
    1905      175712 :   lncsurn = logr_abs(csurn);
    1906             : 
    1907      175712 :   if (r > 1)
    1908             :   {
    1909      175562 :     gel(Z,2) = lncsurn; /* r >= 2 */
    1910      179072 :     for (i = 3; i <= r; i++)
    1911        3510 :       gel(Z,i) = divru(mulrr(gel(Z,i-1), lncsurn), i-1);
    1912             :     /* Z[i] = ln^(i-1)(c1/n) / (i-1)! */
    1913             :   }
    1914             : 
    1915             :   /* i = i0 */
    1916      175712 :     A = gel(aij,i0); t = _addrr(NULL, gel(A,1));
    1917      175712 :     B = gel(bij,i0); s = _addrr(NULL, gel(B,1));
    1918      354784 :     for (j = 2; j <= r; j++)
    1919             :     {
    1920      179072 :       s = _addmulrr(s, gel(Z,j),gel(B,j));
    1921      179072 :       t = _addmulrr(t, gel(Z,j),gel(A,j));
    1922             :     }
    1923    37876550 :   for (i = i0 - 1; i > 1; i--)
    1924             :   {
    1925    37700838 :     A = gel(aij,i); if (t) t = mulrr(t, nsurc);
    1926    37700838 :     B = gel(bij,i); if (s) s = mulrr(s, nsurc);
    1927    56013901 :     for (j = odd(i)? T->rc2: T->rc1; j > 1; j--)
    1928             :     {
    1929    18313063 :       s = _addmulrr(s, gel(Z,j),gel(B,j));
    1930    18313063 :       t = _addmulrr(t, gel(Z,j),gel(A,j));
    1931             :     }
    1932    37700838 :     s = _addrr(s, gel(B,1));
    1933    37700838 :     t = _addrr(t, gel(A,1));
    1934             :   }
    1935             :   /* i = 1 */
    1936      175712 :     A = gel(aij,1); if (t) t = mulrr(t, nsurc);
    1937      175712 :     B = gel(bij,1); if (s) s = mulrr(s, nsurc);
    1938      175712 :     s = _addrr(s, gel(B,1));
    1939      175712 :     t = _addrr(t, gel(A,1));
    1940      354784 :     for (j = 2; j <= r; j++)
    1941             :     {
    1942      179072 :       s = _addmulrr(s, gel(Z,j),gel(B,j));
    1943      179072 :       t = _addmulrr(t, gel(Z,j),gel(A,j));
    1944             :     }
    1945      175712 :   s = _addrr(s, T->b? mulrr(csurn, gel(T->powracpi,T->b+1)): csurn);
    1946      175712 :   if (!s) s = gen_0;
    1947      175712 :   if (!t) t = gen_0;
    1948      175712 :   gel(T->cS,n) = gclone(s);
    1949      175712 :   gel(T->cT,n) = gclone(t); set_avma(av);
    1950             : }
    1951             : 
    1952             : static void
    1953         383 : clear_cScT(ST_t *T, long N)
    1954             : {
    1955         383 :   GEN cS = T->cS, cT = T->cT;
    1956             :   long i;
    1957     1180061 :   for (i=1; i<=N; i++)
    1958     1179678 :     if (cS[i]) {
    1959      175712 :       gunclone(gel(cS,i));
    1960      175712 :       gunclone(gel(cT,i)); gel(cS,i) = gel(cT,i) = NULL;
    1961             :     }
    1962         383 : }
    1963             : 
    1964             : static void
    1965         241 : init_cScT(ST_t *T, GEN dtcr, long N, long prec)
    1966             : {
    1967         241 :   ch_get3(dtcr, &T->a, &T->b, &T->c);
    1968         241 :   T->rc1 = T->a + T->c;
    1969         241 :   T->rc2 = T->b + T->c;
    1970         241 :   T->r   = maxss(T->rc2+1, T->rc1); /* >= 2 */
    1971         241 :   ppgamma(T, prec);
    1972         241 :   clear_cScT(T, N);
    1973         241 : }
    1974             : 
    1975             : /* return a t_REAL */
    1976             : static GEN
    1977         377 : zeta_get_limx(long r1, long r2, long bit)
    1978             : {
    1979         377 :   pari_sp av = avma;
    1980             :   GEN p1, p2, c0, c1, A0;
    1981         377 :   long r = r1 + r2, N = r + r2;
    1982             : 
    1983             :   /* c1 = N 2^(-2r2 / N) */
    1984         377 :   c1 = mulrs(powrfrac(real2n(1, DEFAULTPREC), -2*r2, N), N);
    1985             : 
    1986         377 :   p1 = powru(Pi2n(1, DEFAULTPREC), r - 1);
    1987         377 :   p2 = mulir(powuu(N,r), p1); shiftr_inplace(p2, -r2);
    1988         377 :   c0 = sqrtr( divrr(p2, powru(c1, r+1)) );
    1989             : 
    1990         377 :   A0 = logr_abs( gmul2n(c0, bit) ); p2 = divrr(A0, c1);
    1991         377 :   p1 = divrr(mulur(N*(r+1), logr_abs(p2)), addsr(2*(r+1), gmul2n(A0,2)));
    1992         377 :   return gc_leaf(av, divrr(addrs(p1, 1), powruhalf(p2, N)));
    1993             : }
    1994             : /* N_0 = floor( C_K / limx ). Large */
    1995             : static long
    1996         476 : zeta_get_N0(GEN C,  GEN limx)
    1997             : {
    1998             :   long e;
    1999         476 :   pari_sp av = avma;
    2000         476 :   GEN z = gcvtoi(gdiv(C, limx), &e); /* avoid truncation error */
    2001         476 :   if (e >= 0 || is_bigint(z))
    2002           0 :     pari_err_OVERFLOW("zeta_get_N0 [need too many primes]");
    2003         476 :   return gc_long(av, itos(z));
    2004             : }
    2005             : 
    2006             : static GEN
    2007        1426 : eval_i(long r1, long r2, GEN limx, long i)
    2008             : {
    2009        1426 :   GEN t = powru(limx, i);
    2010        1426 :   if (!r1)      t = mulrr(t, powru(mpfactr(i  , DEFAULTPREC), r2));
    2011        1426 :   else if (!r2) t = mulrr(t, powru(mpfactr(i/2, DEFAULTPREC), r1));
    2012             :   else {
    2013           0 :     GEN u1 = mpfactr(i/2, DEFAULTPREC);
    2014           0 :     GEN u2 = mpfactr(i,   DEFAULTPREC);
    2015           0 :     if (r1 == r2) t = mulrr(t, powru(mulrr(u1,u2), r1));
    2016           0 :     else t = mulrr(t, mulrr(powru(u1,r1), powru(u2,r2)));
    2017             :   }
    2018        1426 :   return t;
    2019             : }
    2020             : 
    2021             : /* "small" even i such that limx^i ( (i\2)! )^r1 ( i! )^r2 > B. */
    2022             : static long
    2023         142 : get_i0(long r1, long r2, GEN B, GEN limx)
    2024             : {
    2025         142 :   long imin = 1, imax = 1400;
    2026         148 :   while (mpcmp(eval_i(r1,r2,limx, imax), B) < 0) { imin = imax; imax *= 2; }
    2027        1420 :   while(imax - imin >= 4)
    2028             :   {
    2029        1278 :     long m = (imax + imin) >> 1;
    2030        1278 :     if (mpcmp(eval_i(r1,r2,limx, m), B) >= 0) imax = m; else imin = m;
    2031             :   }
    2032         142 :   return imax & ~1; /* make it even */
    2033             : }
    2034             : /* limx = zeta_get_limx(r1, r2, bit), a t_REAL */
    2035             : static long
    2036         142 : zeta_get_i0(long r1, long r2, long bit, GEN limx)
    2037             : {
    2038         142 :   pari_sp av = avma;
    2039         142 :   GEN B = gmul(sqrtr( divrr(powrs(mppi(DEFAULTPREC), r2-3), limx) ),
    2040             :                gmul2n(powuu(5, r1), bit + r2));
    2041         142 :   return gc_long(av, get_i0(r1, r2, B, limx));
    2042             : }
    2043             : 
    2044             : static void
    2045         142 : GetST0(GEN bnr, GEN *pS, GEN *pT, GEN CR, long prec)
    2046             : {
    2047             :   pari_sp av, av1, av2;
    2048             :   long n, j, k, jc, n0, prec2, i0, r1, r2, ncond;
    2049         142 :   GEN nf = bnr_get_nf(bnr);
    2050         142 :   GEN vChar = gel(CR,1), dataCR = gel(CR,2), N0, C, an, limx, S, T;
    2051             :   LISTray LIST;
    2052             :   ST_t cScT;
    2053             : 
    2054         142 :   ST_alloc(pS, pT, lg(dataCR), prec); T = *pT; S = *pS;
    2055         142 :   av = avma;
    2056         142 :   nf_get_sign(nf,&r1,&r2);
    2057         142 :   ncond = lg(vChar)-1;
    2058         142 :   C  = cgetg(ncond+1, t_VEC);
    2059         142 :   N0 = cgetg(ncond+1, t_VECSMALL);
    2060         142 :   n0 = 0;
    2061         142 :   limx = zeta_get_limx(r1, r2, prec);
    2062         383 :   for (j = 1; j <= ncond; j++)
    2063             :   {
    2064         241 :     GEN dtcr = gel(dataCR, mael(vChar,j,1)), c = ch_C(dtcr);
    2065         241 :     gel(C,j) = c;
    2066         241 :     N0[j] = zeta_get_N0(c, limx);
    2067         241 :     if (n0 < N0[j]) n0  = N0[j];
    2068             :   }
    2069         142 :   cScT.i0 = i0 = zeta_get_i0(r1, r2, prec, limx);
    2070         142 :   if (DEBUGLEVEL>1) err_printf("i0 = %ld, N0 = %ld\n",i0, n0);
    2071         142 :   InitPrimes(bnr, n0, &LIST);
    2072         142 :   prec2 = precdbl(prec) + EXTRAPREC64;
    2073         142 :   cScT.powracpi = powersr(sqrtr(mppi(prec2)), r1);
    2074         142 :   cScT.cS = cgetg(n0+1, t_VEC);
    2075         142 :   cScT.cT = cgetg(n0+1, t_VEC);
    2076      568817 :   for (j=1; j<=n0; j++) gel(cScT.cS,j) = gel(cScT.cT,j) = NULL;
    2077             : 
    2078         142 :   av1 = avma;
    2079         383 :   for (jc = 1; jc <= ncond; jc++)
    2080             :   {
    2081         241 :     GEN LChar = gel(vChar,jc);
    2082         241 :     long nChar = lg(LChar)-1, N = N0[jc];
    2083             : 
    2084             :     /* Can discard completely a conductor if all chars satisfy L(0,chi) = 0
    2085             :      * Not sure whether this is possible. */
    2086         241 :     if (DEBUGLEVEL>1)
    2087           0 :       err_printf("* conductor no %ld/%ld (N = %ld)\n\tInit: ", jc,ncond,N);
    2088             : 
    2089         241 :     cScT.c1 = gel(C,jc);
    2090         241 :     init_cScT(&cScT, gel(dataCR, LChar[1]), N, prec2);
    2091         241 :     av2 = avma;
    2092         678 :     for (k = 1; k <= nChar; k++)
    2093             :     {
    2094         437 :       long d, c, u = LChar[k];
    2095         437 :       GEN dtcr = gel(dataCR, u), z, s, t;
    2096             :       int **matan;
    2097             : 
    2098         437 :       if (!ch_comp(dtcr)) continue;
    2099         432 :       if (DEBUGLEVEL>1) err_printf("\tchar no: %ld (%ld/%ld)\n", u,k,nChar);
    2100         432 :       z = gel(ch_CHI(dtcr), 2);
    2101         432 :       d = ch_phideg(dtcr); s = t = gen_0;
    2102         432 :       matan = ComputeCoeff(dtcr, &LIST, N, d);
    2103     1503444 :       for (n = 1, c = 0; n <= N; n++)
    2104     1503012 :         if ((an = EvalCoeff(z, matan[n], d)))
    2105             :         {
    2106      374677 :           get_cS_cT(&cScT, n);
    2107      374677 :           s = gadd(s, gmul(an, gel(cScT.cS,n)));
    2108      374677 :           t = gadd(t, gmul(an, gel(cScT.cT,n)));
    2109      374677 :           if (++c == 256) { (void)gc_all(av2,2, &s,&t); c = 0; }
    2110             :         }
    2111         432 :       affgc(s,         gel(S,u));
    2112         432 :       affgc(conj_i(t), gel(T,u));
    2113         432 :       FreeMat(matan, N); set_avma(av2);
    2114             :     }
    2115         241 :     if (DEBUGLEVEL>1) err_printf("\n");
    2116         241 :     set_avma(av1);
    2117             :   }
    2118         142 :   clear_cScT(&cScT, n0);
    2119         142 :   set_avma(av);
    2120         142 : }
    2121             : 
    2122             : static void
    2123         289 : GetST(GEN bnr, GEN *pS, GEN *pT, GEN CR, long prec)
    2124             : {
    2125         289 :   if (nf_get_degree(bnr_get_nf(bnr)) == 2)
    2126         177 :     QuadGetST(bnr, pS, pT, CR, prec);
    2127             :   else
    2128         112 :     GetST0(bnr, pS, pT, CR, prec);
    2129         289 : }
    2130             : 
    2131             : /*******************************************************************/
    2132             : /*                                                                 */
    2133             : /*     Class fields of real quadratic fields using Stark units     */
    2134             : /*                                                                 */
    2135             : /*******************************************************************/
    2136             : /* compute the Hilbert class field using genus class field theory when
    2137             :    the exponent of the class group is exactly 2 (trivial group not covered) */
    2138             : /* Cf Herz, Construction of class fields, LNM 21, Theorem 1 (VII-6) */
    2139             : static GEN
    2140          10 : GenusFieldQuadReal(GEN disc)
    2141             : {
    2142          10 :   GEN T = NULL, p0 = NULL, P = gel(Z_factor(disc), 1);
    2143          10 :   long i, i0 = 0, l = lg(P);
    2144             : 
    2145          30 :   for (i = 1; i < l; i++)
    2146             :   {
    2147          25 :     GEN p = gel(P,i);
    2148          25 :     if (mod4(p) == 3) { p0 = p; i0 = i; break; }
    2149             :   }
    2150          10 :   l--; /* remove last prime */
    2151          10 :   if (i0 == l) l--; /* ... remove p0 and last prime */
    2152          35 :   for (i = 1; i < l; i++)
    2153             :   {
    2154          25 :     GEN p = gel(P,i), d, t;
    2155          25 :     if (i == i0) continue;
    2156          20 :     if (absequaliu(p, 2))
    2157          10 :       switch (mod32(disc))
    2158             :       {
    2159          10 :         case  8: d = gen_2; break;
    2160           0 :         case 24: d = shifti(p0, 1); break;
    2161           0 :         default: d = p0; break;
    2162             :       }
    2163             :     else
    2164          10 :       d = (mod4(p) == 1)? p: mulii(p0, p);
    2165          20 :     t = mkpoln(3, gen_1, gen_0, negi(d)); /* x^2 - d */
    2166          20 :     T = T? ZX_composedsum(T, t): t;
    2167             :   }
    2168          10 :   return polredbest(T, 0);
    2169             : }
    2170             : static GEN
    2171         290 : GenusFieldQuadImag(GEN disc)
    2172             : {
    2173         290 :   GEN T = NULL, P = gel(absZ_factor(disc), 1);
    2174         290 :   long i, n = lg(P)-1;
    2175         845 :   for (i = 1; i < n; i++) /* remove last prime */
    2176             :   {
    2177         555 :     GEN p = gel(P,i), d, t;
    2178         555 :     if (absequaliu(p, 2))
    2179         165 :       switch (mod32(disc))
    2180             :       {
    2181          40 :         case 24: d = gen_2; break;  /* disc = 8 mod 32 */
    2182          30 :         case  8: d = gen_m2; break; /* disc =-8 mod 32 */
    2183          95 :         default: d = gen_m1; break;
    2184             :       }
    2185             :     else
    2186         390 :       d = (mod4(p) == 1)? p: negi(p);
    2187         555 :     t = mkpoln(3, gen_1, gen_0, negi(d)); /* x^2 - d */
    2188         555 :     T = T? ZX_composedsum(T, t): t;
    2189             :   }
    2190         290 :   return polredbest(T, 0);
    2191             : }
    2192             : 
    2193             : /* if flag = 1, computes a fast and crude approximation of the trace of the Stark unit
    2194             :    if flag = 2, computes the Stark unit */
    2195             : static GEN
    2196         470 : AllStark(GEN data, long flag, long newprec)
    2197             : {
    2198         470 :   const long BND = 300;
    2199         470 :   long cl, i, j, cpt = 0, N, h, v, n, r1, r2, den;
    2200             :   pari_sp av, av2;
    2201             :   int **matan;
    2202         470 :   GEN bnr = gel(data,1), nf = bnr_get_nf(bnr), p1, p2, S, T;
    2203         470 :   GEN CR = gel(data,4), dataCR = gel(CR,2);
    2204             :   GEN polrelnum, polrel, Lp, W, vzeta, C, cond1, L1, an;
    2205             :   LISTray LIST;
    2206             :   pari_timer ti;
    2207             : 
    2208         470 :   nf_get_sign(nf, &r1,&r2);
    2209         470 :   N     = nf_get_degree(nf);
    2210         470 :   cond1 = gel(bnr_get_mod(bnr), 2);
    2211             : 
    2212         470 :   v = 1;
    2213         980 :   while (gequal1(gel(cond1,v))) v++;
    2214         470 :   cl = lg(dataCR)-1;
    2215         470 :   h  = itos(ZM_det_triangular(gel(data,2))) >> 1;
    2216             : 
    2217         470 : LABDOUB:
    2218         470 :   if (DEBUGLEVEL) timer_start(&ti);
    2219         470 :   av = avma;
    2220         470 :   W = AllArtinNumbers(CR, newprec);
    2221         470 :   if (DEBUGLEVEL) timer_printf(&ti,"Compute W");
    2222         470 :   Lp = cgetg(cl + 1, t_VEC);
    2223         470 :   if (flag != 1)
    2224             :   {
    2225         235 :     GetST(bnr, &S, &T, CR, newprec);
    2226         235 :     if (DEBUGLEVEL) timer_printf(&ti, "S&T");
    2227        1000 :     for (i = 1; i <= cl; i++)
    2228             :     {
    2229         765 :       GEN chi = gel(dataCR, i), vv = gen_0;
    2230         765 :       if (ch_comp(chi))
    2231         755 :         vv = gel(GetValue(chi, gel(W,i), gel(S,i), gel(T,i), 2, newprec), 2);
    2232         765 :       gel(Lp, i) = vv;
    2233             :     }
    2234             :   }
    2235             :   else
    2236             :   { /* compute a crude approximation of the result */
    2237         235 :     C = cgetg(cl + 1, t_VEC);
    2238        1000 :     for (i = 1; i <= cl; i++) gel(C,i) = ch_C(gel(dataCR, i));
    2239         235 :     n = zeta_get_N0(vecmax(C), zeta_get_limx(r1, r2, newprec));
    2240         235 :     if (n > BND) n = BND;
    2241         235 :     if (DEBUGLEVEL) err_printf("N0 in QuickPol: %ld \n", n);
    2242         235 :     InitPrimes(bnr, n, &LIST);
    2243             : 
    2244         235 :     L1 = cgetg(cl+1, t_VEC);
    2245             :     /* use L(1) = sum (an / n) */
    2246        1000 :     for (i = 1; i <= cl; i++)
    2247             :     {
    2248         765 :       GEN dtcr = gel(dataCR,i);
    2249         765 :       long d = ch_phideg(dtcr);
    2250         765 :       matan = ComputeCoeff(dtcr, &LIST, n, d);
    2251         765 :       av2 = avma;
    2252         765 :       p1 = real_0(newprec); p2 = gel(ch_CHI(dtcr), 2);
    2253      224285 :       for (j = 1; j <= n; j++)
    2254      223520 :         if ( (an = EvalCoeff(p2, matan[j], d)) )
    2255       85840 :           p1 = gadd(p1, gdivgu(an, j));
    2256         765 :       gel(L1,i) = gc_upto(av2, p1);
    2257         765 :       FreeMat(matan, n);
    2258             :     }
    2259         235 :     p1 = gmul2n(powruhalf(mppi(newprec), N-2), 1);
    2260             : 
    2261        1000 :     for (i = 1; i <= cl; i++)
    2262             :     {
    2263             :       long r;
    2264         765 :       GEN WW, A = AChi(gel(dataCR,i), &r, 0, newprec);
    2265         765 :       WW = gmul(gel(C,i), gmul(A, gel(W,i)));
    2266         765 :       gel(Lp,i) = gdiv(gmul(WW, conj_i(gel(L1,i))), p1);
    2267             :     }
    2268             :   }
    2269             : 
    2270         470 :   p1 = gel(data,3);
    2271         470 :   den = (flag == 0) ? 2*h: h;
    2272         470 :   if (flag == 2)
    2273           5 :     vzeta = cgetg(2*h + 1, t_VEC);
    2274             :   else
    2275         465 :     vzeta = cgetg(h + 1,t_VEC);
    2276        3020 :   for (i = 1; i <= h; i++)
    2277             :   {
    2278        2550 :     GEN z = gen_0, sig = gel(p1,i);
    2279       12980 :     for (j = 1; j <= cl; j++)
    2280             :     {
    2281       10430 :       GEN dtcr = gel(dataCR,j), CHI = ch_CHI(dtcr);
    2282       10430 :       GEN t = mulreal(gel(Lp,j), CharEval(CHI, sig));
    2283       10430 :       if (chi_get_deg(CHI) != 2) t = gmul2n(t, 1); /* character not real */
    2284       10430 :       z = gadd(z, t);
    2285             :     }
    2286        2550 :     z = gdivgu(z, den);
    2287        2550 :     if (flag == 2)
    2288             :     {
    2289          55 :       gel(vzeta, i) = gexp(gmul2n(z,1), newprec);
    2290          55 :       gel(vzeta, i+h) = ginv(gel(vzeta,i));
    2291             :     }
    2292             :     /* if flag == 0, we first try with the square-root of the Stark unit */
    2293        2495 :     else gel(vzeta,i) = gmul2n(gcosh(z, newprec), 1);
    2294             :   }
    2295         470 :   polrelnum = roots_to_pol(vzeta, 0);
    2296         470 :   if (DEBUGLEVEL)
    2297             :   {
    2298           0 :     if (DEBUGLEVEL>1) {
    2299           0 :       err_printf("polrelnum = %Ps\n", polrelnum);
    2300           0 :       err_printf("zetavalues = %Ps\n", vzeta);
    2301           0 :       if (flag == 0)
    2302           0 :         err_printf("Checking the square-root of the Stark unit...\n");
    2303             :     }
    2304           0 :     timer_printf(&ti, "Compute %s", flag? "quickpol": "polrelnum");
    2305             :   }
    2306         470 :   if (flag == 1) return gc_GEN(av, polrelnum);
    2307             :   /* try to recognize this polynomial */
    2308         235 :   polrel = RecCoeff(nf, polrelnum, v, flag, newprec);
    2309             : 
    2310             :   /* if this doesn't work, maybe the Stark unit is not a square */
    2311         235 :   if (!polrel && flag == 0)
    2312             :   {
    2313         155 :     if (DEBUGLEVEL)
    2314             :     {
    2315           0 :       if (DEBUGLEVEL>1) {
    2316           0 :         err_printf("It's not a square...\n");
    2317           0 :         err_printf("polrelnum = %Ps\n", polrelnum);
    2318             :       }
    2319           0 :       timer_printf(&ti, "Compute polrelnum");
    2320             :     }
    2321        1080 :     for (j = 1; j <= h; j++)
    2322         925 :       gel(vzeta,j) = gsubgs(gsqr(gel(vzeta,j)), 2);
    2323         155 :     polrelnum = roots_to_pol(vzeta, 0);
    2324             : 
    2325         155 :     polrel = RecCoeff(nf, polrelnum, v, 0, newprec);
    2326             :   }
    2327         235 :   if (!polrel) /* FAILED */
    2328             :   {
    2329           0 :     const long EXTRA_BITS = 64;
    2330             :     long incr_pr;
    2331           0 :     if (++cpt >= 3) pari_err_PREC( "stark (computation impossible)");
    2332             :     /* estimate needed precision */
    2333           0 :     incr_pr = gprecision(polrelnum) - gexpo(polrelnum);
    2334           0 :     if (incr_pr < 0) incr_pr = -incr_pr + EXTRA_BITS;
    2335           0 :     newprec += nbits2extraprec(maxss(3*EXTRA_BITS, cpt*incr_pr));
    2336           0 :     if (DEBUGLEVEL) pari_warn(warnprec, "AllStark", newprec);
    2337           0 :     CharNewPrec(data, newprec);
    2338           0 :     nf = bnr_get_nf(ch_bnr(gel(dataCR,1)));
    2339           0 :     goto LABDOUB;
    2340             :   }
    2341             : 
    2342         235 :   if (DEBUGLEVEL) {
    2343           0 :     if (DEBUGLEVEL>1) err_printf("polrel = %Ps\n", polrel);
    2344           0 :     timer_printf(&ti, "Recpolnum");
    2345             :   }
    2346         235 :   return gc_GEN(av, polrel);
    2347             : }
    2348             : 
    2349             : /********************************************************************/
    2350             : /*                        Main functions                            */
    2351             : /********************************************************************/
    2352             : static GEN
    2353           0 : bnrstark_cyclic(GEN bnr, GEN dtQ, long prec)
    2354             : {
    2355           0 :   GEN v, vH, cyc = gel(dtQ,2), U = gel(dtQ,3), M = ZM_inv(U, NULL);
    2356           0 :   long i, j, l = lg(M);
    2357             : 
    2358             :   /* M = indep. generators of Cl_f/subgp, restrict to cyclic components */
    2359           0 :   vH = cgetg(l, t_VEC);
    2360           0 :   for (i = j = 1; i < l; i++)
    2361             :   {
    2362           0 :     if (is_pm1(gel(cyc,i))) break;
    2363           0 :     gel(vH, j++) = ZM_hnfmodid(vecsplice(M,i), cyc);
    2364             :   }
    2365           0 :   setlg(vH, j); v = cgetg(l, t_VEC);
    2366           0 :   for (i = 1; i < j; i++) gel(v,i) = bnrstark(bnr, gel(vH,i), prec);
    2367           0 :   return v;
    2368             : }
    2369             : GEN
    2370         147 : bnrstark(GEN bnr, GEN subgrp, long prec)
    2371             : {
    2372             :   long newprec;
    2373         147 :   pari_sp av = avma;
    2374             :   GEN nf, data, dtQ;
    2375             : 
    2376             :   /* check the bnr */
    2377         147 :   checkbnr(bnr); nf  = bnr_get_nf(bnr);
    2378         147 :   if (nf_get_degree(nf) == 1) return galoissubcyclo(bnr, subgrp, 0, 0);
    2379         147 :   if (!nf_get_varn(nf))
    2380           0 :     pari_err_PRIORITY("bnrstark", nf_get_pol(nf), "=", 0);
    2381         147 :   if (nf_get_r2(nf)) pari_err_DOMAIN("bnrstark", "r2", "!=", gen_0, nf);
    2382             : 
    2383             :   /* compute bnr(conductor) */
    2384         141 :   bnr_subgroup_sanitize(&bnr, &subgrp);
    2385         141 :   if (gequal1(ZM_det_triangular(subgrp))) { set_avma(av); return pol_x(0); }
    2386             : 
    2387             :   /* check the class field */
    2388         141 :   if (!gequal0(gel(bnr_get_mod(bnr), 2)))
    2389           6 :     pari_err_DOMAIN("bnrstark", "r2(class field)", "!=", gen_0, bnr);
    2390             : 
    2391             :   /* find a suitable extension N */
    2392         135 :   dtQ = InitQuotient(subgrp);
    2393         135 :   data  = FindModulus(bnr, dtQ, &newprec);
    2394         135 :   if (!data) return gc_upto(av, bnrstark_cyclic(bnr, dtQ, prec));
    2395         135 :   if (DEBUGLEVEL>1 && newprec > prec)
    2396           0 :     err_printf("new precision: %ld\n", newprec);
    2397         135 :   return gc_upto(av, AllStark(data, 0, newprec));
    2398             : }
    2399             : 
    2400             : GEN
    2401          10 : bnrstarkunit(GEN bnr, GEN subgrp)
    2402             : {
    2403             :   long newprec, deg;
    2404          10 :   pari_sp av = avma;
    2405             :   GEN map, nf, data, dtQ, bnf, bnrf, Cm, candD, D, QD, CR;
    2406             : 
    2407             :   /* check the input */
    2408          10 :   checkbnr(bnr); bnf = bnr_get_bnf(bnr); nf = bnf_get_nf(bnf);
    2409          10 :   if (!nf_get_varn(nf))
    2410           0 :     pari_err_PRIORITY("bnrstarkunit", nf_get_pol(nf), "=", 0);
    2411          10 :   deg = nf_get_degree(nf);
    2412          10 :   if (deg == 1) pari_err_IMPL("bnrstarkunit for basefield Q");
    2413          10 :   if (nf_get_r2(nf)) pari_err_DOMAIN("bnrstarkunit", "r2", "!=", gen_0, nf);
    2414          10 :   bnr_subgroup_sanitize(&bnr, &subgrp);
    2415          10 :   if (lg(bid_get_archp(bnr_get_bid(bnr)))-1 != deg-1)
    2416           5 :     pari_err_DOMAIN("bnrstarkunit", "# unramified places", "!=", gen_1, bnr);
    2417           5 :   bnrf = Buchray(bnf, gel(bnr_get_mod(bnr), 1), nf_INIT);
    2418           5 :   map = bnrsurjection(bnr, bnrf);
    2419           5 :   subgrp = abmap_subgroup_image(map, subgrp);
    2420           5 :   dtQ  = InitQuotient(subgrp);
    2421           5 :   Cm   = ComputeKernel(map, dtQ);
    2422           5 :   candD = subgrouplist_cond_sub(bnr, Cm, mkvec(gen_2));
    2423           5 :   if (lg(candD) != 2) pari_err(e_MISC, "incorrect modulus in bnrstark");
    2424           5 :   D    = gel(candD, 1);
    2425           5 :   QD   = InitQuotient(D);
    2426           5 :   CR   = InitChar(bnr, AllChars(bnr, QD, 1), 0, DEFAULTPREC);
    2427           5 :   data = mkvec4(bnr, D, subgroup_classes(Cm), CR);
    2428           5 :   CplxModulus(data, &newprec);
    2429           5 :   return gc_upto(av, AllStark(data, 2, newprec));
    2430             : }
    2431             : 
    2432             : /* For each character of Cl(bnr)/subgp, compute L(1, chi) (or equivalently
    2433             :  * the first nonzero term c(chi) of the expansion at s = 0).
    2434             :  * If flag & 1: compute the value at s = 1 (for nontrivial characters),
    2435             :  * else compute the term c(chi) and return [r(chi), c(chi)] where r(chi) is
    2436             :  *   the order of L(s, chi) at s = 0.
    2437             :  * If flag & 2: compute the value of the L-function L_S(s, chi) where S is the
    2438             :  *   set of places dividing the modulus of bnr (and the infinite places),
    2439             :  * else
    2440             :  *   compute the value of the primitive L-function attached to chi,
    2441             :  * If flag & 4: return also the character */
    2442             : GEN
    2443          60 : bnrL1(GEN bnr, GEN subgp, long flag, long prec)
    2444             : {
    2445             :   GEN L1, ch, Qt, z;
    2446             :   long l, h;
    2447          60 :   pari_sp av = avma;
    2448             : 
    2449          60 :   checkbnr(bnr);
    2450          60 :   if (flag < 0 || flag > 8) pari_err_FLAG("bnrL1");
    2451             : 
    2452          60 :   subgp = bnr_subgroup_check(bnr, subgp, NULL);
    2453          60 :   if (!subgp) subgp = diagonal_shallow(bnr_get_cyc(bnr));
    2454             : 
    2455          60 :   Qt = InitQuotient(subgp);
    2456          60 :   ch = AllChars(bnr, Qt, 0); l = lg(ch);
    2457          60 :   h = itou(gel(Qt,1));
    2458          60 :   L1 = cgetg((flag&1)? h: h+1, t_VEC);
    2459          60 :   if (l > 1)
    2460             :   {
    2461          54 :     GEN W, S, T, CR = InitChar(bnr, ch, 1, prec), dataCR = gel(CR,2);
    2462             :     long i, j;
    2463             : 
    2464          54 :     GetST(bnr, &S, &T, CR, prec);
    2465          54 :     W = AllArtinNumbers(CR, prec);
    2466         258 :     for (i = j = 1; i < l; i++)
    2467             :     {
    2468         204 :       GEN chi = gel(ch,i);
    2469         204 :       z = GetValue(gel(dataCR,i), gel(W,i), gel(S,i), gel(T,i), flag, prec);
    2470         204 :       gel(L1,j++) = (flag & 4)? mkvec2(gel(chi,1), z): z;
    2471         204 :       if (lg(chi) == 4)
    2472             :       { /* nonreal */
    2473         114 :         z = conj_i(z);
    2474         114 :         gel(L1, j++) = (flag & 4)? mkvec2(gel(chi,3), z): z;
    2475             :       }
    2476             :     }
    2477             :   }
    2478          60 :   if (!(flag & 1))
    2479             :   { /* trivial character */
    2480          30 :     z = GetValue1(bnr, flag & 2, prec);
    2481          30 :     if (flag & 4)
    2482             :     {
    2483           0 :       GEN chi = zerovec(lg(bnr_get_cyc(bnr))-1);
    2484           0 :       z = mkvec2(chi, z);
    2485             :     }
    2486          30 :     gel(L1,h) = z;
    2487             :   }
    2488          60 :   return gc_GEN(av, L1);
    2489             : }
    2490             : 
    2491             : /*******************************************************************/
    2492             : /*                                                                 */
    2493             : /*       Hilbert and Ray Class field using Stark                   */
    2494             : /*                                                                 */
    2495             : /*******************************************************************/
    2496             : /* P in A[x,y], deg_y P < 2, return P0 and P1 in A[x] such that P = P0 + P1 y */
    2497             : static void
    2498          95 : split_pol_quad(GEN P, GEN *gP0, GEN *gP1)
    2499             : {
    2500          95 :   long i, l = lg(P);
    2501          95 :   GEN P0 = cgetg(l, t_POL), P1 = cgetg(l, t_POL);
    2502          95 :   P0[1] = P1[1] = P[1];
    2503         865 :   for (i = 2; i < l; i++)
    2504             :   {
    2505         770 :     GEN c = gel(P,i), c0 = c, c1 = gen_0;
    2506         770 :     if (typ(c) == t_POL) /* write c = c1 y + c0 */
    2507         675 :       switch(degpol(c))
    2508             :       {
    2509           0 :         case -1: c0 = gen_0; break;
    2510         675 :         default: c1 = gel(c,3); /* fall through */
    2511         675 :         case  0: c0 = gel(c,2); break;
    2512             :       }
    2513         770 :     gel(P0,i) = c0; gel(P1,i) = c1;
    2514             :   }
    2515          95 :   *gP0 = normalizepol_lg(P0, l);
    2516          95 :   *gP1 = normalizepol_lg(P1, l);
    2517          95 : }
    2518             : 
    2519             : /* k = nf quadratic field, P relative equation of H_k (Hilbert class field)
    2520             :  * return T in Z[X], such that H_k / Q is the compositum of Q[X]/(T) and k */
    2521             : static GEN
    2522          95 : makescind(GEN nf, GEN P)
    2523             : {
    2524          95 :   GEN Pp, p, pol, G, L, a, roo, P0,P1, Ny,Try, nfpol = nf_get_pol(nf);
    2525             :   long i, is_P;
    2526             : 
    2527          95 :   P = lift_shallow(P);
    2528          95 :   split_pol_quad(P, &P0, &P1);
    2529             :   /* P = P0 + y P1, Norm_{k/Q}(P) = P0^2 + Tr y P0P1 + Ny P1^2, irreducible/Q */
    2530          95 :   Ny = gel(nfpol, 2);
    2531          95 :   Try = negi(gel(nfpol, 3));
    2532          95 :   pol = RgX_add(RgX_sqr(P0), RgX_Rg_mul(RgX_sqr(P1), Ny));
    2533          95 :   if (signe(Try)) pol = RgX_add(pol, RgX_Rg_mul(RgX_mul(P0,P1), Try));
    2534             :   /* pol = rnfequation(nf, P); */
    2535          95 :   G = galoisinit(pol, NULL);
    2536          95 :   L = gal_get_group(G);
    2537          95 :   p = gal_get_p(G);
    2538          95 :   a = FpX_oneroot(nfpol, p);
    2539             :   /* P mod a prime \wp above p (which splits) */
    2540          95 :   Pp = FpXY_evalx(P, a, p);
    2541          95 :   roo = gal_get_roots(G);
    2542          95 :   is_P = gequal0( FpX_eval(Pp, remii(gel(roo,1),p), p) );
    2543             :   /* each roo[i] mod p is a root of P or (exclusive) tau(P) mod \wp */
    2544             :   /* record whether roo[1] is a root of P or tau(P) */
    2545             : 
    2546         730 :   for (i = 1; i < lg(L); i++)
    2547             :   {
    2548         730 :     GEN perm = gel(L,i);
    2549         730 :     long k = perm[1]; if (k == 1) continue;
    2550         635 :     k = gequal0( FpX_eval(Pp, remii(gel(roo,k),p), p) );
    2551             :     /* roo[k] is a root of the other polynomial */
    2552         635 :     if (k != is_P)
    2553             :     {
    2554          95 :       ulong o = perm_orderu(perm);
    2555          95 :       if (o != 2) perm = perm_powu(perm, o >> 1);
    2556             :       /* perm has order two and doesn't belong to Gal(H_k/k) */
    2557          95 :       return polredbest(galoisfixedfield(G, perm, 1, varn(P)), 0);
    2558             :     }
    2559             :   }
    2560           0 :   pari_err_BUG("makescind");
    2561             :   return NULL; /*LCOV_EXCL_LINE*/
    2562             : }
    2563             : 
    2564             : /* pbnf = NULL if no bnf is needed, f = NULL may be passed for a trivial
    2565             :  * conductor */
    2566             : static void
    2567         607 : quadray_init(GEN *pD, GEN *pbnf, long prec)
    2568             : {
    2569         607 :   GEN D = *pD, nf, bnf = NULL;
    2570         607 :   if (typ(D) == t_INT)
    2571             :   {
    2572             :     int isfund;
    2573         581 :     if (pbnf) {
    2574         180 :       bnf = Buchall(quadpoly0(D, 1), nf_FORCE, prec);
    2575         180 :       nf = bnf_get_nf(bnf);
    2576         180 :       isfund = equalii(D, nf_get_disc(nf));
    2577             :     }
    2578             :     else
    2579         401 :       isfund = Z_isfundamental(D);
    2580         581 :     if (!isfund) pari_err_DOMAIN("quadray", "isfundamental(D)", "=",gen_0, D);
    2581             :   }
    2582             :   else
    2583             :   {
    2584          26 :     bnf = checkbnf(D);
    2585          26 :     nf = bnf_get_nf(bnf);
    2586          26 :     if (nf_get_degree(nf) != 2)
    2587           6 :       pari_err_DOMAIN("quadray", "degree", "!=", gen_2, nf_get_pol(nf));
    2588          20 :     D = nf_get_disc(nf);
    2589             :   }
    2590         595 :   if (pbnf) *pbnf = bnf;
    2591         595 :   *pD = D;
    2592         595 : }
    2593             : 
    2594             : /* compute the polynomial over Q of the Hilbert class field of
    2595             :    Q(sqrt(D)) where D is a positive fundamental discriminant */
    2596             : static GEN
    2597         105 : quadhilbertreal(GEN D, long prec)
    2598             : {
    2599             :   GEN bnf, bnr, dtQ, data, M;
    2600             :   long newprec;
    2601             :   pari_timer T;
    2602             : 
    2603         105 :   quadray_init(&D, &bnf, prec);
    2604         105 :   switch(itou_or_0(cyc_get_expo(bnf_get_cyc(bnf))))
    2605             :   {
    2606           0 :     case 1: return pol_x(0);
    2607          10 :     case 2: return GenusFieldQuadReal(D);
    2608             :   }
    2609          95 :   bnr  = Buchray(bnf, gen_1, nf_INIT);
    2610          95 :   M = diagonal_shallow(bnr_get_cyc(bnr));
    2611          95 :   dtQ = InitQuotient(M);
    2612             : 
    2613          95 :   if (DEBUGLEVEL) timer_start(&T);
    2614          95 :   data = FindModulus(bnr, dtQ, &newprec);
    2615          95 :   if (DEBUGLEVEL) timer_printf(&T,"FindModulus");
    2616          95 :   if (!data) return bnrstark_cyclic(bnr, dtQ, prec);
    2617          95 :   return makescind(bnf_get_nf(bnf), AllStark(data, 0, newprec));
    2618             : }
    2619             : 
    2620             : /*******************************************************************/
    2621             : /*                                                                 */
    2622             : /*       Hilbert and Ray Class field using CM (Schertz)            */
    2623             : /*                                                                 */
    2624             : /*******************************************************************/
    2625             : /* form^2 = 1 ? */
    2626             : static int
    2627         581 : hasexp2(GEN form)
    2628             : {
    2629         581 :   GEN a = gel(form,1), b = gel(form,2), c = gel(form,3);
    2630         581 :   return !signe(b) || absequalii(a,b) || equalii(a,c);
    2631             : }
    2632             : static int
    2633         945 : uhasexp2(GEN form)
    2634             : {
    2635         945 :   long a = form[1], b = form[2], c = form[3];
    2636         945 :   return !b || a == labs(b) || a == c;
    2637             : }
    2638             : 
    2639             : GEN
    2640         325 : qfbforms(GEN D)
    2641             : {
    2642         325 :   ulong d = itou(D), dover3 = d/3, t, b2, a, b, c, h;
    2643         325 :   GEN L = cgetg((long)(sqrt((double)d) * log2(d)), t_VEC);
    2644         325 :   b2 = b = (d&1); h = 0;
    2645         325 :   if (!b) /* b = 0 treated separately to avoid special cases */
    2646             :   {
    2647         180 :     t = d >> 2; /* (b^2 - D) / 4*/
    2648        2110 :     for (a=1; a*a<=t; a++)
    2649        1930 :       if (c = t/a, t == c*a) gel(L,++h) = mkvecsmall3(a,0,c);
    2650         180 :     b = 2; b2 = 4;
    2651             :   }
    2652             :   /* now b > 0, b = D (mod 2) */
    2653        5770 :   for ( ; b2 <= dover3; b += 2, b2 = b*b)
    2654             :   {
    2655        5445 :     t = (b2 + d) >> 2; /* (b^2 - D) / 4*/
    2656             :     /* b = a */
    2657        5445 :     if (c = t/b, t == c*b) gel(L,++h) = mkvecsmall3(b,b,c);
    2658             :     /* b < a < c */
    2659     1365735 :     for (a = b+1; a*a < t; a++)
    2660     1360290 :       if (c = t/a, t == c*a)
    2661             :       {
    2662         755 :         gel(L,++h) = mkvecsmall3(a, b,c);
    2663         755 :         gel(L,++h) = mkvecsmall3(a,-b,c);
    2664             :       }
    2665             :     /* a = c */
    2666        5445 :     if (a * a == t) gel(L,++h) = mkvecsmall3(a,b,a);
    2667             :   }
    2668         325 :   setlg(L,h+1); return L;
    2669             : }
    2670             : 
    2671             : /* gcd(n, 24) */
    2672             : static long
    2673         581 : GCD24(long n)
    2674             : {
    2675         581 :   switch(n % 24)
    2676             :   {
    2677          25 :     case 0: return 24;
    2678          25 :     case 1: return 1;
    2679          20 :     case 2: return 2;
    2680           0 :     case 3: return 3;
    2681          85 :     case 4: return 4;
    2682           0 :     case 5: return 1;
    2683          75 :     case 6: return 6;
    2684           0 :     case 7: return 1;
    2685           0 :     case 8: return 8;
    2686           0 :     case 9: return 3;
    2687          65 :     case 10: return 2;
    2688           0 :     case 11: return 1;
    2689          85 :     case 12: return 12;
    2690           0 :     case 13: return 1;
    2691           0 :     case 14: return 2;
    2692           0 :     case 15: return 3;
    2693          65 :     case 16: return 8;
    2694           0 :     case 17: return 1;
    2695          66 :     case 18: return 6;
    2696           0 :     case 19: return 1;
    2697           0 :     case 20: return 4;
    2698           0 :     case 21: return 3;
    2699          70 :     case 22: return 2;
    2700           0 :     case 23: return 1;
    2701           0 :     default: return 0;
    2702             :   }
    2703             : }
    2704             : 
    2705             : struct gpq_data {
    2706             :   long p, q;
    2707             :   GEN sqd; /* sqrt(D), t_REAL */
    2708             :   GEN u, D;
    2709             :   GEN pq, pq2; /* p*q, 2*p*q */
    2710             :   GEN qfpq ; /* class of \P * \Q */
    2711             : };
    2712             : 
    2713             : /* find P and Q two non principal prime ideals (above p <= q) such that
    2714             :  *   cl(P) = cl(Q) if P,Q have order 2 in Cl(K).
    2715             :  *   Ensure that e = 24 / gcd(24, (p-1)(q-1)) = 1 */
    2716             : /* D t_INT, discriminant */
    2717             : static void
    2718          35 : init_pq(GEN D, struct gpq_data *T)
    2719             : {
    2720          35 :   const long Np = 6547; /* N.B. primepi(50000) = 5133 */
    2721          35 :   const ulong maxq = 50000;
    2722          35 :   GEN listp = cgetg(Np + 1, t_VECSMALL); /* primes p */
    2723          35 :   GEN listP = cgetg(Np + 1, t_VEC); /* primeform(p) if of order 2, else NULL */
    2724          35 :   GEN gcd24 = cgetg(Np + 1, t_VECSMALL); /* gcd(p-1, 24) */
    2725             :   forprime_t S;
    2726          35 :   long l = 1;
    2727          35 :   double best = 0.;
    2728             :   ulong q;
    2729             : 
    2730          35 :   u_forprime_init(&S, 2, ULONG_MAX);
    2731          35 :   T->D = D;
    2732          35 :   T->p = T->q = 0;
    2733             :   for(;;)
    2734        1271 :   {
    2735             :     GEN Q;
    2736             :     long i, gcdq, mod;
    2737             :     int order2, store;
    2738             :     double t;
    2739             : 
    2740        1306 :     q = u_forprime_next(&S);
    2741        1306 :     if (best > 0 && q >= maxq)
    2742             :     {
    2743           0 :       if (DEBUGLEVEL)
    2744           0 :         pari_warn(warner,"possibly suboptimal (p,q) for D = %Ps", D);
    2745           0 :       break;
    2746             :     }
    2747        1306 :     if (kroiu(D, q) < 0) continue; /* inert */
    2748         636 :     Q = qfi_red(primeform_u(D, q));
    2749         636 :     if (is_pm1(gel(Q,1))) continue; /* Q | q is principal */
    2750             : 
    2751         581 :     store = 1;
    2752         581 :     order2 = hasexp2(Q);
    2753         581 :     gcd24[l] = gcdq = GCD24(q-1);
    2754         581 :     mod = 24 / gcdq; /* mod must divide p-1 otherwise e > 1 */
    2755         581 :     listp[l] = q;
    2756         581 :     gel(listP,l) = order2 ? Q : NULL;
    2757         581 :     t = (q+1)/(double)(q-1);
    2758        1521 :     for (i = 1; i < l; i++) /* try all (p, q), p < q in listp */
    2759             :     {
    2760        1186 :       long p = listp[i], gcdp = gcd24[i];
    2761             :       double b;
    2762             :       /* P,Q order 2 => cl(Q) = cl(P) */
    2763        1186 :       if (order2 && gel(listP,i) && !gequal(gel(listP,i), Q)) continue;
    2764        1181 :       if (gcdp % gcdq == 0) store = 0; /* already a better one in the list */
    2765        1181 :       if ((p-1) % mod) continue;
    2766             : 
    2767         246 :       b = (t*(p+1)) / (p-1); /* (p+1)(q+1) / (p-1)(q-1) */
    2768         246 :       if (b > best) {
    2769          70 :         store = 0; /* (p,q) always better than (q,r) for r >= q */
    2770          70 :         best = b; T->q = q; T->p = p;
    2771          70 :         if (DEBUGLEVEL>2) err_printf("p,q = %ld,%ld\n", p, q);
    2772             :       }
    2773             :       /* won't improve with this q as largest member */
    2774         246 :       if (best > 0) break;
    2775             :     }
    2776             :     /* if !store or (q,r) won't improve on current best pair, forget that q */
    2777         581 :     if (store && t*t > best)
    2778          85 :       if (++l >= Np) pari_err_BUG("quadhilbert (not enough primes)");
    2779         581 :     if (!best) /* (p,q) with p < q always better than (q,q) */
    2780             :     { /* try (q,q) */
    2781         100 :       if (gcdq >= 12 && umodiu(D, q)) /* e = 1 and unramified */
    2782             :       {
    2783           5 :         double b = (t*q) / (q-1); /* q(q+1) / (q-1)^2 */
    2784           5 :         if (b > best) {
    2785           5 :           best = b; T->q = T->p = q;
    2786           5 :           if (DEBUGLEVEL>2) err_printf("p,q = %ld,%ld\n", q, q);
    2787             :         }
    2788             :       }
    2789             :     }
    2790             :     /* If (p1+1)(q+1) / (p1-1)(q-1) <= best, we can no longer improve
    2791             :      * even with best p : stop */
    2792         581 :     if ((listp[1]+1)*t <= (listp[1]-1)*best) break;
    2793             :   }
    2794          35 :   if (DEBUGLEVEL>1)
    2795           0 :     err_printf("(p, q) = %ld, %ld; gain = %f\n", T->p, T->q, 12*best);
    2796          35 : }
    2797             : 
    2798             : static GEN
    2799        2930 : gpq(GEN form, struct gpq_data *T)
    2800             : {
    2801        2930 :   pari_sp av = avma;
    2802        2930 :   long a = form[1], b = form[2], c = form[3], p = T->p, q = T->q;
    2803             :   GEN form2, w, z;
    2804        2930 :   int fl, real = 0;
    2805             : 
    2806        2930 :   form2 = qfbcomp_i(T->qfpq, mkqfb(stoi(a), stoi(-b), stoi(c), T->D));
    2807             :   /* form2 and form yield complex conjugate roots : only compute for the
    2808             :    * lexicographically smallest of the 2 */
    2809        2930 :   fl = cmpis(gel(form2,1), a);
    2810        2930 :   if (fl <= 0)
    2811             :   {
    2812        1540 :     if (fl < 0) return NULL;
    2813         150 :     fl = cmpis(gel(form2,2), b);
    2814         150 :     if (fl <= 0)
    2815             :     {
    2816         105 :       if (fl < 0) return NULL;
    2817             :       /* form == form2 : real root */
    2818          60 :       real = 1;
    2819             :     }
    2820             :   }
    2821             : 
    2822        1495 :   if (p == 2) { /* (a,b,c) = (1,1,0) mod 2 ? */
    2823         145 :     if (a % q == 0 && (a & b & 1) && !(c & 1))
    2824             :     { /* apply S : make sure that (a,b,c) represents odd values */
    2825           0 :       lswap(a,c); b = -b;
    2826             :     }
    2827             :   }
    2828        1495 :   if (a % p == 0 || a % q == 0)
    2829             :   { /* apply T^k, look for c' = a k^2 + b k + c coprime to N */
    2830         425 :     while (c % p == 0 || c % q == 0)
    2831             :     {
    2832          70 :       c += a + b;
    2833          70 :       b += a << 1;
    2834             :     }
    2835         355 :     lswap(a, c); b = -b; /* apply S */
    2836             :   }
    2837             :   /* now (a,b,c) ~ form and (a,pq) = 1 */
    2838             : 
    2839             :   /* gcd(2a, u) = 2,  w = u mod 2pq, -b mod 2a */
    2840        1495 :   w = Z_chinese(T->u, stoi(-b), T->pq2, utoipos(a << 1));
    2841        1495 :   z = double_eta_quotient(utoipos(a), w, T->D, T->p, T->q, T->pq, T->sqd);
    2842        1495 :   if (real && typ(z) == t_COMPLEX) z = gcopy(gel(z, 1));
    2843        1495 :   return gc_upto(av, z);
    2844             : }
    2845             : 
    2846             : /* returns an equation for the Hilbert class field of Q(sqrt(D)), D < 0
    2847             :  * fundamental discriminant */
    2848             : static GEN
    2849         330 : quadhilbertimag(GEN D)
    2850             : {
    2851             :   GEN L, P, Pi, Pr, qfp, u;
    2852             :   long h, i, prec;
    2853             :   struct gpq_data T;
    2854             :   pari_timer ti;
    2855             : 
    2856         330 :   if (DEBUGLEVEL>1) timer_start(&ti);
    2857         330 :   if (lgefint(D) == 3)
    2858         330 :     switch (D[2]) { /* = |D|; special cases where e > 1 */
    2859           5 :       case 3:
    2860             :       case 4:
    2861             :       case 7:
    2862             :       case 8:
    2863             :       case 11:
    2864             :       case 19:
    2865             :       case 43:
    2866             :       case 67:
    2867           5 :       case 163: return pol_x(0);
    2868             :     }
    2869         325 :   L = qfbforms(D);
    2870         325 :   h = lg(L)-1;
    2871         325 :   if (! (h & (h - 1))) /* power of 2 */
    2872             :   { /* check whether > |Cl|/2 elements have order <= 2 ==> 2-elementary */
    2873         295 :     long lim = (h>>1) + 1;
    2874        1235 :     for (i=1; i <= lim; i++)
    2875         945 :       if (!uhasexp2(gel(L,i))) break;
    2876         295 :     if (i > lim) return GenusFieldQuadImag(D);
    2877             :   }
    2878          35 :   if (DEBUGLEVEL>1) timer_printf(&ti,"class number = %ld",h);
    2879          35 :   init_pq(D, &T);
    2880          35 :   qfp = primeform_u(D, T.p);
    2881          35 :   T.pq =  muluu(T.p, T.q);
    2882          35 :   T.pq2 = shifti(T.pq,1);
    2883          35 :   if (T.p == T.q)
    2884             :   {
    2885           0 :     GEN qfbp2 = qfbcompraw(qfp, qfp);
    2886           0 :     u = gel(qfbp2,2);
    2887           0 :     T.u = modii(u, T.pq2);
    2888           0 :     T.qfpq = qfi_red(qfbp2);
    2889             :   }
    2890             :   else
    2891             :   {
    2892          35 :     GEN qfq = primeform_u(D, T.q), bp = gel(qfp,2), bq = gel(qfq,2);
    2893          35 :     T.u = Z_chinese(bp, bq, utoipos(T.p << 1), utoipos(T.q << 1));
    2894             :     /* T.u = bp (mod 2p), T.u = bq (mod 2q) */
    2895          35 :     T.qfpq = qfbcomp_i(qfp, qfq);
    2896             :   }
    2897             :   /* u modulo 2pq */
    2898          35 :   prec = LOWDEFAULTPREC;
    2899          35 :   Pr = cgetg(h+1,t_VEC);
    2900          35 :   Pi = cgetg(h+1,t_VEC);
    2901             :   for(;;)
    2902          10 :   {
    2903          45 :     long ex, exmax = 0, r1 = 0, r2 = 0;
    2904          45 :     pari_sp av0 = avma;
    2905          45 :     T.sqd = sqrtr_abs(itor(D, prec));
    2906        2975 :     for (i=1; i<=h; i++)
    2907             :     {
    2908        2930 :       GEN s = gpq(gel(L,i), &T);
    2909        2930 :       if (DEBUGLEVEL>3) err_printf("%ld ", i);
    2910        2930 :       if (!s) continue;
    2911        1495 :       if (typ(s) != t_COMPLEX) gel(Pr, ++r1) = s; /* real root */
    2912        1435 :       else                     gel(Pi, ++r2) = s;
    2913        1495 :       ex = gexpo(s); if (ex > 0) exmax += ex;
    2914             :     }
    2915          45 :     if (DEBUGLEVEL>1) timer_printf(&ti,"roots");
    2916          45 :     setlg(Pr, r1+1);
    2917          45 :     setlg(Pi, r2+1);
    2918          45 :     P = roots_to_pol_r1(shallowconcat(Pr,Pi), 0, r1);
    2919          45 :     P = grndtoi(P,&exmax);
    2920          45 :     if (DEBUGLEVEL>1) timer_printf(&ti,"product, error bits = %ld",exmax);
    2921          45 :     if (exmax <= -10) break;
    2922          10 :     set_avma(av0); prec += nbits2extraprec(DEFAULTPREC + exmax);
    2923          10 :     if (DEBUGLEVEL) pari_warn(warnprec,"quadhilbertimag",prec);
    2924             :   }
    2925          35 :   return P;
    2926             : }
    2927             : 
    2928             : GEN
    2929         412 : quadhilbert(GEN D, long prec)
    2930             : {
    2931         412 :   pari_sp av = avma;
    2932         412 :   GEN d = D;
    2933         412 :   quadray_init(&d, NULL, 0);
    2934         695 :   return gc_upto(av, signe(d)>0? quadhilbertreal(D,prec)
    2935         295 :                                     : quadhilbertimag(d));
    2936             : }
    2937             : 
    2938             : /* return a vector of all roots of 1 in bnf [not necessarily quadratic] */
    2939             : static GEN
    2940          50 : getallrootsof1(GEN bnf)
    2941             : {
    2942          50 :   GEN T, u, nf = bnf_get_nf(bnf), tu;
    2943          50 :   long i, n = bnf_get_tuN(bnf);
    2944             : 
    2945          50 :   if (n == 2) {
    2946          40 :     long N = nf_get_degree(nf);
    2947          40 :     return mkvec2(scalarcol_shallow(gen_m1, N),
    2948             :                   scalarcol_shallow(gen_1, N));
    2949             :   }
    2950          10 :   tu = poltobasis(nf, bnf_get_tuU(bnf));
    2951          10 :   T = zk_multable(nf, tu);
    2952          10 :   u = cgetg(n+1, t_VEC); gel(u,1) = tu;
    2953          40 :   for (i=2; i <= n; i++) gel(u,i) = ZM_ZC_mul(T, gel(u,i-1));
    2954          10 :   return u;
    2955             : }
    2956             : static int
    2957          90 : _ZV_cmp(void *E, GEN x, GEN y) { (void)E; return ZV_cmp(x,y); }
    2958             : 
    2959             : /* assume bnr has the right conductor. Find la != root of 1 mod F such that
    2960             :  * Norm la = 1 mod 2*(F \cap Z) */
    2961             : static GEN
    2962          50 : get_lambda(GEN bnr)
    2963             : {
    2964          50 :   GEN bnf = bnr_get_bnf(bnr), nf = bnf_get_nf(bnf), T = nf_get_pol(nf);
    2965          50 :   GEN U, la, F = gel(bnr_get_mod(bnr), 1);
    2966          50 :   long a, b, i, lu, f = itos(gcoeff(F,1,1)), f2 = 2*f, v = varn(T);
    2967             : 
    2968          50 :   U = getallrootsof1(bnf); lu = lg(U);
    2969         170 :   for (i=1; i<lu; i++) gel(U,i) = ZC_modhnf(gel(U,i), F, NULL);
    2970          50 :   U = gen_sort_uniq(U, NULL, &_ZV_cmp); /* set of roots of 1 mod F */
    2971         200 :   for (b = 3; b < f; b += 2)
    2972             :   { /* a = 0: shortcut */
    2973         150 :     if (Fl_sqr(b, f2) != 1) continue;
    2974           0 :     la = scalarcol_shallow(stoi(b), 2); /* unit mod F ? */
    2975           0 :     if (!tablesearch(U, la, &ZV_cmp)) return la;
    2976             :   }
    2977          70 :   for (a = 1; a < f2; a++)
    2978         940 :     for (b = 0; b < f2; b++)
    2979             :     {
    2980         920 :       if (!odd(a) && !odd(b)) continue;
    2981         840 :       la = deg1pol_shallow(stoi(a), stoi(b), v); /* ax + b */
    2982         840 :       if (umodiu(ZX_resultant(la,T), f2) != 1) continue;
    2983          60 :       la = poltobasis(nf, la); /* unit mod F ? */
    2984          60 :       if (!tablesearch(U, ZC_modhnf(la, F, NULL), &ZV_cmp)) return la;
    2985             :     }
    2986           0 :   pari_err_BUG("get_lambda");
    2987             :   return NULL;/*LCOV_EXCL_LINE*/
    2988             : }
    2989             : 
    2990             : static GEN
    2991        6270 : to_approx(GEN nf, GEN a)
    2992             : {
    2993        6270 :   GEN M = nf_get_M(nf);
    2994        6270 :   return gadd(gel(a,1), gmul(gcoeff(M,1,2),gel(a,2)));
    2995             : }
    2996             : /* Z-basis for a (over C) */
    2997             : static GEN
    2998        3110 : get_om(GEN nf, GEN a) {
    2999        3110 :   return mkvec2(to_approx(nf,gel(a,2)),
    3000        3110 :                 to_approx(nf,gel(a,1)));
    3001             : }
    3002             : 
    3003             : /* Compute all elts in class group G = [|G|,c,g], c=cyclic factors, g=gens.
    3004             :  * Set list[j + 1] = g1^e1...gk^ek where j is the integer
    3005             :  *   ek + ck [ e(k-1) + c(k-1) [... + c2 [e1]]...] */
    3006             : static GEN
    3007          50 : getallelts(GEN bnr)
    3008             : {
    3009             :   GEN nf, C, c, g, list, pows, gk;
    3010             :   long lc, i, j, no;
    3011             : 
    3012          50 :   nf = bnr_get_nf(bnr);
    3013          50 :   no = itos( bnr_get_no(bnr) );
    3014          50 :   c = bnr_get_cyc(bnr);
    3015          50 :   g = bnr_get_gen_nocheck(bnr); lc = lg(c)-1;
    3016          50 :   list = cgetg(no+1,t_VEC);
    3017          50 :   gel(list,1) = matid(nf_get_degree(nf)); /* (1) */
    3018          50 :   if (!no) return list;
    3019             : 
    3020          50 :   pows = cgetg(lc+1,t_VEC);
    3021          50 :   c = leafcopy(c); settyp(c, t_VECSMALL);
    3022         100 :   for (i=1; i<=lc; i++)
    3023             :   {
    3024          50 :     long k = itos(gel(c,i));
    3025          50 :     c[i] = k;
    3026          50 :     gk = cgetg(k, t_VEC); gel(gk,1) = gel(g,i);
    3027        3060 :     for (j=2; j<k; j++)
    3028        3010 :       gel(gk,j) = idealmoddivisor(bnr, idealmul(nf, gel(gk,j-1), gel(gk,1)));
    3029          50 :     gel(pows,i) = gk; /* powers of g[i] */
    3030             :   }
    3031             : 
    3032          50 :   C = cgetg(lc+1, t_VECSMALL); C[1] = c[lc];
    3033          50 :   for (i=2; i<=lc; i++) C[i] = C[i-1] * c[lc-i+1];
    3034             :   /* C[i] = c(k-i+1) * ... * ck */
    3035             :   /* j < C[i+1] <==> j only involves g(k-i)...gk */
    3036          50 :   i = 1;
    3037        3110 :   for (j=1; j < C[1]; j++)
    3038        3060 :     gel(list, j+1) = gmael(pows,lc,j);
    3039          50 :   while(j<no)
    3040             :   {
    3041             :     long k;
    3042             :     GEN a;
    3043           0 :     if (j == C[i+1]) i++;
    3044           0 :     a = gmael(pows,lc-i,j/C[i]);
    3045           0 :     k = j%C[i] + 1;
    3046           0 :     if (k > 1) a = idealmoddivisor(bnr, idealmul(nf, a, gel(list,k)));
    3047           0 :     gel(list, ++j) = a;
    3048             :   }
    3049          50 :   return list;
    3050             : }
    3051             : 
    3052             : /* x quadratic integer (approximate), recognize it. If error return NULL */
    3053             : static GEN
    3054        3160 : findbezk(GEN nf, GEN x)
    3055             : {
    3056        3160 :   GEN a,b, M = nf_get_M(nf), u = gcoeff(M,1,2);
    3057             :   long ea, eb;
    3058             : 
    3059             :   /* u t_COMPLEX generator of nf.zk, write x ~ a + b u, a,b in Z */
    3060        3160 :   b = grndtoi(mpdiv(imag_i(x), gel(u,2)), &eb);
    3061        3160 :   if (eb > -20) return NULL;
    3062        3160 :   a = grndtoi(mpsub(real_i(x), mpmul(b,gel(u,1))), &ea);
    3063        3160 :   if (ea > -20) return NULL;
    3064        3160 :   return nf_to_scalar_or_polmod(nf, mkcol2(a,b));
    3065             : }
    3066             : 
    3067             : static GEN
    3068          50 : findbezk_pol(GEN nf, GEN x)
    3069             : {
    3070          50 :   long i, lx = lg(x);
    3071          50 :   GEN y = cgetg(lx,t_POL);
    3072        3210 :   for (i=2; i<lx; i++)
    3073        3160 :     if (! (gel(y,i) = findbezk(nf,gel(x,i))) ) return NULL;
    3074          50 :   y[1] = x[1]; return y;
    3075             : }
    3076             : 
    3077             : /* P approximation computed at initial precision prec. Compute needed prec
    3078             :  * to know P with 1 word worth of trailing decimals */
    3079             : static long
    3080           0 : get_prec(GEN P, long prec)
    3081             : {
    3082           0 :   long k = gprecision(P);
    3083           0 :   if (k == 3) return precdbl(prec); /* approximation not trustworthy */
    3084           0 :   k = prec - k; /* lost precision when computing P */
    3085           0 :   if (k < 0) k = 0;
    3086           0 :   k += nbits2prec(gexpo(P) + 128);
    3087           0 :   if (k <= prec) k = precdbl(prec); /* dubious: old prec should have worked */
    3088           0 :   return k;
    3089             : }
    3090             : 
    3091             : /* Compute data for ellphist */
    3092             : static GEN
    3093        3110 : ellphistinit(GEN om, long prec)
    3094             : {
    3095        3110 :   GEN res,om1b,om2b, om1 = gel(om,1), om2 = gel(om,2);
    3096             : 
    3097        3110 :   if (gsigne(mulimag(om1,conj_i(om2))) < 0) { swap(om1,om2); om = mkvec2(om1,om2); }
    3098        3110 :   om1b = conj_i(om1);
    3099        3110 :   om2b = conj_i(om2); res = cgetg(4,t_VEC);
    3100        3110 :   gel(res,1) = gdivgu(elleisnum(om,2,prec),12);
    3101        3110 :   gel(res,2) = gdiv(PiI2(prec), gmul(om2, mulimag(om1b,om2)));
    3102        3110 :   gel(res,3) = om2b; return res;
    3103             : }
    3104             : 
    3105             : /* Computes log(phi^*(z,om)), using res computed by ellphistinit */
    3106             : static GEN
    3107        6220 : ellphist(GEN om, GEN res, GEN z, long prec)
    3108             : {
    3109        6220 :   GEN u = mulimag(z, gel(res,3));
    3110        6220 :   GEN zst = gsub(gmul(u, gel(res,2)), gmul(z,gel(res,1)));
    3111        6220 :   return gsub(ellsigma(om,z,1,prec),gmul2n(gmul(z,zst),-1));
    3112             : }
    3113             : 
    3114             : /* Computes phi^*(la,om)/phi^*(1,om) where (1,om) is an oriented basis of the
    3115             :    ideal gf*gc^{-1} */
    3116             : static GEN
    3117        3110 : computeth2(GEN om, GEN la, long prec)
    3118             : {
    3119        3110 :   GEN p1,p2,res = ellphistinit(om,prec);
    3120             : 
    3121        3110 :   p1 = gsub(ellphist(om,res,la,prec), ellphist(om,res,gen_1,prec));
    3122        3110 :   p2 = imag_i(p1);
    3123        3110 :   if (gexpo(real_i(p1)) > 20 || gexpo(p2) > minss(prec,realprec(p2)) - 10)
    3124           0 :     return NULL;
    3125        3110 :   return gexp(p1,prec);
    3126             : }
    3127             : 
    3128             : /* Computes P_2(X)=polynomial in Z_K[X] closest to prod_gc(X-th2(gc)) where
    3129             :    the product is over the ray class group bnr.*/
    3130             : static GEN
    3131          50 : computeP2(GEN bnr, long prec)
    3132             : {
    3133          50 :   long clrayno, i, first = 1;
    3134          50 :   pari_sp av=avma, av2;
    3135          50 :   GEN listray, P0, P, lanum, la = get_lambda(bnr);
    3136          50 :   GEN nf = bnr_get_nf(bnr), f = gel(bnr_get_mod(bnr), 1);
    3137          50 :   if (DEBUGLEVEL) err_printf("lambda = %Ps\n",la);
    3138          50 :   listray = getallelts(bnr);
    3139          50 :   clrayno = lg(listray)-1; av2 = avma;
    3140          50 : PRECPB:
    3141          50 :   if (!first)
    3142             :   {
    3143           0 :     if (DEBUGLEVEL) pari_warn(warnprec,"computeP2",prec);
    3144           0 :     nf = gc_GEN(av2, nfnewprec_shallow(bnr_get_nf(bnr),prec));
    3145             :   }
    3146          50 :   first = 0; lanum = to_approx(nf,la);
    3147          50 :   P = cgetg(clrayno+1,t_VEC);
    3148        3160 :   for (i=1; i<=clrayno; i++)
    3149             :   {
    3150        3110 :     GEN om = get_om(nf, idealdiv(nf,f,gel(listray,i)));
    3151        3110 :     GEN s = computeth2(om,lanum,prec);
    3152        3110 :     if (!s) { prec = precdbl(prec); goto PRECPB; }
    3153        3110 :     gel(P,i) = s;
    3154             :   }
    3155          50 :   P0 = roots_to_pol(P, 0);
    3156          50 :   P = findbezk_pol(nf, P0);
    3157          50 :   if (!P) { prec = get_prec(P0, prec); goto PRECPB; }
    3158          50 :   return gc_GEN(av, P);
    3159             : }
    3160             : 
    3161             : #define nexta(a) (a>0 ? -a : 1-a)
    3162             : static GEN
    3163          35 : do_compo(GEN A0, GEN B)
    3164             : {
    3165          35 :   long a, i, l = lg(B), v = fetch_var_higher();
    3166             :   GEN A, z;
    3167             :   /* now v > x = pol_x(0) > nf variable */
    3168          35 :   B = leafcopy(B); setvarn(B, v);
    3169         150 :   for (i = 2; i < l; i++) gel(B,i) = monomial(gel(B,i), l-i-1, 0);
    3170             :   /* B := x^deg(B) B(v/x) */
    3171          35 :   A = A0 = leafcopy(A0); setvarn(A0, v);
    3172          40 :   for  (a = 0;; a = nexta(a))
    3173             :   {
    3174          40 :     if (a) A = RgX_Rg_translate(A0, stoi(a));
    3175          40 :     z = resultant(A,B); /* in variable 0 */
    3176          40 :     if (issquarefree(z)) break;
    3177             :   }
    3178          35 :   (void)delete_var(); return z;
    3179             : }
    3180             : #undef nexta
    3181             : 
    3182             : static GEN
    3183          10 : galoisapplypol(GEN nf, GEN s, GEN x)
    3184             : {
    3185          10 :   long i, lx = lg(x);
    3186          10 :   GEN y = cgetg(lx,t_POL);
    3187             : 
    3188          40 :   for (i=2; i<lx; i++) gel(y,i) = galoisapply(nf,s,gel(x,i));
    3189          10 :   y[1] = x[1]; return y;
    3190             : }
    3191             : /* x quadratic, write it as ua + v, u,v rational */
    3192             : static GEN
    3193          50 : findquad(GEN a, GEN x, GEN T)
    3194             : {
    3195          50 :   pari_sp av = avma;
    3196             :   GEN u, v;
    3197          50 :   if (typ(x) == t_POLMOD) x = gel(x,2);
    3198          50 :   if (typ(a) == t_POLMOD) a = gel(a,2);
    3199          50 :   u = poldivrem(x, a, &v);
    3200          50 :   u = simplify_shallow(u);
    3201          50 :   v = simplify_shallow(v);
    3202          50 :   if (!is_scalar_t(typ(u))) pari_err_TYPE("findquad", u);
    3203          50 :   if (!is_scalar_t(typ(v))) pari_err_TYPE("findquad", v);
    3204          50 :   if (isintzero(u)) return gc_GEN(av, v);
    3205          20 :   return gc_upto(av, gmodulo(deg1pol_shallow(u, v, varn(a)), T));
    3206             : }
    3207             : static GEN
    3208          10 : RgX_findquad(GEN x, GEN a, GEN T)
    3209          60 : { pari_APPLY_pol_normalized(findquad(a, gel(x,i), T)); }
    3210             : /* m is 3, 4 or 12 */
    3211             : static GEN
    3212          25 : compocyclo(GEN D, long m)
    3213          25 : { return do_compo(quadhilbertimag(D), polcyclo(m,0)); }
    3214             : /* m is prime or 4 * prime */
    3215             : static GEN
    3216          10 : compocyclop(GEN nf, long m)
    3217             : {
    3218          10 :   GEN sb, a, b, s, p1, p2, R, res, polL, polLK, nfL, D = nf_get_disc(nf);
    3219             :   long ell, vx;
    3220             : 
    3221          10 :   p1 = quadhilbertimag(D);
    3222          10 :   p2 = polcyclo(m,0);
    3223          10 :   ell = odd(m)? m: (m>>2); /* prime */
    3224          10 :   if (absequalui(ell,D)) /* ell = |D| */
    3225             :   {
    3226           0 :     p2 = gcoeff(nffactor(nf,p2),1,1);
    3227           0 :     return do_compo(p1,p2);
    3228             :   }
    3229          10 :   if (ell%4 == 3) ell = -ell;
    3230             :   /* nf = K = Q(a), L = K(b) quadratic extension = Q(t) */
    3231          10 :   polLK = quadpoly_i(stoi(ell)); /* relative polynomial */
    3232          10 :   res = rnfequation2(nf, polLK);
    3233          10 :   vx = nf_get_varn(nf);
    3234          10 :   polL = gsubst(gel(res,1),0,pol_x(vx)); /* = charpoly(t) */
    3235          10 :   a = gsubst(lift_shallow(gel(res,2)), 0,pol_x(vx));
    3236          10 :   b = gsub(pol_x(vx), gmul(gel(res,3), a));
    3237          10 :   nfL = nfinit(polL, DEFAULTPREC);
    3238          10 :   p1 = gcoeff(nffactor(nfL,p1),1,1);
    3239          10 :   p2 = gcoeff(nffactor(nfL,p2),1,1);
    3240          10 :   R = do_compo(p1,p2); /* relative equation over L */
    3241             :   /* compute non trivial s in Gal(L / K) */
    3242          10 :   sb= gneg(gadd(b, RgX_coeff(polLK,1))); /* s(b) = Tr(b) - b */
    3243          10 :   s = gadd(pol_x(vx), gsub(sb, b)); /* s(t) = t + s(b) - b */
    3244          10 :   R = gmul(R, galoisapplypol(nfL, s, R));
    3245          10 :   return RgX_findquad(R, a, nf_get_pol(nf));
    3246             : }
    3247             : 
    3248             : /* I integral ideal in HNF. (x) = I, x small in Z ? */
    3249             : static long
    3250          85 : isZ(GEN I)
    3251             : {
    3252          85 :   GEN x = gcoeff(I,1,1);
    3253          85 :   if (signe(gcoeff(I,1,2)) || !equalii(x, gcoeff(I,2,2))) return 0;
    3254          75 :   return is_bigint(x)? -1: itos(x);
    3255             : }
    3256             : 
    3257             : /* Treat special cases directly. return NULL if not special case */
    3258             : static GEN
    3259          85 : treatspecialsigma(GEN bnr)
    3260             : {
    3261          85 :   GEN bnf = bnr_get_bnf(bnr), nf = bnf_get_nf(bnf);
    3262          85 :   GEN f = gel(bnr_get_mod(bnr), 1),  D = nf_get_disc(nf);
    3263             :   GEN p1, p2;
    3264          85 :   long Ds, fl, tryf, i = isZ(f);
    3265             : 
    3266          85 :   if (i == 1) return quadhilbertimag(D); /* f = 1 */
    3267             : 
    3268          85 :   if (absequaliu(D,3)) /* Q(j) */
    3269             :   {
    3270           0 :     if (i == 4 || i == 5 || i == 7) return polcyclo(i,0);
    3271           0 :     if (!absequaliu(gcoeff(f,1,1),9) || !absequaliu(Z_content(f),3)) return NULL;
    3272             :     /* f = P_3^3 */
    3273           0 :     p1 = mkpolmod(bnf_get_tuU(bnf), nf_get_pol(nf));
    3274           0 :     return gadd(pol_xn(3,0), p1); /* x^3+j */
    3275             :   }
    3276          85 :   if (absequaliu(D,4)) /* Q(i) */
    3277             :   {
    3278          10 :     if (i == 3 || i == 5) return polcyclo(i,0);
    3279          10 :     if (i != 4) return NULL;
    3280           0 :     p1 = mkpolmod(bnf_get_tuU(bnf), nf_get_pol(nf));
    3281           0 :     return gadd(pol_xn(2,0), p1); /* x^2+i */
    3282             :   }
    3283          75 :   Ds = smodis(D,48);
    3284          75 :   if (i)
    3285             :   {
    3286          65 :     if (i==2 && Ds%16== 8) return compocyclo(D, 4);
    3287          60 :     if (i==3 && Ds% 3== 1) return compocyclo(D, 3);
    3288          50 :     if (i==4 && Ds% 8== 1) return compocyclo(D, 4);
    3289          45 :     if (i==6 && Ds   ==40) return compocyclo(D,12);
    3290          40 :     return NULL;
    3291             :   }
    3292             : 
    3293          10 :   p1 = gcoeff(f,1,1); /* integer > 0 */
    3294          10 :   tryf = itou_or_0(p1); if (!tryf) return NULL;
    3295          10 :   p2 = gcoeff(f,2,2); /* integer > 0 */
    3296          10 :   if (is_pm1(p2)) fl = 0;
    3297             :   else {
    3298           0 :     if (Ds % 16 != 8 || !absequaliu(Z_content(f),2)) return NULL;
    3299           0 :     fl = 1; tryf >>= 1;
    3300             :   }
    3301          10 :   if (tryf <= 3 || umodiu(D, tryf) || !uisprime(tryf)) return NULL;
    3302          10 :   if (fl) tryf <<= 2;
    3303          10 :   return compocyclop(nf, tryf);
    3304             : }
    3305             : 
    3306             : GEN
    3307         117 : quadray(GEN D, GEN f, long prec)
    3308             : {
    3309         117 :   GEN bnr, y, bnf, H = NULL;
    3310         117 :   pari_sp av = avma;
    3311             : 
    3312         117 :   if (isint1(f)) return quadhilbert(D, prec);
    3313          90 :   if (typ(D) == t_INT && typ(f) != t_INT)
    3314           0 :     pari_err_TYPE("quadray [conductor]", f);
    3315          90 :   quadray_init(&D, &bnf, prec);
    3316          90 :   bnr = Buchray(bnf, f, nf_INIT|nf_GEN);
    3317          90 :   if (is_pm1(bnr_get_no(bnr))) { set_avma(av); return pol_x(0); }
    3318          90 :   if (signe(D) > 0)
    3319           5 :     y = bnrstark(bnr, H, prec);
    3320             :   else
    3321             :   {
    3322          85 :     bnr_subgroup_sanitize(&bnr, &H);
    3323          85 :     y = treatspecialsigma(bnr);
    3324          85 :     if (!y) y = computeP2(bnr, prec);
    3325             :   }
    3326          90 :   return gc_upto(av, y);
    3327             : }

Generated by: LCOV version 1.16