Code coverage tests

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

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

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

LCOV - code coverage report
Current view: top level - basemath - buch1.c (source / functions) Hit Total Coverage
Test: PARI/GP v2.18.0 lcov report (development 29732-95c6201d93) Lines: 679 711 95.5 %
Date: 2024-11-21 09:08:54 Functions: 51 55 92.7 %
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             : #include "pari.h"
      15             : #include "paripriv.h"
      16             : #include "mpqs.h"
      17             : 
      18             : #define DEBUGLEVEL DEBUGLEVEL_quadclassunit
      19             : 
      20             : /*******************************************************************/
      21             : /*                                                                 */
      22             : /*         CLASS GROUP AND REGULATOR (McCURLEY, BUCHMANN)          */
      23             : /*                   QUADRATIC FIELDS                              */
      24             : /*                                                                 */
      25             : /*******************************************************************/
      26             : /* For largeprime() hashtable. Note that hashed pseudoprimes are odd (unless
      27             :  * 2 | index), hence the low order bit is not useful. So we hash
      28             :  * HASHBITS bits starting at bit 1, not bit 0 */
      29             : #define HASHBITS 10
      30             : static const long HASHT = 1L << HASHBITS;
      31             : 
      32             : static long
      33     2565259 : hash(long q) { return (q & ((1L << (HASHBITS+1)) - 1)) >> 1; }
      34             : #undef HASHBITS
      35             : 
      36             : /* See buch2.c:
      37             :  * B->subFB contains split p such that \prod p > sqrt(B->Disc)
      38             :  * B->powsubFB contains powers of forms in B->subFB */
      39             : #define RANDOM_BITS 4
      40             : static const long CBUCH = (1L<<RANDOM_BITS)-1;
      41             : 
      42             : struct buch_quad
      43             : {
      44             :   ulong limhash;
      45             :   long KC, KC2, PRECREG;
      46             :   long *primfact, *exprimfact, **hashtab;
      47             :   GEN FB, numFB, prodFB;
      48             :   GEN powsubFB, vperm, subFB, badprim;
      49             :   struct qfr_data *q;
      50             : };
      51             : 
      52             : /*******************************************************************/
      53             : /*                                                                 */
      54             : /*  Routines related to binary quadratic forms (for internal use)  */
      55             : /*                                                                 */
      56             : /*******************************************************************/
      57             : /* output canonical representative wrt projection Cl^+ --> Cl (a > 0) */
      58             : static GEN
      59     1166963 : qfr3_canon(GEN x, struct qfr_data *S)
      60             : {
      61     1166963 :   GEN a = gel(x,1), c = gel(x,3);
      62     1166963 :   if (signe(a) < 0) {
      63      403298 :     if (absequalii(a,c)) return qfr3_rho(x, S);
      64      403284 :     setsigne(a, 1);
      65      403284 :     setsigne(c,-1);
      66             :   }
      67     1166949 :   return x;
      68             : }
      69             : static GEN
      70        3710 : qfr3_canon_safe(GEN x, struct qfr_data *S)
      71             : {
      72        3710 :   GEN a = gel(x,1), c = gel(x,3);
      73        3710 :   if (signe(a) < 0) {
      74         224 :     if (absequalii(a,c)) return qfr3_rho(x, S);
      75         224 :     gel(x,1) = negi(a);
      76         224 :     gel(x,3) = negi(c);
      77             :   }
      78        3710 :   return x;
      79             : }
      80             : static GEN
      81     1952748 : qfr5_canon(GEN x, struct qfr_data *S)
      82             : {
      83     1952748 :   GEN a = gel(x,1), c = gel(x,3);
      84     1952748 :   if (signe(a) < 0) {
      85      723268 :     if (absequalii(a,c)) return qfr5_rho(x, S);
      86      723254 :     setsigne(a, 1);
      87      723254 :     setsigne(c,-1);
      88             :   }
      89     1952734 :   return x;
      90             : }
      91             : static GEN
      92     1735160 : QFR5_comp(GEN x,GEN y, struct qfr_data *S)
      93     1735160 : { return qfr5_canon(qfr5_comp(x,y,S), S); }
      94             : static GEN
      95     1007545 : QFR3_comp(GEN x, GEN y, struct qfr_data *S)
      96     1007545 : { return qfr3_canon(qfr3_comp(x,y,S), S); }
      97             : 
      98             : /* compute rho^n(x) */
      99             : static GEN
     100      220157 : qfr5_rho_pow(GEN x, long n, struct qfr_data *S)
     101             : {
     102             :   long i;
     103      220157 :   pari_sp av = avma;
     104     2194360 :   for (i=1; i<=n; i++)
     105             :   {
     106     1974203 :     x = qfr5_rho(x,S);
     107     1974203 :     if (gc_needed(av,1))
     108             :     {
     109           0 :       if(DEBUGMEM>1) pari_warn(warnmem,"qfr5_rho_pow");
     110           0 :       x = gerepilecopy(av, x);
     111             :     }
     112             :   }
     113      220157 :   return gerepilecopy(av, x);
     114             : }
     115             : 
     116             : static GEN
     117      217588 : qfr5_pf(struct qfr_data *S, long p, long prec)
     118             : {
     119      217588 :   GEN y = primeform_u(S->D,p);
     120      217588 :   return qfr5_canon(qfr5_red(qfr_to_qfr5(y,prec), S), S);
     121             : }
     122             : 
     123             : static GEN
     124      159390 : qfr3_pf(struct qfr_data *S, long p)
     125             : {
     126      159390 :   GEN y = primeform_u(S->D,p);
     127      159390 :   return qfr3_canon(qfr3_red(y, S), S);
     128             : }
     129             : 
     130             : #define qfi_pf primeform_u
     131             : 
     132             : /* Warning: ex[0] not set in general */
     133             : static GEN
     134     4061435 : init_form(struct buch_quad *B, GEN ex,
     135             :           GEN (*comp)(GEN,GEN,struct qfr_data *S))
     136             : {
     137     4061435 :   long i, l = lg(B->powsubFB);
     138     4061435 :   GEN F = NULL;
     139    22906817 :   for (i=1; i<l; i++)
     140    18856493 :     if (ex[i])
     141             :     {
     142    17679988 :       GEN t = gmael(B->powsubFB,i,ex[i]);
     143    17679988 :       F = F? comp(F, t, B->q): t;
     144             :     }
     145     4050324 :   return F;
     146             : }
     147             : static GEN
     148      197015 : qfr5_factorback(struct buch_quad *B, GEN ex) { return init_form(B, ex, &QFR5_comp); }
     149             : 
     150             : static GEN
     151    11739092 : QFI_comp(GEN x, GEN y, struct qfr_data *S) { (void)S; return qfbcomp_i(x,y); }
     152             : 
     153             : static GEN
     154      160618 : qfi_factorback(struct buch_quad *B, GEN ex) { return init_form(B, ex, &QFI_comp); }
     155             : 
     156             : static GEN
     157     3703409 : random_form(struct buch_quad *B, GEN ex,
     158             :             GEN (*comp)(GEN,GEN, struct qfr_data *S))
     159             : {
     160     3703409 :   long i, l = lg(ex);
     161     3703409 :   pari_sp av = avma;
     162             :   GEN F;
     163             :   for(;;)
     164             :   {
     165    20591116 :     for (i=1; i<l; i++) ex[i] = random_bits(RANDOM_BITS);
     166     3703695 :     if ((F = init_form(B, ex, comp))) return F;
     167         459 :     set_avma(av);
     168             :   }
     169             : }
     170             : static GEN
     171      161567 : qfr3_random(struct buch_quad *B,GEN ex){ return random_form(B, ex, &QFR3_comp); }
     172             : static GEN
     173     3541868 : qfi_random(struct buch_quad *B,GEN ex) { return random_form(B, ex, &QFI_comp); }
     174             : 
     175             : /*******************************************************************/
     176             : /*                                                                 */
     177             : /*                     Common subroutines                          */
     178             : /*                                                                 */
     179             : /*******************************************************************/
     180             : /* We assume prime ideals of norm < D generate Cl(K); failure with
     181             :  * a factor base of primes with norm < LIMC <= D. Suggest new value.
     182             :  * All values of the form c * log^2 (disc K) [where D has c = 4
     183             :  * (Grenie-Molteni, under GRH)]. A value of c in [0.3, 1] should be OK for most
     184             :  * fields. No example is known where c >= 2 is necessary. */
     185             : ulong
     186        2392 : bnf_increase_LIMC(ulong LIMC, ulong D)
     187             : {
     188        2392 :   if (LIMC >= D) pari_err_BUG("Buchmann's algorithm");
     189        2392 :   if (LIMC <= D / 13.333)
     190         236 :     LIMC *= 2; /* tiny c <= 0.3 : double it */
     191             :   else
     192        2156 :     LIMC += maxuu(1, D / 20); /* large c, add 0.2 to it */
     193        2392 :   if (LIMC > D) LIMC = D;
     194        2392 :   return LIMC;
     195             : }
     196             : 
     197             : /* Is |q| <= p ? */
     198             : static int
     199       10416 : isless_iu(GEN q, ulong p) {
     200       10416 :   long l = lgefint(q);
     201       10416 :   return l==2 || (l == 3 && uel(q,2) <= p);
     202             : }
     203             : 
     204             : static GEN
     205     5029157 : Z_isquasismooth_prod(GEN N, GEN P)
     206             : {
     207     5029157 :   P = gcdii(P,N);
     208    10032952 :   while (!is_pm1(P))
     209             :   {
     210     4998126 :     N = diviiexact(N, P);
     211     5006155 :     P = gcdii(N, P);
     212             :   }
     213     5024670 :   return N;
     214             : }
     215             : 
     216             : static long
     217     5035836 : factorquad(struct buch_quad *B, GEN f, long nFB, ulong limp)
     218             : {
     219             :   ulong X;
     220     5035836 :   long i, lo = 0;
     221     5035836 :   GEN F, x = gel(f,1), FB = B->FB, P = B->primfact, E = B->exprimfact;
     222     5035836 :   if (B->badprim && !is_pm1(gcdii(x, B->badprim))) return 0;
     223     5028899 :   F =  Z_isquasismooth_prod(x, B->prodFB);
     224     5024531 :   if (cmpiu(F, B->limhash) > 0) return 0;
     225     4400848 :   for (i=1; lgefint(x) > 3; i++)
     226             :   {
     227       10416 :     ulong p = uel(FB,i), r;
     228       10416 :     GEN q = absdiviu_rem(x, p, &r);
     229       10416 :     if (!r)
     230             :     {
     231        1530 :       long k = 0;
     232        2633 :       do { k++; x = q; q = absdiviu_rem(x, p, &r); } while (!r);
     233        1530 :       lo++; P[lo] = p; E[lo] = k;
     234             :     }
     235       10416 :     if (isless_iu(q,p)) {
     236           1 :       if (lgefint(x) == 3) { X = uel(x,2); goto END; }
     237          34 :       return 0;
     238             :     }
     239       10415 :     if (i == nFB) return 0;
     240             :   }
     241     4390432 :   X = uel(x,2);
     242     4390432 :   if (X == 1) { P[0] = 0; return 1; }
     243    68604574 :   for (;; i++)
     244    68604574 :   { /* single precision affair, split for efficiency */
     245    72976338 :     ulong p = uel(FB,i);
     246    72976338 :     ulong q = X / p, r = X % p; /* gcc makes a single div */
     247    72976338 :     if (!r)
     248             :     {
     249     5620521 :       long k = 0;
     250     6787854 :       do { k++; X = q; q = X / p; r = X % p; } while (!r);
     251     5620521 :       lo++; P[lo] = p; E[lo] = k;
     252             :     }
     253    72976338 :     if (q <= p) break;
     254    68637539 :     if (i == nFB) return 0;
     255             :   }
     256     4338800 : END:
     257     4338800 :   if (X > B->limhash) return 0;
     258     4338800 :   if (X != 1 && X <= limp) { lo++; P[lo] = X; E[lo] = 1; X = 1; }
     259     4338800 :   P[0] = lo; return X;
     260             : }
     261             : 
     262             : /* Check for a "large prime relation" involving q; q may not be prime */
     263             : static long *
     264     2565221 : largeprime(struct buch_quad *B, long q, GEN ex, long np, long nrho)
     265             : {
     266     2565221 :   const long hashv = hash(q);
     267     2565261 :   long *pt, i, l = lg(B->subFB);
     268             : 
     269     2786998 :   for (pt = B->hashtab[hashv]; ; pt = (long*) pt[0])
     270             :   {
     271     2786998 :     if (!pt)
     272             :     {
     273     2365893 :       pt = (long*) pari_malloc((l+3) * sizeof(long));
     274     2366091 :       *pt++ = nrho; /* nrho = pt[-3] */
     275     2366091 :       *pt++ = np;   /* np   = pt[-2] */
     276     2366091 :       *pt++ = q;    /* q    = pt[-1] */
     277     2366091 :       pt[0] = (long)B->hashtab[hashv];
     278    14987910 :       for (i=1; i<l; i++) pt[i]=ex[i];
     279     2366091 :       B->hashtab[hashv]=pt; return NULL;
     280             :     }
     281      421105 :     if (pt[-1] == q) break;
     282             :   }
     283      237900 :   for(i=1; i<l; i++)
     284      234610 :     if (pt[i] != ex[i]) return pt;
     285        3290 :   return (pt[-2]==np)? NULL: pt;
     286             : }
     287             : 
     288             : static void
     289      169724 : clearhash(long **hash)
     290             : {
     291             :   long *pt;
     292             :   long i;
     293   173939863 :   for (i=0; i<HASHT; i++) {
     294   176136253 :     for (pt = hash[i]; pt; ) {
     295     2366114 :       void *z = (void*)(pt - 3);
     296     2366114 :       pt = (long*) pt[0]; pari_free(z);
     297             :     }
     298   173770139 :     hash[i] = NULL;
     299             :   }
     300      169385 : }
     301             : 
     302             : /* last prime stored */
     303             : ulong
     304           0 : GRH_last_prime(GRHcheck_t *S) { return (S->primes + S->nprimes-1)->p; }
     305             : /* ensure that S->primes can hold at least nb primes */
     306             : void
     307      401059 : GRH_ensure(GRHcheck_t *S, long nb)
     308             : {
     309      401059 :   if (S->maxprimes <= nb)
     310             :   {
     311           0 :     do S->maxprimes *= 2; while (S->maxprimes <= nb);
     312           0 :     pari_realloc_ip((void**)&S->primes, S->maxprimes*sizeof(*S->primes));
     313             :   }
     314      401059 : }
     315             : /* cache data for all primes up to the LIM */
     316             : static void
     317     1174824 : cache_prime_quad(GRHcheck_t *S, ulong LIM, GEN D)
     318             : {
     319             :   GRHprime_t *pr;
     320             :   long nb;
     321             : 
     322     1174824 :   if (S->limp >= LIM) return;
     323       72499 :   nb = (long)primepi_upper_bound((double)LIM); /* #{p <= LIM} <= nb */
     324       72507 :   GRH_ensure(S, nb+1); /* room for one extra prime */
     325       72507 :   for (pr = S->primes + S->nprimes;;)
     326    12670703 :   {
     327    12743210 :     ulong p = u_forprime_next(&(S->P));
     328    12742858 :     pr->p = p;
     329    12742858 :     pr->logp = log((double)p);
     330    12742858 :     pr->dec = (GEN)kroiu(D,p);
     331    12743209 :     S->nprimes++;
     332    12743209 :     pr++;
     333             :     /* store up to nextprime(LIM) included */
     334    12743209 :     if (p >= LIM) { S->limp = p; break; }
     335             :   }
     336             : }
     337             : 
     338             : static GEN
     339       70374 : compute_invresquad(GRHcheck_t *S, long LIMC)
     340             : {
     341       70374 :   pari_sp av = avma;
     342       70374 :   GEN invres = real_1(DEFAULTPREC);
     343       70374 :   double limp = log((double)LIMC) / 2;
     344             :   GRHprime_t *pr;
     345             :   long i;
     346    12726808 :   for (pr = S->primes, i = S->nprimes; i > 0; pr++, i--)
     347             :   {
     348    12658682 :     long s = (long)pr->dec;
     349    12658682 :     if (s)
     350             :     {
     351    12546526 :       ulong p = pr->p;
     352    12546526 :       if (s > 0 || pr->logp <= limp)
     353             :         /* Both p and P contribute */
     354     6377894 :         invres = mulur(p - s, divru(invres, p));
     355     6168632 :       else if (s<0)
     356             :         /* Only p contributes */
     357     6168609 :         invres = mulur(p, divru(invres, p - 1));
     358             :     }
     359             :   }
     360       68126 :   return gerepileuptoleaf(av, invres);
     361             : }
     362             : 
     363             : /* p | conductor of order of disc D ? */
     364             : static int
     365      392089 : is_bad(GEN D, ulong p)
     366             : {
     367      392089 :   pari_sp av = avma;
     368      392089 :   if (p == 2)
     369             :   {
     370       89692 :     long r = mod16(D) >> 1;
     371       89692 :     if (r && signe(D) < 0) r = 8-r;
     372       89692 :     return (r < 4);
     373             :   }
     374      302397 :   return gc_bool(av, dvdii(D, sqru(p))); /* p^2 | D ? */
     375             : }
     376             : 
     377             : /* returns the n-th suitable ideal for the factorbase */
     378             : static long
     379       70363 : nthidealquad(GEN D, long n)
     380             : {
     381       70363 :   pari_sp av = avma;
     382             :   forprime_t S;
     383             :   ulong p;
     384       70363 :   (void)u_forprime_init(&S, 2, ULONG_MAX);
     385      357602 :   while ((p = u_forprime_next(&S)) && n > 0)
     386      287242 :     if (!is_bad(D, p) && kroiu(D, p) >= 0) n--;
     387       70361 :   return gc_long(av, p);
     388             : }
     389             : 
     390             : static int
     391     1033961 : quadGRHchk(GEN D, GRHcheck_t *S, ulong LIMC)
     392             : {
     393     1033961 :   double logC = log((double)LIMC), SA = 0, SB = 0;
     394             :   long i;
     395             : 
     396     1033961 :   cache_prime_quad(S, LIMC, D);
     397     1033957 :   for (i = 0;; i++)
     398    29532780 :   {
     399    30566737 :     GRHprime_t *pr = S->primes+i;
     400    30566737 :     ulong p = pr->p;
     401             :     long M;
     402             :     double logNP, q, A, B;
     403    30566737 :     if (p > LIMC) break;
     404    29622097 :     if ((long)pr->dec < 0)
     405             :     {
     406    14789717 :       logNP = 2 * pr->logp;
     407    14789717 :       q = 1/(double)p;
     408             :     }
     409             :     else
     410             :     {
     411    14832380 :       logNP = pr->logp;
     412    14832380 :       q = 1/sqrt((double)p);
     413             :     }
     414    29622097 :     A = logNP * q; B = logNP * A; M = (long)(logC/logNP);
     415    29622097 :     if (M > 1)
     416             :     {
     417     2502445 :       double inv1_q = 1 / (1-q);
     418     2502445 :       A *= (1 - pow(q, (double) M)) * inv1_q;
     419     2502445 :       B *= (1 - pow(q, (double) M)*(M+1 - M*q)) * inv1_q * inv1_q;
     420             :     }
     421    29622097 :     if ((long)pr->dec>0) { SA += 2*A;SB += 2*B; } else { SA += A; SB += B; }
     422    29622097 :     if (p == LIMC) break;
     423             :   }
     424     1033957 :   return GRHok(S, logC, SA, SB);
     425             : }
     426             : 
     427             : /* C2 >= C1; create B->FB, B->numFB; set B->badprim. Return L(kro_D, 1) */
     428             : static void
     429       70510 : FBquad(struct buch_quad *B, ulong C2, ulong C1, GRHcheck_t *S)
     430             : {
     431       70510 :   GEN D = B->q->D;
     432             :   long i;
     433             :   pari_sp av;
     434             :   GRHprime_t *pr;
     435             : 
     436       70510 :   cache_prime_quad(S, C2, D);
     437       70510 :   pr = S->primes;
     438       70510 :   B->numFB = cgetg(C2+1, t_VECSMALL);
     439       70510 :   B->FB    = cgetg(C2+1, t_VECSMALL);
     440       70478 :   av = avma;
     441       70478 :   B->KC = 0; i = 0;
     442       70478 :   B->badprim = gen_1;
     443     2813577 :   for (;; pr++) /* p <= C2 */
     444     2813577 :   {
     445     2884055 :     ulong p = pr->p;
     446     2884055 :     if (!B->KC && p > C1) B->KC = i;
     447     2884055 :     if (p > C2) break;
     448     2821062 :     switch ((long)pr->dec)
     449             :     {
     450     1390192 :       case -1: break; /* inert */
     451      104850 :       case  0: /* ramified */
     452      104850 :         if (is_bad(D, p)) { B->badprim = muliu(B->badprim, p); break; }
     453             :         /* fall through */
     454             :       default:  /* split */
     455     1430839 :         i++; B->numFB[p] = i; B->FB[i] = p; break;
     456             :     }
     457     2821094 :     if (p == C2)
     458             :     {
     459        7517 :       if (!B->KC) B->KC = i;
     460        7517 :       break;
     461             :     }
     462             :   }
     463       70510 :   B->KC2 = i;
     464       70510 :   setlg(B->FB, B->KC2+1);
     465       70509 :   if (B->badprim != gen_1)
     466          42 :     B->badprim = gerepileuptoint(av, B->badprim);
     467             :   else
     468             :   {
     469       70467 :     B->badprim = NULL; set_avma(av);
     470             :   }
     471       70509 :   B->prodFB = zv_prod_Z(B->FB);
     472       70506 : }
     473             : 
     474             : /* create B->vperm, return B->subFB */
     475             : static GEN
     476       70507 : subFBquad(struct buch_quad *B, GEN D, double PROD, long minSFB)
     477             : {
     478       70507 :   long i, j, lgsub = 1, ino = 1, lv = B->KC+1;
     479       70507 :   double prod = 1.;
     480             :   pari_sp av;
     481             :   GEN no;
     482             : 
     483       70507 :   B->vperm = cgetg(lv, t_VECSMALL);
     484       70507 :   av = avma;
     485       70507 :   no    = cgetg(lv, t_VECSMALL);
     486      335115 :   for (j = 1; j < lv; j++)
     487             :   {
     488      334975 :     ulong p = uel(B->FB,j);
     489      334975 :     if (!umodiu(D, p)) no[ino++] = j; /* ramified */
     490             :     else
     491             :     {
     492      255037 :       B->vperm[lgsub++] = j;
     493      255037 :       prod *= p;
     494      255037 :       if (lgsub > minSFB && prod > PROD) break;
     495             :     }
     496             :   }
     497             :   /* lgsub >= 1 otherwise quadGRHchk is false */
     498       70509 :   i = lgsub;
     499      150448 :   for (j = 1; j < ino;i++,j++) B->vperm[i] = no[j];
     500     1165981 :   for (     ; i < lv; i++)     B->vperm[i] = i;
     501       70509 :   no = gclone(vecslice(B->vperm, 1, lgsub-1));
     502       70508 :   set_avma(av); return no;
     503             : }
     504             : 
     505             : /* assume n >= 1, x[i][j] = B->subFB[i]^j, for j = 1..n */
     506             : static GEN
     507       99239 : powsubFBquad(struct buch_quad *B, long n)
     508             : {
     509       99239 :   pari_sp av = avma;
     510       99239 :   long i,j, l = lg(B->subFB);
     511       99239 :   GEN F, y, x = cgetg(l, t_VEC), D = B->q->D;
     512             : 
     513       99239 :   if (B->PRECREG) /* real */
     514             :   {
     515       39627 :     for (i=1; i<l; i++)
     516             :     {
     517       34510 :       F = qfr5_pf(B->q, B->FB[B->subFB[i]], B->PRECREG);
     518       34510 :       y = cgetg(n+1, t_VEC); gel(x,i) = y;
     519       34510 :       gel(y,1) = F;
     520      552160 :       for (j=2; j<=n; j++) gel(y,j) = QFR5_comp(gel(y,j-1), F, B->q);
     521             :     }
     522             :   }
     523             :   else /* imaginary */
     524             :   {
     525      427148 :     for (i=1; i<l; i++)
     526             :     {
     527      333061 :       F = qfi_pf(D, B->FB[B->subFB[i]]);
     528      333020 :       y = cgetg(n+1, t_VEC); gel(x,i) = y;
     529      333935 :       gel(y,1) = F;
     530     5318775 :       for (j=2; j<=n; j++) gel(y,j) = qfbcomp_i(gel(y,j-1), F);
     531             :     }
     532             :   }
     533       99204 :   x = gclone(x); set_avma(av); return x;
     534             : }
     535             : 
     536             : static void
     537       98298 : sub_fact(struct buch_quad *B, GEN col, GEN F)
     538             : {
     539       98298 :   GEN b = gel(F,2);
     540             :   long i;
     541      208018 :   for (i=1; i<=B->primfact[0]; i++)
     542             :   {
     543      109720 :     ulong p = B->primfact[i], k = B->numFB[p];
     544      109720 :     long e = B->exprimfact[i];
     545      109720 :     if (umodiu(b, p<<1) > p) e = -e;
     546      109720 :     col[k] -= e;
     547             :   }
     548       98298 : }
     549             : 
     550             : #if 0
     551             : static void
     552             : dbg_fact(struct buch_quad *B)
     553             : {
     554             :   long i;
     555             :   for (i=1; i<=B->primfact[0]; i++)
     556             :   {
     557             :     ulong p = B->primfact[i];
     558             :     long e = B->exprimfact[i];
     559             :     err_printf("%lu^%ld ",p,e);
     560             :   }
     561             : }
     562             : 
     563             : static void
     564             : chk_fact(struct buch_quad *B, GEN col)
     565             : {
     566             :   long i, l = lg(col);
     567             :   GEN Q = qfi_pf(B->q->D, 1);
     568             :   for (i=1; i< l; i++)
     569             :   {
     570             :     ulong p = B->FB[i];
     571             :     long k = col[i];
     572             :     Q = qfbcomp(qfbpowraw(qfi_pf(B->q->D, p),k),Q);
     573             :   }
     574             :   if (!gequal1(gel(Q,1))) pari_err_BUG("chk_fact");
     575             : }
     576             : #endif
     577             : 
     578             : static void
     579     1892427 : add_fact(struct buch_quad *B, GEN col, GEN F)
     580             : {
     581     1892427 :   GEN b = gel(F,2);
     582             :   long i;
     583     5930690 :   for (i=1; i<=B->primfact[0]; i++)
     584             :   {
     585     4038145 :     ulong p = B->primfact[i], k = B->numFB[p];
     586     4038145 :     long e = B->exprimfact[i];
     587     4038145 :     if (umodiu(b, p<<1) > p) e = -e;
     588     4038263 :     col[k] += e;
     589             :   }
     590     1892545 : }
     591             : 
     592             : static GEN
     593       70363 : get_clgp(struct buch_quad *B, GEN W, GEN *ptD)
     594             : {
     595       70363 :   GEN res, init, u1, D = ZM_snf_group(W,NULL,&u1);
     596       70363 :   long i, j, l = lg(W), c = lg(D);
     597             : 
     598       70363 :   res=cgetg(c,t_VEC); init = cgetg(l,t_VEC);
     599      216240 :   for (i=1; i<l; i++) gel(init,i) = primeform_u(B->q->D, B->FB[B->vperm[i]]);
     600      195869 :   for (j=1; j<c; j++)
     601             :   {
     602      125510 :     GEN g = NULL;
     603      125510 :     if (signe(B->q->D) > 0)
     604             :     {
     605       13328 :       for (i=1; i<l; i++)
     606             :       {
     607        9618 :         GEN t, u = gcoeff(u1,i,j);
     608        9618 :         if (!signe(u)) continue;
     609        4543 :         t = qfr3_pow(gel(init,i), u, B->q);
     610        4543 :         g = g? qfr3_comp(g, t, B->q): t;
     611             :       }
     612        3710 :       g = qfr3_to_qfr(qfr3_canon_safe(qfr3_red(g, B->q), B->q), B->q->D);
     613             :     }
     614             :     else
     615             :     {
     616      423812 :       for (i=1; i<l; i++)
     617             :       {
     618      302014 :         GEN t, u = gcoeff(u1,i,j);
     619      302014 :         if (!signe(u)) continue;
     620      203592 :         t = powgi(gel(init,i), u);
     621      203597 :         g = g? qfbcomp_i(g, t): t;
     622             :       }
     623             :     }
     624      125508 :     gel(res,j) = g;
     625             :   }
     626       70359 :   *ptD = D; return res;
     627             : }
     628             : 
     629             : static long
     630       70377 : trivial_relations(struct buch_quad *B, GEN mat, GEN C)
     631             : {
     632       70377 :   long i, j = 0;
     633       70377 :   GEN col, D = B->q->D;
     634     1500393 :   for (i = 1; i <= B->KC; i++)
     635             :   { /* ramified prime ==> trivial relation */
     636     1430016 :     if (umodiu(D, B->FB[i])) continue;
     637      104434 :     col = zero_zv(B->KC);
     638      104434 :     col[i] = 2; j++;
     639      104434 :     gel(mat,j) = col;
     640      104434 :     gel(C,j) = gen_0;
     641             :   }
     642       70377 :   return j;
     643             : }
     644             : 
     645             : static void
     646           0 : dbg_all(pari_timer *T, const char *phase, long s, long n)
     647             : {
     648           0 :   err_printf("\n");
     649           0 :   timer_printf(T, "%s rel [#rel/#test = %ld/%ld]", phase,s,n);
     650           0 : }
     651             : 
     652             : /* Imaginary Quadratic fields */
     653             : 
     654             : static void
     655        8813 : rel_to_col(struct buch_quad *B, GEN col, GEN rel, GEN b)
     656             : {
     657        8813 :   GEN P = gel(rel, 1), E = gel(rel, 2);
     658        8813 :   long i, lP = lg(P);
     659      105784 :   for (i=1; i<lP; i++)
     660             :   {
     661       96971 :     ulong p = uel(P, i), e = uel(E, i);
     662       96971 :     col[B->numFB[p]] += umodiu(b, p<<1) > p ? -e :e;
     663             :   }
     664        8813 : }
     665             : 
     666             : static void
     667       99047 : imag_relations(struct buch_quad *B, long need, long *pc, ulong LIMC, GEN mat)
     668             : {
     669             :   pari_timer T;
     670       99047 :   long lgsub = lg(B->subFB), current = *pc, nbtest = 0, s = 0;
     671             :   long i, fpc;
     672             :   pari_sp av;
     673       99047 :   GEN col, form, ex = cgetg(lgsub, t_VECSMALL);
     674             : 
     675       99047 :   if (!current) current = 1;
     676       99047 :   if (DEBUGLEVEL>2) timer_start(&T);
     677       99047 :   av = avma;
     678             :   for(;;)
     679             :   {
     680     3640393 :     if (s >= need) break;
     681     3541347 :     set_avma(av);
     682     3541374 :     form = qfi_random(B,ex);
     683     3538503 :     form = qfbcomp_i(form, qfi_pf(B->q->D, B->FB[current]));
     684     3538055 :     nbtest++; fpc = factorquad(B,form,B->KC,LIMC);
     685     3540686 :     if (!fpc)
     686             :     {
     687      291564 :       if (DEBUGLEVEL>3) err_printf(".");
     688      291564 :       if ((nbtest & 0xff) == 0 && ++current > B->KC) current = 1;
     689      291564 :       continue;
     690             :     }
     691     3249122 :     if (fpc > 1)
     692             :     {
     693     1800115 :       long *fpd = largeprime(B,fpc,ex,current,0);
     694             :       ulong b1, b2, p;
     695             :       GEN form2;
     696     1800542 :       if (!fpd)
     697             :       {
     698     1639948 :         if (DEBUGLEVEL>3) err_printf(".");
     699     1639948 :         continue;
     700             :       }
     701      160594 :       form2 = qfbcomp_i(qfi_factorback(B,fpd), qfi_pf(B->q->D, B->FB[fpd[-2]]));
     702      160623 :       p = fpc << 1;
     703      160623 :       b1 = umodiu(gel(form2,2), p);
     704      160624 :       b2 = umodiu(gel(form,2),  p);
     705      160624 :       if (b1 != b2 && b1+b2 != p) continue;
     706             : 
     707      160624 :       col = gel(mat,++s);
     708      160624 :       add_fact(B,col, form);
     709      160627 :       (void)factorquad(B,form2,B->KC,LIMC);
     710      160628 :       if (b1==b2)
     711             :       {
     712      414949 :         for (i=1; i<lgsub; i++) col[B->subFB[i]] += fpd[i]-ex[i];
     713       80671 :         sub_fact(B, col, form2); col[fpd[-2]]++;
     714             :       }
     715             :       else
     716             :       {
     717      410316 :         for (i=1; i<lgsub; i++) col[B->subFB[i]] += -fpd[i]-ex[i];
     718       79957 :         add_fact(B, col, form2); col[fpd[-2]]--;
     719             :       }
     720      160631 :       if (DEBUGLEVEL>2) err_printf(" %ldP",s);
     721             :     }
     722             :     else
     723             :     {
     724     1449007 :       col = gel(mat,++s);
     725     6959363 :       for (i=1; i<lgsub; i++) col[B->subFB[i]] = -ex[i];
     726     1449007 :       add_fact(B, col, form);
     727     1449433 :       if (DEBUGLEVEL>2) err_printf(" %ld",s);
     728             :     }
     729     1609834 :     col[current]--;
     730     1609834 :     if (++current > B->KC) current = 1;
     731             :   }
     732       99046 :   if (DEBUGLEVEL>2) dbg_all(&T, "random", s, nbtest);
     733       99046 :   *pc = current;
     734       99046 : }
     735             : 
     736             : static void
     737         336 : mpqs_relations(struct buch_quad *B, long need, long *pc, ulong LIMC, GEN mat, mpqs_handle_t *H, GEN missing_primes)
     738             : {
     739             :   pari_timer T;
     740             :   long i, lV;
     741             :   GEN V;
     742         336 :   if (DEBUGLEVEL>2) timer_start(&T);
     743         336 :   V = mpqs_class_rels(H, need, missing_primes);
     744         336 :   if (!V) { imag_relations(B, need, pc, LIMC, mat); return; }
     745         336 :   lV = lg(V);
     746        9149 :   for (i = 1; i < lV && i <= need; i++)
     747             :   {
     748        8813 :     GEN formA = gel(V,i), rel = gel(formA,2), b = gel(formA,1);
     749        8813 :     GEN col = gel(mat,i);
     750        8813 :     rel_to_col(B, col, rel, b);
     751             :   }
     752         336 :   if (DEBUGLEVEL>2) timer_printf(&T, "MPQS rel [#rel = %ld]", i-1);
     753             : }
     754             : 
     755             : static int
     756           7 : imag_be_honest(struct buch_quad *B)
     757             : {
     758           7 :   long p, fpc, s = B->KC, nbtest = 0;
     759           7 :   GEN F, ex = cgetg(lg(B->subFB), t_VECSMALL);
     760           7 :   pari_sp av = avma;
     761             : 
     762         525 :   while (s<B->KC2)
     763             :   {
     764         518 :     p = B->FB[s+1]; if (DEBUGLEVEL>2) err_printf(" %ld",p);
     765         518 :     F = qfbcomp_i(qfi_pf(B->q->D, p), qfi_random(B, ex));
     766         518 :     fpc = factorquad(B,F,s,p-1);
     767         518 :     if (fpc == 1) { nbtest=0; s++; }
     768             :     else
     769         392 :       if (++nbtest > 40) return 0;
     770         518 :     set_avma(av);
     771             :   }
     772           7 :   return 1;
     773             : }
     774             : 
     775             : static GEN
     776      184688 : dist(GEN e, GEN d, long prec)
     777             : {
     778      184688 :   GEN t = qfr5_dist(e, d, prec);
     779      184688 :   return signe(d) < 0 ? mkcomplex(t, gen_1): t;
     780             : }
     781             : 
     782             : /* Real Quadratic fields */
     783             : 
     784             : static void
     785        5145 : real_relations(struct buch_quad *B, long need, long *pc, long lim, ulong LIMC, GEN mat, GEN C)
     786             : {
     787             :   pari_timer T;
     788        5145 :   long lgsub = lg(B->subFB), prec = B->PRECREG, current = *pc, nbtest=0, s=0;
     789             :   long i, fpc, endcycle, rhoacc, rho;
     790             :   /* in a 2nd phase, don't include FB[current] but run along the cyle
     791             :    * ==> get more units */
     792        5145 :   int first = (current == 0);
     793             :   pari_sp av, av1;
     794        5145 :   GEN d, col, form, form0, form1, ex = cgetg(lgsub, t_VECSMALL);
     795             : 
     796        5145 :   if (DEBUGLEVEL>2) timer_start(&T);
     797        5145 :   if (!current) current = 1;
     798        5145 :   if (lim > need) lim = need;
     799        5145 :   av = avma;
     800             :   for(;;)
     801             :   {
     802      166691 :     if (s >= need) break;
     803      161546 :     if (first && s >= lim) {
     804        2142 :       first = 0;
     805        2142 :       if (DEBUGLEVEL>2) dbg_all(&T, "initial", s, nbtest);
     806             :     }
     807      161546 :     set_avma(av); form = qfr3_random(B, ex);
     808      161546 :     if (!first)
     809      159369 :       form = QFR3_comp(form, qfr3_pf(B->q, B->FB[current]), B->q);
     810      161546 :     av1 = avma;
     811      161546 :     form0 = form; form1 = NULL;
     812      161546 :     endcycle = rhoacc = 0;
     813      161546 :     rho = -1;
     814             : 
     815     1300964 : CYCLE:
     816     1300964 :     if (endcycle || rho > 5000)
     817             :     {
     818          21 :       if (++current > B->KC) current = 1;
     819          21 :       continue;
     820             :     }
     821     1300943 :     if (gc_needed(av,1))
     822             :     {
     823           0 :       if(DEBUGMEM>1) pari_warn(warnmem,"real_relations");
     824           0 :       gerepileall(av1, form1? 2: 1, &form, &form1);
     825             :     }
     826     1300943 :     if (rho < 0) rho = 0; /* first time in */
     827             :     else
     828             :     {
     829     1139397 :       form = qfr3_rho(form, B->q); rho++;
     830     1139397 :       rhoacc++;
     831     1139397 :       if (first)
     832      442526 :         endcycle = (absequalii(gel(form,1),gel(form0,1))
     833      221263 :              && equalii(gel(form,2),gel(form0,2)));
     834             :       else
     835             :       {
     836      918134 :         if (absequalii(gel(form,1), gel(form,3))) /* a = -c */
     837             :         {
     838           0 :           if (absequalii(gel(form,1),gel(form0,1)) &&
     839           0 :                   equalii(gel(form,2),gel(form0,2))) continue;
     840           0 :           form = qfr3_rho(form, B->q); rho++;
     841           0 :           rhoacc++;
     842             :         }
     843             :         else
     844      918134 :           { setsigne(form[1],1); setsigne(form[3],-1); }
     845      918176 :         if (equalii(gel(form,1),gel(form0,1)) &&
     846          42 :             equalii(gel(form,2),gel(form0,2))) continue;
     847             :       }
     848             :     }
     849     1300943 :     nbtest++; fpc = factorquad(B,form,B->KC,LIMC);
     850     1300943 :     if (!fpc)
     851             :     {
     852      386778 :       if (DEBUGLEVEL>3) err_printf(".");
     853      386778 :       goto CYCLE;
     854             :     }
     855      914165 :     if (fpc > 1)
     856             :     { /* look for Large Prime relation */
     857      764946 :       long *fpd = largeprime(B,fpc,ex,first? 0: current,rhoacc);
     858             :       ulong b1, b2, p;
     859             :       GEN form2;
     860      764946 :       if (!fpd)
     861             :       {
     862      729477 :         if (DEBUGLEVEL>3) err_printf(".");
     863      729477 :         goto CYCLE;
     864             :       }
     865       35469 :       if (!form1)
     866             :       {
     867       35469 :         form1 = qfr5_factorback(B,ex);
     868       35469 :         if (!first)
     869       35469 :           form1 = QFR5_comp(form1, qfr5_pf(B->q, B->FB[current], prec), B->q);
     870             :       }
     871       35469 :       form1 = qfr5_rho_pow(form1, rho, B->q);
     872       35469 :       rho = 0;
     873             : 
     874       35469 :       form2 = qfr5_factorback(B,fpd);
     875       35469 :       if (fpd[-2])
     876       23709 :         form2 = QFR5_comp(form2, qfr5_pf(B->q, B->FB[fpd[-2]], prec), B->q);
     877       35469 :       form2 = qfr5_rho_pow(form2, fpd[-3], B->q);
     878       35469 :       if (!absequalii(gel(form2,1),gel(form2,3)))
     879             :       {
     880       35469 :         setsigne(form2[1], 1);
     881       35469 :         setsigne(form2[3],-1);
     882             :       }
     883       35469 :       p = fpc << 1;
     884       35469 :       b1 = umodiu(gel(form2,2), p);
     885       35469 :       b2 = umodiu(gel(form1,2), p);
     886       35469 :       if (b1 != b2 && b1+b2 != p) goto CYCLE;
     887             : 
     888       35469 :       col = gel(mat,++s);
     889       35469 :       add_fact(B, col, form1);
     890       35469 :       (void)factorquad(B,form2,B->KC,LIMC);
     891       35469 :       if (b1==b2)
     892             :       {
     893      135100 :         for (i=1; i<lgsub; i++) col[B->subFB[i]] += fpd[i]-ex[i];
     894       17626 :         sub_fact(B,col, form2);
     895       17626 :         if (fpd[-2]) col[fpd[-2]]++;
     896       17626 :         d = dist(subii(gel(form1,4),gel(form2,4)),
     897       17626 :                       divrr(gel(form1,5),gel(form2,5)), prec);
     898             :       }
     899             :       else
     900             :       {
     901      136780 :         for (i=1; i<lgsub; i++) col[B->subFB[i]] += -fpd[i]-ex[i];
     902       17843 :         add_fact(B, col, form2);
     903       17843 :         if (fpd[-2]) col[fpd[-2]]--;
     904       17843 :         d = dist(addii(gel(form1,4),gel(form2,4)),
     905       17843 :                       mulrr(gel(form1,5),gel(form2,5)), prec);
     906             :       }
     907       35469 :       if (DEBUGLEVEL>2) err_printf(" %ldP",s);
     908             :     }
     909             :     else
     910             :     { /* standard relation */
     911      149219 :       if (!form1)
     912             :       {
     913      126077 :         form1 = qfr5_factorback(B, ex);
     914      126077 :         if (!first)
     915      123900 :           form1 = QFR5_comp(form1, qfr5_pf(B->q, B->FB[current], prec), B->q);
     916             :       }
     917      149219 :       form1 = qfr5_rho_pow(form1, rho, B->q);
     918      149219 :       rho = 0;
     919             : 
     920      149219 :       col = gel(mat,++s);
     921     1147489 :       for (i=1; i<lgsub; i++) col[B->subFB[i]] = -ex[i];
     922      149219 :       add_fact(B, col, form1);
     923      149219 :       d = dist(gel(form1,4), gel(form1,5), prec);
     924      149219 :       if (DEBUGLEVEL>2) err_printf(" %ld",s);
     925             :     }
     926      184688 :     gaffect(d, gel(C,s));
     927      184688 :     if (first)
     928             :     {
     929       25319 :       if (s >= lim) continue;
     930       23163 :       goto CYCLE;
     931             :     }
     932             :     else
     933             :     {
     934      159369 :       col[current]--;
     935      159369 :       if (++current > B->KC) current = 1;
     936             :     }
     937             :   }
     938        5145 :   if (DEBUGLEVEL>2) dbg_all(&T, "random", s, nbtest);
     939        5145 :   *pc = current;
     940        5145 : }
     941             : 
     942             : static int
     943           7 : real_be_honest(struct buch_quad *B)
     944             : {
     945           7 :   long p, fpc, s = B->KC, nbtest = 0;
     946           7 :   GEN F,F0, ex = cgetg(lg(B->subFB), t_VECSMALL);
     947           7 :   pari_sp av = avma;
     948             : 
     949          28 :   while (s<B->KC2)
     950             :   {
     951          21 :     p = B->FB[s+1]; if (DEBUGLEVEL>2) err_printf(" %ld",p);
     952          21 :     F = QFR3_comp(qfr3_random(B, ex), qfr3_pf(B->q, p), B->q);
     953          21 :     for (F0 = F;;)
     954             :     {
     955          49 :       fpc = factorquad(B,F,s,p-1);
     956          49 :       if (fpc == 1) { nbtest=0; s++; break; }
     957          28 :       if (++nbtest > 40) return 0;
     958          28 :       F = qfr3_canon(qfr3_rho(F, B->q), B->q);
     959          28 :       if (equalii(gel(F,1),gel(F0,1))
     960           0 :        && equalii(gel(F,2),gel(F0,2))) break;
     961             :     }
     962          21 :     set_avma(av);
     963             :   }
     964           7 :   return 1;
     965             : }
     966             : 
     967             : static GEN
     968       55146 : crabs(GEN a)
     969             : {
     970       55146 :   return signe(real_i(a)) < 0 ? gneg(a): a;
     971             : }
     972             : 
     973             : static GEN
     974       33978 : gcdreal(GEN a,GEN b)
     975             : {
     976       33978 :   if (!signe(real_i(a))) return crabs(b);
     977       33012 :   if (!signe(real_i(b))) return crabs(a);
     978       32897 :   if (expo(real_i(a))<-5) return crabs(b);
     979       12040 :   if (expo(real_i(b))<-5) return crabs(a);
     980        9093 :   a = crabs(a); b = crabs(b);
     981       20019 :   while (expo(real_i(b)) >= -5  && signe(real_i(b)))
     982             :   {
     983             :     long e;
     984       10926 :     GEN r, q = gcvtoi(divrr(real_i(a),real_i(b)),&e);
     985       10926 :     if (e > 0) return NULL;
     986       10926 :     r = gsub(a, gmul(q,b)); a = b; b = r;
     987             :   }
     988        9093 :   return crabs(a);
     989             : }
     990             : 
     991             : static int
     992       91616 : get_R(struct buch_quad *B, GEN C, long sreg, GEN z, GEN *ptR)
     993             : {
     994       91616 :   GEN R = gen_1;
     995             :   double c;
     996             :   long i;
     997       91616 :   if (B->PRECREG)
     998             :   {
     999        2982 :     R = crabs(gel(C,1));
    1000       36960 :     for (i=2; i<=sreg; i++)
    1001             :     {
    1002       33978 :       R = gcdreal(gel(C,i), R);
    1003       33978 :       if (!R) return fupb_PRECI;
    1004             :     }
    1005        2982 :     if (gexpo(real_i(R)) <= -3)
    1006             :     {
    1007           0 :       if (DEBUGLEVEL>2) err_printf("regulator is zero.\n");
    1008           0 :       return fupb_RELAT;
    1009             :     }
    1010        2982 :     if (DEBUGLEVEL>2) err_printf("#### Tentative regulator: %Ps\n",R);
    1011             :   }
    1012       91616 :   c = gtodouble(gmul(z, real_i(R)));
    1013       91618 :   if (c < 0.8 || c > 1.3) return fupb_RELAT;
    1014       70357 :   *ptR = R; return fupb_NONE;
    1015             : }
    1016             : 
    1017             : static int
    1018       70357 : quad_be_honest(struct buch_quad *B)
    1019             : {
    1020             :   int r;
    1021       70357 :   if (B->KC2 <= B->KC) return 1;
    1022          14 :   if (DEBUGLEVEL>2)
    1023           0 :     err_printf("be honest for primes from %ld to %ld\n", B->FB[B->KC+1],B->FB[B->KC2]);
    1024          14 :   r = B->PRECREG? real_be_honest(B): imag_be_honest(B);
    1025          14 :   if (DEBUGLEVEL>2) err_printf("\n");
    1026          14 :   return r;
    1027             : }
    1028             : 
    1029             : static GEN
    1030       70536 : Buchquad_i(GEN D, double cbach, double cbach2, long prec)
    1031             : {
    1032       70536 :   const long MAXRELSUP = 20, SFB_MAX = 3, MPQS_THRESHOLD = 60;
    1033             :   pari_timer T;
    1034             :   pari_sp av, av2;
    1035       70536 :   const long RELSUP = 5;
    1036             :   long i, s, current, triv, sfb_trials, nrelsup, nreldep, need, nsubFB, minSFB;
    1037             :   ulong low, high, LIMC0, LIMC, LIMC2, LIMCMAX, cp;
    1038       70536 :   GEN W, cyc, gen, dep, mat, C, extraC, B, R, invhr, h = NULL; /*-Wall*/
    1039             :   double drc, sdrc, lim, LOGD, LOGD2;
    1040             :   GRHcheck_t GRHcheck;
    1041             :   struct qfr_data q;
    1042             :   struct buch_quad BQ;
    1043       70536 :   int FIRST = 1, use_mpqs = 0;
    1044             :   mpqs_handle_t H;
    1045             :   GEN missing_primes;
    1046             : 
    1047       70536 :   check_quaddisc(D, &s, /*junk*/&i, "Buchquad");
    1048       70538 :   R = NULL; /* -Wall */
    1049       70538 :   BQ.q = &q;
    1050       70538 :   q.D = D;
    1051       70538 :   if (s < 0)
    1052             :   {
    1053       68382 :     if (abscmpiu(q.D,4) <= 0)
    1054         175 :       retmkvec4(gen_1, cgetg(1,t_VEC), cgetg(1,t_VEC), gen_1);
    1055       68207 :     prec = BQ.PRECREG = 0;
    1056       68207 :     if (expi(D) >= MPQS_THRESHOLD)
    1057          28 :       use_mpqs = 1;
    1058             :   } else {
    1059        2156 :     BQ.PRECREG = maxss(prec+EXTRAPREC64, nbits2prec(2*expi(q.D) + 128));
    1060             :   }
    1061       70363 :   if (DEBUGLEVEL>2) timer_start(&T);
    1062       70363 :   BQ.primfact   = new_chunk(100);
    1063       70363 :   BQ.exprimfact = new_chunk(100);
    1064       70363 :   BQ.hashtab = (long**) new_chunk(HASHT);
    1065    72117518 :   for (i=0; i<HASHT; i++) BQ.hashtab[i] = NULL;
    1066             : 
    1067       70363 :   drc = fabs(gtodouble(q.D));
    1068       70362 :   LOGD = log(drc);
    1069       70362 :   LOGD2 = LOGD * LOGD;
    1070             : 
    1071       70362 :   sdrc = lim = sqrt(drc);
    1072       70362 :   if (!BQ.PRECREG) lim /= sqrt(3.);
    1073       70362 :   cp = (ulong)exp(sqrt(LOGD*log(LOGD)/8.0));
    1074       70362 :   if (cp < 20) cp = 20;
    1075       70362 :   if (cbach > 6.) {
    1076           0 :     if (cbach2 < cbach) cbach2 = cbach;
    1077           0 :     cbach = 6.;
    1078             :   }
    1079       70362 :   if (cbach < 0.)
    1080           0 :     pari_err_DOMAIN("Buchquad","Bach constant","<",gen_0,dbltor(cbach));
    1081       70362 :   av = avma;
    1082       70362 :   BQ.powsubFB = BQ.subFB = NULL;
    1083       70362 :   minSFB = (expi(D) > 15)? 3: 2;
    1084       70363 :   init_GRHcheck(&GRHcheck, 2, BQ.PRECREG? 2: 0, LOGD);
    1085       70363 :   high = low = LIMC0 = maxss((long)(cbach2*LOGD2), 1);
    1086       70363 :   LIMCMAX = (long)(4.*LOGD2);
    1087             :   /* 97/1223 below to ensure a good enough approximation of residue */
    1088       70363 :   cache_prime_quad(&GRHcheck, expi(D) < 16 ? 97: 1223, D);
    1089      587313 :   while (!quadGRHchk(D, &GRHcheck, high))
    1090             :   {
    1091      516951 :     low = high;
    1092      516951 :     high *= 2;
    1093             :   }
    1094      517027 :   while (high - low > 1)
    1095             :   {
    1096      446666 :     long test = (low+high)/2;
    1097      446666 :     if (quadGRHchk(D, &GRHcheck, test))
    1098      233960 :       high = test;
    1099             :     else
    1100      212706 :       low = test;
    1101             :   }
    1102       70361 :   if (high == LIMC0+1 && quadGRHchk(D, &GRHcheck, LIMC0))
    1103           0 :     LIMC2 = LIMC0;
    1104             :   else
    1105       70361 :     LIMC2 = high;
    1106       70361 :   if (LIMC2 > LIMCMAX) LIMC2 = LIMCMAX;
    1107       70361 :   LIMC0 = (long)(cbach*LOGD2);
    1108       70361 :   LIMC = cbach ? LIMC0 : LIMC2;
    1109       70361 :   LIMC = maxss(LIMC, nthidealquad(D, 2));
    1110             : 
    1111             : /* LIMC = Max(cbach*(log D)^2, exp(sqrt(log D loglog D) / 8)) */
    1112       70377 : START:
    1113       70377 :   missing_primes = NULL;
    1114             :   do
    1115             :   {
    1116       70510 :     if (!FIRST) LIMC = bnf_increase_LIMC(LIMC,LIMCMAX);
    1117       70510 :     if (DEBUGLEVEL>2 && LIMC > LIMC0)
    1118           0 :       err_printf("%s*** Bach constant: %f\n", FIRST?"":"\n", LIMC/LOGD2);
    1119       70510 :     FIRST = 0; set_avma(av);
    1120       70510 :     guncloneNULL(BQ.subFB);
    1121       70510 :     guncloneNULL(BQ.powsubFB);
    1122       70510 :     clearhash(BQ.hashtab);
    1123       70510 :     if (LIMC < cp) LIMC = cp;
    1124       70510 :     if (LIMC2 < LIMC) LIMC2 = LIMC;
    1125       70510 :     if (BQ.PRECREG) qfr_data_init(q.D, BQ.PRECREG, &q);
    1126             : 
    1127       70510 :     FBquad(&BQ, LIMC2, LIMC, &GRHcheck);
    1128       70506 :     if (DEBUGLEVEL>2) timer_printf(&T, "factor base");
    1129       70506 :     BQ.subFB = subFBquad(&BQ, q.D, lim + 0.5, minSFB);
    1130       70507 :     if (DEBUGLEVEL>2) timer_printf(&T, "subFBquad = %Ps",
    1131             :                                  vecpermute(BQ.FB, BQ.subFB));
    1132       70507 :     nsubFB = lg(BQ.subFB) - 1;
    1133             :   }
    1134       70507 :   while (nsubFB < (expi(D) > 15 ? 3 : 2));
    1135             :   /* invhr = 2^r1 (2pi)^r2 / sqrt(D) w ~ L(chi,1) / hR */
    1136       70374 :   invhr = gmul(dbltor((BQ.PRECREG?2.:M_PI)/sdrc),
    1137             :                compute_invresquad(&GRHcheck, LIMC));
    1138       70373 :   BQ.powsubFB = powsubFBquad(&BQ,CBUCH+1);
    1139       70377 :   if (DEBUGLEVEL>2) timer_printf(&T, "powsubFBquad");
    1140       70377 :   BQ.limhash = (LIMC & HIGHMASK)? (HIGHBIT>>1): LIMC*LIMC;
    1141             : 
    1142       70377 :   need = BQ.KC + RELSUP - 2;
    1143       70377 :   current = 0;
    1144       70377 :   W = NULL;
    1145       70377 :   sfb_trials = nreldep = nrelsup = 0;
    1146       70377 :   s = nsubFB + RELSUP;
    1147       70377 :   if (use_mpqs)
    1148          28 :     use_mpqs = mpqs_class_init(&H, D, BQ.KC);
    1149       70377 :   av2 = avma;
    1150             : 
    1151             :   do
    1152             :   {
    1153      104529 :     if ((nreldep & 3) == 1 || (nrelsup & 7) == 1) {
    1154       28864 :       if (DEBUGLEVEL>2) err_printf("*** Changing sub factor base\n");
    1155       28864 :       gunclone(BQ.subFB);
    1156       28864 :       gunclone(BQ.powsubFB);
    1157       28864 :       BQ.subFB = gclone(vecslice(BQ.vperm, 1, nsubFB));
    1158       28864 :       BQ.powsubFB = powsubFBquad(&BQ,CBUCH+1);
    1159       28864 :       if (DEBUGLEVEL>2) timer_printf(&T, "powsubFBquad");
    1160       28864 :       clearhash(BQ.hashtab);
    1161             :     }
    1162      104529 :     if (!use_mpqs) need += 2;
    1163      104529 :     mat    = cgetg(need+1, t_MAT);
    1164      104528 :     extraC = cgetg(need+1, t_VEC);
    1165      104529 :     if (!W) { /* first time */
    1166       70377 :       C = extraC;
    1167       70377 :       triv = trivial_relations(&BQ, mat, C);
    1168       70377 :       if (DEBUGLEVEL>2) err_printf("KC = %ld, need %ld relations\n", BQ.KC, need);
    1169             :     } else {
    1170       34152 :       triv = 0;
    1171       34152 :       if (DEBUGLEVEL>2) err_printf("...need %ld more relations\n", need);
    1172             :     }
    1173      104529 :     if (BQ.PRECREG) {
    1174      189833 :       for (i = triv+1; i<=need; i++) {
    1175      184688 :         gel(mat,i) = zero_zv(BQ.KC);
    1176      184688 :         gel(extraC,i) = mkcomplex(cgetr(BQ.PRECREG), cgeti(3));
    1177             :       }
    1178        5145 :       real_relations(&BQ, need - triv, &current, s,LIMC,mat + triv,extraC + triv);
    1179             :     } else {
    1180     1718528 :       for (i = triv+1; i<=need; i++) {
    1181     1619150 :         gel(mat,i) = zero_zv(BQ.KC);
    1182     1619144 :         gel(extraC,i) = gen_0;
    1183             :       }
    1184       99378 :       if (use_mpqs)
    1185         336 :         mpqs_relations(&BQ, need - triv, &current, LIMC,mat + triv, &H, missing_primes);
    1186             :       else
    1187       99042 :         imag_relations(&BQ, need - triv, &current, LIMC,mat + triv);
    1188             :     }
    1189      104527 :     if (DEBUGLEVEL>2) timer_start(&T);
    1190      104528 :     if (!W)
    1191       70376 :       W = hnfspec_i(mat,BQ.vperm,&dep,&B,&C,nsubFB);
    1192             :     else
    1193       34152 :       W = hnfadd_i(W,BQ.vperm,&dep,&B,&C, mat,extraC);
    1194      104529 :     if (DEBUGLEVEL>2) timer_printf(&T, "hnfadd");
    1195      104529 :     need = BQ.KC - (lg(W)-1) - (lg(B)-1);
    1196      104529 :     gerepileall(av2, 4, &W,&C,&B,&dep);
    1197      104528 :     missing_primes = vecslice(BQ.vperm,1,need);
    1198      104528 :     if (need)
    1199             :     {
    1200       12905 :       if (++nreldep > 15 && cbach < 1) goto START;
    1201       12905 :       continue;
    1202             :     }
    1203             : 
    1204       91623 :     h = ZM_det_triangular(W);
    1205       91616 :     if (DEBUGLEVEL>2) err_printf("\n#### Tentative class number: %Ps\n", h);
    1206             : 
    1207       91616 :     switch(get_R(&BQ, C, (lg(C)-1) - (lg(B)-1) - (lg(W)-1), mulir(h,invhr), &R))
    1208             :     {
    1209           0 :       case fupb_PRECI:
    1210           0 :         BQ.PRECREG = precdbl(BQ.PRECREG);
    1211           0 :         FIRST = 1; goto START;
    1212             : 
    1213       21261 :       case fupb_RELAT:
    1214       21261 :         if (++nrelsup > MAXRELSUP)
    1215             :         {
    1216          63 :           if (++sfb_trials > SFB_MAX && cbach <= 1) goto START;
    1217          49 :           if (nsubFB < minss(10,BQ.KC)) nsubFB++;
    1218             :         }
    1219       21247 :         need = minss(BQ.KC, nrelsup);
    1220             :     }
    1221             :   }
    1222      104509 :   while (need);
    1223             :   /* DONE */
    1224       70357 :   if (!quad_be_honest(&BQ)) goto START;
    1225       70356 :   if (DEBUGLEVEL>2) timer_printf(&T, "be honest");
    1226       70356 :   clearhash(BQ.hashtab);
    1227       70363 :   free_GRHcheck(&GRHcheck);
    1228             : 
    1229       70363 :   gen = get_clgp(&BQ,W,&cyc);
    1230       70361 :   gunclone(BQ.subFB);
    1231       70362 :   gunclone(BQ.powsubFB);
    1232       70363 :   if (BQ.PRECREG)
    1233        2156 :     return mkvec5(h, cyc, gen, real_i(R), mpodd(imag_i(R)) ? gen_m1:gen_1);
    1234             :   else
    1235       68207 :     return mkvec4(h, cyc, gen, real_i(R));
    1236             : }
    1237             : GEN
    1238        4424 : Buchquad(GEN D, double c, double c2, long prec)
    1239             : {
    1240        4424 :   pari_sp av = avma;
    1241        4424 :   GEN z = Buchquad_i(D, c, c2, prec);
    1242        4424 :   return gerepilecopy(av, z);
    1243             : }
    1244             : 
    1245             : GEN
    1246           0 : buchimag(GEN D, GEN c, GEN c2, GEN REL)
    1247           0 : { (void)REL; return Buchquad(D,gtodouble(c),gtodouble(c2),0); }
    1248             : 
    1249             : GEN
    1250           0 : buchreal(GEN D, GEN flag, GEN c, GEN c2, GEN REL, long prec) {
    1251           0 :   if (signe(flag)) pari_err_IMPL("narrow class group");
    1252           0 :   (void)REL; return Buchquad(D,gtodouble(c),gtodouble(c2),prec);
    1253             : }
    1254             : 
    1255             : GEN
    1256        4424 : quadclassunit0(GEN x, long flag, GEN data, long prec)
    1257             : {
    1258             :   long lx;
    1259        4424 :   double c1 = 0.0, c2 = 0.0;
    1260             : 
    1261        4424 :   if (!data) lx=1;
    1262             :   else
    1263             :   {
    1264          28 :     lx = lg(data);
    1265          28 :     if (typ(data)!=t_VEC) pari_err_TYPE("quadclassunit", data);
    1266          28 :     if (lx > 7) pari_err_DIM("quadclassunit [tech vector]");
    1267          28 :     if (lx > 3) lx = 3;
    1268             :   }
    1269        4424 :   switch(lx)
    1270             :   {
    1271          21 :     case 3: c2 = gtodouble(gel(data,2));
    1272          28 :     case 2: c1 = gtodouble(gel(data,1));
    1273             :   }
    1274        4424 :   if (flag) pari_err_IMPL("narrow class group");
    1275        4424 :   return Buchquad(x,c1,c2,prec);
    1276             : }
    1277             : GEN
    1278       61173 : quadclassno(GEN D)
    1279             : {
    1280       61173 :   pari_sp av = avma;
    1281       61173 :   GEN h = abgrp_get_no(Buchquad_i(D, 0, 0, 0));
    1282       61175 :   return icopy_avma(h, av);
    1283             : }
    1284             : long
    1285        4938 : quadclassnos(long D)
    1286             : {
    1287        4938 :   pari_sp av = avma;
    1288        4938 :   long h = itos(abgrp_get_no(Buchquad_i(stoi(D), 0, 0, 0)));
    1289        4938 :   return gc_long(av, h);
    1290             : }

Generated by: LCOV version 1.16