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 - buch2.c (source / functions) Hit Total Coverage
Test: PARI/GP v2.12.1 lcov report (development 25819-e703fe1174) Lines: 2178 2369 91.9 %
Date: 2020-09-18 06:10:04 Functions: 151 163 92.6 %
Legend: Lines: hit not hit

          Line data    Source code
       1             : /* Copyright (C) 2000  The PARI group.
       2             : 
       3             : This file is part of the PARI/GP package.
       4             : 
       5             : PARI/GP is free software; you can redistribute it and/or modify it under the
       6             : terms of the GNU General Public License as published by the Free Software
       7             : Foundation. It is distributed in the hope that it will be useful, but WITHOUT
       8             : ANY WARRANTY WHATSOEVER.
       9             : 
      10             : Check the License for details. You should have received a copy of it, along
      11             : with the package; see the file 'COPYING'. If not, write to the Free Software
      12             : Foundation, Inc., 51 Franklin Street, Fifth Floor, Boston, MA 02110-1301 USA. */
      13             : #include "pari.h"
      14             : #include "paripriv.h"
      15             : /*******************************************************************/
      16             : /*                                                                 */
      17             : /*         CLASS GROUP AND REGULATOR (McCURLEY, BUCHMANN)          */
      18             : /*                    GENERAL NUMBER FIELDS                        */
      19             : /*                                                                 */
      20             : /*******************************************************************/
      21             : /* get_random_ideal */
      22             : static const long RANDOM_BITS = 4;
      23             : /* Buchall */
      24             : static const long RELSUP = 5;
      25             : static const long FAIL_DIVISOR = 32;
      26             : static const long MINFAIL = 10;
      27             : /* small_norm */
      28             : static const long BNF_RELPID = 4;
      29             : static const long BMULT = 8;
      30             : static const long maxtry_ELEMENT = 1000*1000;
      31             : static const long maxtry_DEP = 20;
      32             : static const long maxtry_FACT = 500;
      33             : /* rnd_rel */
      34             : static const long RND_REL_RELPID = 1;
      35             : /* random relations */
      36             : static const long MINSFB = 3;
      37             : static const long SFB_MAX = 3;
      38             : static const long DEPSIZESFBMULT = 16;
      39             : static const long DEPSFBDIV = 10;
      40             : /* add_rel_i */
      41             : static const ulong mod_p = 27449UL;
      42             : /* be_honest */
      43             : static const long maxtry_HONEST = 50;
      44             : 
      45             : typedef struct FACT {
      46             :     long pr, ex;
      47             : } FACT;
      48             : 
      49             : typedef struct subFB_t {
      50             :   GEN subFB;
      51             :   struct subFB_t *old;
      52             : } subFB_t;
      53             : 
      54             : /* a factor base contains only non-inert primes
      55             :  * KC = # of P in factor base (p <= n, NP <= n2)
      56             :  * KC2= # of P assumed to generate class group (NP <= n2)
      57             :  *
      58             :  * KCZ = # of rational primes under ideals counted by KC
      59             :  * KCZ2= same for KC2 */
      60             : 
      61             : typedef struct FB_t {
      62             :   GEN FB; /* FB[i] = i-th rational prime used in factor base */
      63             :   GEN LP; /* vector of all prime ideals in FB */
      64             :   GEN *LV; /* LV[p] = vector of P|p, NP <= n2
      65             :             * isclone() is set for LV[p] iff all P|p are in FB
      66             :             * LV[i], i not prime or i > n2, is undefined! */
      67             :   GEN iLP; /* iLP[p] = i such that LV[p] = [LP[i],...] */
      68             :   GEN L_jid; /* indexes of "useful" prime ideals for rnd_rel */
      69             :   long KC, KCZ, KCZ2;
      70             :   GEN subFB; /* LP o subFB =  part of FB used to build random relations */
      71             :   int sfb_chg; /* need to change subFB ? */
      72             :   GEN perm; /* permutation of LP used to represent relations [updated by
      73             :                hnfspec/hnfadd: dense rows come first] */
      74             :   GEN idealperm; /* permutation of ideals under field automorphisms */
      75             :   GEN minidx; /* minidx[i] min ideal in orbit of LP[i] under field autom */
      76             :   subFB_t *allsubFB; /* all subFB's used */
      77             :   GEN embperm; /* permutations of the complex embeddings */
      78             :   long MAXDEPSIZESFB; /* # trials before increasing subFB */
      79             :   long MAXDEPSFB; /* MAXDEPSIZESFB / DEPSFBDIV, # trials befor rotating subFB */
      80             : } FB_t;
      81             : 
      82             : enum { sfb_CHANGE = 1, sfb_INCREASE = 2 };
      83             : 
      84             : typedef struct REL_t {
      85             :   GEN R; /* relation vector as t_VECSMALL; clone */
      86             :   long nz; /* index of first non-zero elt in R (hash) */
      87             :   GEN m; /* pseudo-minimum yielding the relation; clone */
      88             :   long relorig; /* relation this one is an image of */
      89             :   long relaut; /* automorphim used to compute this relation from the original */
      90             :   GEN emb; /* archimedean embeddings */
      91             :   GEN junk[2]; /*make sure sizeof(struct) is a power of two.*/
      92             : } REL_t;
      93             : 
      94             : typedef struct RELCACHE_t {
      95             :   REL_t *chk; /* last checkpoint */
      96             :   REL_t *base; /* first rel found */
      97             :   REL_t *last; /* last rel found so far */
      98             :   REL_t *end; /* target for last relation. base <= last <= end */
      99             :   size_t len; /* number of rels pre-allocated in base */
     100             :   long relsup; /* how many linearly dependent relations we allow */
     101             :   GEN basis; /* mod p basis (generating family actually) */
     102             :   ulong missing; /* missing vectors in generating family above */
     103             : } RELCACHE_t;
     104             : 
     105             : typedef struct FP_t {
     106             :   double **q;
     107             :   GEN x;
     108             :   double *y;
     109             :   double *z;
     110             :   double *v;
     111             : } FP_t;
     112             : 
     113             : typedef struct RNDREL_t {
     114             :   long jid;
     115             :   GEN ex;
     116             : } RNDREL_t;
     117             : 
     118             : static void
     119           0 : wr_rel(GEN e)
     120             : {
     121           0 :   long i, l = lg(e);
     122           0 :   for (i = 1; i < l; i++)
     123           0 :     if (e[i]) err_printf("%ld^%ld ",i,e[i]);
     124           0 : }
     125             : static void
     126           0 : dbg_newrel(RELCACHE_t *cache)
     127             : {
     128           0 :   if (DEBUGLEVEL > 1)
     129             :   {
     130           0 :     err_printf("\n++++ cglob = %ld\nrel = ", cache->last - cache->base);
     131           0 :     wr_rel(cache->last->R);
     132           0 :     err_printf("\n");
     133             :   }
     134             :   else
     135           0 :     err_printf("%ld ", cache->last - cache->base);
     136           0 : }
     137             : 
     138             : static void
     139       12082 : delete_cache(RELCACHE_t *M)
     140             : {
     141             :   REL_t *rel;
     142      225547 :   for (rel = M->base+1; rel <= M->last; rel++)
     143             :   {
     144      213465 :     gunclone(rel->R);
     145      213465 :     if (rel->m) gunclone(rel->m);
     146             :   }
     147       12082 :   pari_free((void*)M->base); M->base = NULL;
     148       12082 : }
     149             : 
     150             : static void
     151       12705 : delete_FB(FB_t *F)
     152             : {
     153             :   subFB_t *s, *sold;
     154       26199 :   for (s = F->allsubFB; s; s = sold) { sold = s->old; pari_free(s); }
     155       12705 :   gunclone(F->minidx);
     156       12705 :   gunclone(F->idealperm);
     157       12705 : }
     158             : 
     159             : static void
     160       12124 : reallocate(RELCACHE_t *M, long len)
     161             : {
     162       12124 :   REL_t *old = M->base;
     163       12124 :   M->len = len;
     164       12124 :   pari_realloc_ip((void**)&M->base, (len+1) * sizeof(REL_t));
     165       12124 :   if (old)
     166             :   {
     167          42 :     size_t last = M->last - old, chk = M->chk - old, end = M->end - old;
     168          42 :     M->last = M->base + last;
     169          42 :     M->chk  = M->base + chk;
     170          42 :     M->end  = M->base + end;
     171             :   }
     172       12124 : }
     173             : 
     174             : #define pr_get_smallp(pr) gel(pr,1)[2]
     175             : 
     176             : /* don't take P|p all other Q|p are already there */
     177             : static int
     178       53564 : bad_subFB(FB_t *F, long t)
     179             : {
     180       53564 :   GEN LP, P = gel(F->LP,t);
     181       53564 :   long p = pr_get_smallp(P);
     182       53564 :   LP = F->LV[p];
     183       53564 :   return (isclone(LP) && t == F->iLP[p] + lg(LP)-1);
     184             : }
     185             : 
     186             : static void
     187       13494 : assign_subFB(FB_t *F, GEN yes, long iyes)
     188             : {
     189       13494 :   long i, lv = sizeof(subFB_t) + iyes*sizeof(long); /* for struct + GEN */
     190       13494 :   subFB_t *s = (subFB_t *)pari_malloc(lv);
     191       13494 :   s->subFB = (GEN)&s[1];
     192       13494 :   s->old = F->allsubFB; F->allsubFB = s;
     193       54311 :   for (i = 0; i < iyes; i++) s->subFB[i] = yes[i];
     194       13494 :   F->subFB = s->subFB;
     195       13494 :   F->MAXDEPSIZESFB = (iyes-1) * DEPSIZESFBMULT;
     196       13494 :   F->MAXDEPSFB = F->MAXDEPSIZESFB / DEPSFBDIV;
     197       13494 : }
     198             : 
     199             : /* Determine the permutation of the ideals made by each field automorphism */
     200             : static GEN
     201       12705 : FB_aut_perm(FB_t *F, GEN auts, GEN cyclic)
     202             : {
     203       12705 :   long i, j, m, KC = F->KC, nauts = lg(auts)-1;
     204       12705 :   GEN minidx, perm = zero_Flm_copy(KC, nauts);
     205             : 
     206       12705 :   if (!nauts) { F->minidx = gclone(identity_zv(KC)); return cgetg(1,t_MAT); }
     207       12173 :   minidx = zero_Flv(KC);
     208       26404 :   for (m = 1; m < lg(cyclic); m++)
     209             :   {
     210       14231 :     GEN thiscyc = gel(cyclic, m);
     211       14231 :     long k0 = thiscyc[1];
     212       14231 :     GEN aut = gel(auts, k0), permk0 = gel(perm, k0), ppermk;
     213       14231 :     i = 1;
     214       75278 :     while (i <= KC)
     215             :     {
     216       61047 :       pari_sp av2 = avma;
     217       61047 :       GEN seen = zero_Flv(KC), P = gel(F->LP, i);
     218       61047 :       long imin = i, p, f, l;
     219       61047 :       p = pr_get_smallp(P);
     220       61047 :       f = pr_get_f(P);
     221             :       do
     222             :       {
     223      179473 :         if (++i > KC) break;
     224      165242 :         P = gel(F->LP, i);
     225             :       }
     226      165242 :       while (p == pr_get_smallp(P) && f == pr_get_f(P));
     227      240520 :       for (j = imin; j < i; j++)
     228             :       {
     229      179473 :         GEN img = ZM_ZC_mul(aut, pr_get_gen(gel(F->LP, j)));
     230      608559 :         for (l = imin; l < i; l++)
     231      608559 :           if (!seen[l] && ZC_prdvd(img, gel(F->LP, l)))
     232             :           {
     233      179473 :             seen[l] = 1; permk0[j] = l; break;
     234             :           }
     235             :       }
     236       61047 :       set_avma(av2);
     237             :     }
     238       15897 :     for (ppermk = permk0, i = 2; i < lg(thiscyc); i++)
     239             :     {
     240        1666 :       GEN permk = gel(perm, thiscyc[i]);
     241      111937 :       for (j = 1; j <= KC; j++) permk[j] = permk0[ppermk[j]];
     242        1666 :       ppermk = permk;
     243             :     }
     244             :   }
     245      111013 :   for (j = 1; j <= KC; j++)
     246             :   {
     247       98840 :     if (minidx[j]) continue;
     248       46851 :     minidx[j] = j;
     249      120428 :     for (i = 1; i <= nauts; i++) minidx[coeff(perm, j, i)] = j;
     250             :   }
     251       12173 :   F->minidx = gclone(minidx); return perm;
     252             : }
     253             : 
     254             : /* set subFB.
     255             :  * Fill F->perm (if != NULL): primes ideals sorted by increasing norm (except
     256             :  * the ones in subFB come first [dense rows for hnfspec]) */
     257             : static void
     258       12705 : subFBgen(FB_t *F, GEN auts, GEN cyclic, double PROD, long minsFB)
     259             : {
     260             :   GEN y, perm, yes, no;
     261       12705 :   long i, j, k, iyes, ino, lv = F->KC + 1;
     262             :   double prod;
     263             :   pari_sp av;
     264             : 
     265       12705 :   F->LP   = cgetg(lv, t_VEC);
     266       12705 :   F->L_jid = F->perm = cgetg(lv, t_VECSMALL);
     267       12705 :   av = avma;
     268       12705 :   y = cgetg(lv,t_COL); /* Norm P */
     269       72037 :   for (k=0, i=1; i <= F->KCZ; i++)
     270             :   {
     271       59332 :     GEN LP = F->LV[F->FB[i]];
     272       59332 :     long l = lg(LP);
     273      179438 :     for (j = 1; j < l; j++)
     274             :     {
     275      120106 :       GEN P = gel(LP,j);
     276      120106 :       k++;
     277      120106 :       gel(y,k) = pr_norm(P);
     278      120106 :       gel(F->LP,k) = P;
     279             :     }
     280             :   }
     281             :   /* perm sorts LP by increasing norm */
     282       12705 :   perm = indexsort(y);
     283       12705 :   no  = cgetg(lv, t_VECSMALL); ino  = 1;
     284       12705 :   yes = cgetg(lv, t_VECSMALL); iyes = 1;
     285       12705 :   prod = 1.0;
     286       62622 :   for (i = 1; i < lv; i++)
     287             :   {
     288       53564 :     long t = perm[i];
     289       53564 :     if (bad_subFB(F, t)) { no[ino++] = t; continue; }
     290             : 
     291       24101 :     yes[iyes++] = t;
     292       24101 :     prod *= (double)itos(gel(y,t));
     293       24101 :     if (iyes > minsFB && prod > PROD) break;
     294             :   }
     295       12705 :   setlg(yes, iyes);
     296       36806 :   for (j=1; j<iyes; j++)     F->perm[j] = yes[j];
     297       42168 :   for (i=1; i<ino; i++, j++) F->perm[j] =  no[i];
     298       79247 :   for (   ; j<lv; j++)       F->perm[j] =  perm[j];
     299       12705 :   F->allsubFB = NULL;
     300       12705 :   F->idealperm = gclone(FB_aut_perm(F, auts, cyclic));
     301       12705 :   if (iyes) assign_subFB(F, yes, iyes);
     302       12705 :   set_avma(av);
     303       12705 : }
     304             : static int
     305        3166 : subFB_change(FB_t *F)
     306             : {
     307        3166 :   long i, iyes, minsFB, lv = F->KC + 1, l = lg(F->subFB)-1;
     308        3166 :   pari_sp av = avma;
     309        3166 :   GEN yes, L_jid = F->L_jid, present = zero_zv(lv-1);
     310             : 
     311        3166 :   switch (F->sfb_chg)
     312             :   {
     313          99 :     case sfb_INCREASE: minsFB = l + 1; break;
     314        3067 :     default: minsFB = l; break;
     315             :   }
     316             : 
     317        3166 :   yes = cgetg(minsFB+1, t_VECSMALL); iyes = 1;
     318        3166 :   if (L_jid)
     319             :   {
     320       12485 :     for (i = 1; i < lg(L_jid); i++)
     321             :     {
     322       11631 :       long l = L_jid[i];
     323       11631 :       yes[iyes++] = l;
     324       11631 :       present[l] = 1;
     325       11631 :       if (iyes > minsFB) break;
     326             :     }
     327             :   }
     328           0 :   else i = 1;
     329        3166 :   if (iyes <= minsFB)
     330             :   {
     331        2099 :     for ( ; i < lv; i++)
     332             :     {
     333        2099 :       long l = F->perm[i];
     334        2099 :       if (present[l]) continue;
     335        2099 :       yes[iyes++] = l;
     336        2099 :       if (iyes > minsFB) break;
     337             :     }
     338         854 :     if (i == lv) return 0;
     339             :   }
     340        3166 :   if (zv_equal(F->subFB, yes))
     341             :   {
     342        2377 :     if (DEBUGLEVEL) err_printf("\n*** NOT Changing sub factor base\n");
     343             :   }
     344             :   else
     345             :   {
     346         789 :     if (DEBUGLEVEL) err_printf("\n*** Changing sub factor base\n");
     347         789 :     assign_subFB(F, yes, iyes);
     348             :   }
     349        3166 :   F->sfb_chg = 0; return gc_bool(av, 1);
     350             : }
     351             : 
     352             : /* make sure enough room to store n more relations */
     353             : static void
     354       78163 : pre_allocate(RELCACHE_t *cache, size_t n)
     355             : {
     356       78163 :   size_t len = (cache->last - cache->base) + n;
     357       78163 :   if (len >= cache->len) reallocate(cache, len << 1);
     358       78163 : }
     359             : 
     360             : void
     361       28433 : init_GRHcheck(GRHcheck_t *S, long N, long R1, double LOGD)
     362             : {
     363       28433 :   const double c1 = M_PI*M_PI/2;
     364       28433 :   const double c2 = 3.663862376709;
     365       28433 :   const double c3 = 3.801387092431; /* Euler + log(8*Pi)*/
     366       28433 :   S->clone = 0;
     367       28433 :   S->cN = R1*c2 + N*c1;
     368       28433 :   S->cD = LOGD - N*c3 - R1*M_PI/2;
     369       28433 :   S->maxprimes = 16000; /* sufficient for LIMC=176081*/
     370       28433 :   S->primes = (GRHprime_t*)pari_malloc(S->maxprimes*sizeof(*S->primes));
     371       28433 :   S->nprimes = 0;
     372       28433 :   S->limp = 0;
     373       28433 :   u_forprime_init(&S->P, 2, ULONG_MAX);
     374       28433 : }
     375             : 
     376             : void
     377       28433 : free_GRHcheck(GRHcheck_t *S)
     378             : {
     379       28433 :   if (S->clone)
     380             :   {
     381       12040 :     long i = S->nprimes;
     382             :     GRHprime_t *pr;
     383     1505385 :     for (pr = S->primes, i = S->nprimes; i > 0; pr++, i--) gunclone(pr->dec);
     384             :   }
     385       28433 :   pari_free(S->primes);
     386       28433 : }
     387             : 
     388             : int
     389      310961 : GRHok(GRHcheck_t *S, double L, double SA, double SB)
     390             : {
     391      310961 :   return (S->cD + (S->cN + 2*SB) / L - 2*SA < -1e-8);
     392             : }
     393             : 
     394             : /* Return factorization pattern of p: [f,n], where n[i] primes of
     395             :  * residue degree f[i] */
     396             : static GEN
     397     1493345 : get_fs(GEN nf, GEN P, GEN index, ulong p)
     398             : {
     399             :   long j, k, f, n, l;
     400             :   GEN fs, ns;
     401             : 
     402     1493345 :   if (umodiu(index, p))
     403             :   { /* easy case: p does not divide index */
     404     1490013 :     GEN F = Flx_degfact(ZX_to_Flx(P,p), p);
     405     1490013 :     fs = gel(F,1); l = lg(fs);
     406             :   }
     407             :   else
     408             :   {
     409        3332 :     GEN F = idealprimedec(nf, utoipos(p));
     410        3332 :     l = lg(F);
     411        3332 :     fs = cgetg(l, t_VECSMALL);
     412       13734 :     for (j = 1; j < l; j++) fs[j] = pr_get_f(gel(F,j));
     413             :   }
     414     1493345 :   ns = cgetg(l, t_VECSMALL);
     415     1493345 :   f = fs[1]; n = 1;
     416     2702679 :   for (j = 2, k = 1; j < l; j++)
     417     1209334 :     if (fs[j] == f)
     418     1112622 :       n++;
     419             :     else
     420             :     {
     421       96712 :       ns[k] = n; fs[k] = f; k++;
     422       96712 :       f = fs[j]; n = 1;
     423             :     }
     424     1493345 :   ns[k] = n; fs[k] = f; k++;
     425     1493345 :   setlg(fs, k);
     426     1493345 :   setlg(ns, k); return mkvec2(fs,ns);
     427             : }
     428             : 
     429             : /* cache data for all rational primes up to the LIM */
     430             : static void
     431      164416 : cache_prime_dec(GRHcheck_t *S, ulong LIM, GEN nf)
     432             : {
     433      164416 :   pari_sp av = avma;
     434             :   GRHprime_t *pr;
     435             :   GEN index, P;
     436             :   double nb;
     437             : 
     438      164416 :   if (S->limp >= LIM) return;
     439       59927 :   S->clone = 1;
     440       59927 :   nb = primepi_upper_bound((double)LIM); /* #{p <= LIM} <= nb */
     441       59927 :   GRH_ensure(S, nb+1); /* room for one extra prime */
     442       59927 :   P = nf_get_pol(nf);
     443       59927 :   index = nf_get_index(nf);
     444       59927 :   for (pr = S->primes + S->nprimes;;)
     445     1433418 :   {
     446     1493345 :     ulong p = u_forprime_next(&(S->P));
     447     1493345 :     pr->p = p;
     448     1493345 :     pr->logp = log((double)p);
     449     1493345 :     pr->dec = gclone(get_fs(nf, P, index, p));
     450     1493345 :     S->nprimes++;
     451     1493345 :     pr++;
     452     1493345 :     set_avma(av);
     453             :     /* store up to nextprime(LIM) included */
     454     1493345 :     if (p >= LIM) { S->limp = p; break; }
     455             :   }
     456             : }
     457             : 
     458             : static double
     459      424242 : tailresback(long R1, long R2, double rK, long C, double C2, double C3, double r1K, double r2K, double logC, double logC2, double logC3)
     460             : {
     461      424242 :   const double  rQ = 1.83787706641;
     462      424242 :   const double r1Q = 1.98505372441;
     463      424242 :   const double r2Q = 1.07991541347;
     464     1272726 :   return fabs((R1+R2-1)*(12*logC3+4*logC2-9*logC-6)/(2*C*logC3)
     465      424242 :          + (rK-rQ)*(6*logC2 + 5*logC + 2)/(C*logC3)
     466      424242 :          - R2*(6*logC2+11*logC+6)/(C2*logC2)
     467      424242 :          - 2*(r1K-r1Q)*(3*logC2 + 4*logC + 2)/(C2*logC3)
     468      424242 :          + (R1+R2-1)*(12*logC3+40*logC2+45*logC+18)/(6*C3*logC3)
     469      424242 :          + (r2K-r2Q)*(2*logC2 + 3*logC + 2)/(C3*logC3));
     470             : }
     471             : 
     472             : static double
     473      212121 : tailres(long R1, long R2, double al2K, double rKm, double rKM, double r1Km,
     474             :         double r1KM, double r2Km, double r2KM, double C, long i)
     475             : {
     476             :   /* C >= 3*2^i, lower bound for eint1(log(C)/2) */
     477             :   /* for(i=0,30,print(eint1(log(3*2^i)/2))) */
     478             :   static double tab[] = {
     479             :     0.50409264803,
     480             :     0.26205336997,
     481             :     0.14815491171,
     482             :     0.08770540561,
     483             :     0.05347651832,
     484             :     0.03328934284,
     485             :     0.02104510690,
     486             :     0.01346475900,
     487             :     0.00869778586,
     488             :     0.00566279855,
     489             :     0.00371111950,
     490             :     0.00244567837,
     491             :     0.00161948049,
     492             :     0.00107686891,
     493             :     0.00071868750,
     494             :     0.00048119961,
     495             :     0.00032312188,
     496             :     0.00021753772,
     497             :     0.00014679818,
     498             :     9.9272855581E-5,
     499             :     6.7263969995E-5,
     500             :     4.5656812967E-5,
     501             :     3.1041124593E-5,
     502             :     2.1136011590E-5,
     503             :     1.4411645381E-5,
     504             :     9.8393304088E-6,
     505             :     6.7257395409E-6,
     506             :     4.6025878272E-6,
     507             :     3.1529719271E-6,
     508             :     2.1620490021E-6,
     509             :     1.4839266071E-6
     510             :   };
     511      212121 :   const double logC = log(C), logC2 = logC*logC, logC3 = logC*logC2;
     512      212121 :   const double C2 = C*C, C3 = C*C2;
     513      212121 :   double E1 = i >30? 0: tab[i];
     514      212121 :   return al2K*((33*logC2+22*logC+8)/(8*logC3*sqrt(C))+15*E1/16)
     515      212121 :     + maxdd(tailresback(rKm,r1KM,r2Km, C,C2,C3,R1,R2,logC,logC2,logC3),
     516      212121 :             tailresback(rKM,r1Km,r2KM, C,C2,C3,R1,R2,logC,logC2,logC3))/2
     517      212121 :     + ((R1+R2-1)*4*C+R2)*(C2+6*logC)/(4*C2*C2*logC2);
     518             : }
     519             : 
     520             : static long
     521       12040 : primeneeded(long N, long R1, long R2, double LOGD)
     522             : {
     523       12040 :   const double lim = 0.25; /* should be log(2)/2 == 0.34657... */
     524       12040 :   const double al2K =  0.3526*LOGD - 0.8212*N + 4.5007;
     525       12040 :   const double  rKm = -1.0155*LOGD + 2.1042*N - 8.3419;
     526       12040 :   const double  rKM = -0.5   *LOGD + 1.2076*N + 1;
     527       12040 :   const double r1Km = -       LOGD + 1.4150*N;
     528       12040 :   const double r1KM = -       LOGD + 1.9851*N;
     529       12040 :   const double r2Km = -       LOGD + 0.9151*N;
     530       12040 :   const double r2KM = -       LOGD + 1.0800*N;
     531       12040 :   long Cmin = 3, Cmax = 3, i = 0;
     532      107667 :   while (tailres(R1, R2, al2K, rKm, rKM, r1Km, r1KM, r2Km, r2KM, Cmax, i) > lim)
     533             :   {
     534       95627 :     Cmin = Cmax;
     535       95627 :     Cmax *= 2;
     536       95627 :     i++;
     537             :   }
     538       12040 :   i--;
     539      116494 :   while (Cmax - Cmin > 1)
     540             :   {
     541      104454 :     long t = (Cmin + Cmax)/2;
     542      104454 :     if (tailres(R1, R2, al2K, rKm, rKM, r1Km, r1KM, r2Km, r2KM, t, i) > lim)
     543       69405 :       Cmin = t;
     544             :     else
     545       35049 :       Cmax = t;
     546             :   }
     547       12040 :   return Cmax;
     548             : }
     549             : 
     550             : /* ~ 1 / Res(s = 1, zeta_K) */
     551             : static GEN
     552       12040 : compute_invres(GRHcheck_t *S, long LIMC)
     553             : {
     554       12040 :   pari_sp av = avma;
     555       12040 :   double loginvres = 0.;
     556             :   GRHprime_t *pr;
     557             :   long i;
     558       12040 :   double logLIMC = log((double)LIMC);
     559       12040 :   double logLIMC2 = logLIMC*logLIMC, denc;
     560             :   double c0, c1, c2;
     561       12040 :   denc = 1/(pow((double)LIMC, 3.) * logLIMC * logLIMC2);
     562       12040 :   c2 = (    logLIMC2 + 3 * logLIMC / 2 + 1) * denc;
     563       12040 :   denc *= LIMC;
     564       12040 :   c1 = (3 * logLIMC2 + 4 * logLIMC     + 2) * denc;
     565       12040 :   denc *= LIMC;
     566       12040 :   c0 = (3 * logLIMC2 + 5 * logLIMC / 2 + 1) * denc;
     567     1494493 :   for (pr = S->primes, i = S->nprimes; i > 0; pr++, i--)
     568             :   {
     569             :     GEN dec, fs, ns;
     570             :     long addpsi;
     571             :     double addpsi1, addpsi2;
     572     1493345 :     double logp = pr->logp, NPk;
     573     1493345 :     long j, k, limp = logLIMC/logp;
     574     1493345 :     ulong p = pr->p, p2 = p*p;
     575     1493345 :     if (limp < 1) break;
     576     1482453 :     dec = pr->dec;
     577     1482453 :     fs = gel(dec, 1); ns = gel(dec, 2);
     578     1482453 :     loginvres += 1./p;
     579             :     /* NB: limp = 1 nearly always and limp > 2 for very few primes */
     580     1741789 :     for (k=2, NPk = p; k <= limp; k++) { NPk *= p; loginvres += 1/(k * NPk); }
     581     1482453 :     addpsi = limp;
     582     1482453 :     addpsi1 = p *(pow((double)p , (double)limp)-1)/(p -1);
     583     1482453 :     addpsi2 = p2*(pow((double)p2, (double)limp)-1)/(p2-1);
     584     1482453 :     j = lg(fs);
     585     3061205 :     while (--j > 0)
     586             :     {
     587             :       long f, nb, kmax;
     588             :       double NP, NP2, addinvres;
     589     1578752 :       f = fs[j]; if (f > limp) continue;
     590      726418 :       nb = ns[j];
     591      726418 :       NP = pow((double)p, (double)f);
     592      726418 :       addinvres = 1/NP;
     593      726418 :       kmax = limp / f;
     594      893263 :       for (k=2, NPk = NP; k <= kmax; k++) { NPk *= NP; addinvres += 1/(k*NPk); }
     595      726418 :       NP2 = NP*NP;
     596      726418 :       loginvres -= nb * addinvres;
     597      726418 :       addpsi -= nb * f * kmax;
     598      726418 :       addpsi1 -= nb*(f*NP *(pow(NP ,(double)kmax)-1)/(NP -1));
     599      726418 :       addpsi2 -= nb*(f*NP2*(pow(NP2,(double)kmax)-1)/(NP2-1));
     600             :     }
     601     1482453 :     loginvres -= (addpsi*c0 - addpsi1*c1 + addpsi2*c2)*logp;
     602             :   }
     603       12040 :   return gerepileuptoleaf(av, mpexp(dbltor(loginvres)));
     604             : }
     605             : 
     606             : static long
     607       12040 : nthideal(GRHcheck_t *S, GEN nf, long n)
     608             : {
     609       12040 :   pari_sp av = avma;
     610       12040 :   GEN P = nf_get_pol(nf);
     611       12040 :   ulong p = 0, *vecN = (ulong*)const_vecsmall(n, LONG_MAX);
     612       12040 :   long i, N = poldegree(P, -1);
     613       12040 :   for (i = 0; ; i++)
     614       35112 :   {
     615             :     GRHprime_t *pr;
     616             :     GEN fs;
     617       47152 :     cache_prime_dec(S, p+1, nf);
     618       47152 :     pr = S->primes + i;
     619       47152 :     fs = gel(pr->dec, 1);
     620       47152 :     p = pr->p;
     621       47152 :     if (fs[1] != N)
     622             :     {
     623       32025 :       GEN ns = gel(pr->dec, 2);
     624       32025 :       long k, l, j = lg(fs);
     625       66808 :       while (--j > 0)
     626             :       {
     627       34783 :         ulong NP = upowuu(p, fs[j]);
     628             :         long nf;
     629       34783 :         if (!NP) continue;
     630      117755 :         for (k = 1; k <= n; k++) if (vecN[k] > NP) break;
     631       34473 :         if (k > n) continue;
     632             :         /* vecN[k] <= NP */
     633       21923 :         nf = ns[j]; /*#{primes of norme NP} = nf, insert them here*/
     634       45890 :         for (l = k+nf; l <= n; l++) vecN[l] = vecN[l-nf];
     635       56046 :         for (l = 0; l < nf && k+l <= n; l++) vecN[k+l] = NP;
     636       48132 :         while (l <= k) vecN[l++] = NP;
     637             :       }
     638             :     }
     639       47152 :     if (p > vecN[n]) break;
     640             :   }
     641       12040 :   return gc_long(av, vecN[n]);
     642             : }
     643             : 
     644             : /* Compute FB, LV, iLP + KC*. Reset perm
     645             :  * C2: bound for norm of tested prime ideals (includes be_honest())
     646             :  * C1: bound for p, such that P|p (NP <= C2) used to build relations */
     647             : static void
     648       12705 : FBgen(FB_t *F, GEN nf, long N, ulong C1, ulong C2, GRHcheck_t *S)
     649             : {
     650             :   GRHprime_t *pr;
     651             :   long i, ip;
     652             :   GEN prim;
     653       12705 :   const double L = log((double)C2 + 0.5);
     654             : 
     655       12705 :   cache_prime_dec(S, C2, nf);
     656       12705 :   pr = S->primes;
     657       12705 :   F->sfb_chg = 0;
     658       12705 :   F->FB  = cgetg(C2+1, t_VECSMALL);
     659       12705 :   F->iLP = cgetg(C2+1, t_VECSMALL);
     660       12705 :   F->LV = (GEN*)const_vec(C2, NULL);
     661             : 
     662       12705 :   prim = icopy(gen_1);
     663       12705 :   i = ip = 0;
     664       12705 :   F->KC = F->KCZ = 0;
     665      119455 :   for (;; pr++) /* p <= C2 */
     666      119455 :   {
     667      132160 :     ulong p = pr->p;
     668             :     long k, l, m;
     669             :     GEN LP, nb, f;
     670             : 
     671      132160 :     if (!F->KC && p > C1) { F->KCZ = i; F->KC = ip; }
     672      132160 :     if (p > C2) break;
     673             : 
     674      125083 :     if (DEBUGLEVEL>1) err_printf(" %ld",p);
     675             : 
     676      125083 :     f = gel(pr->dec, 1); nb = gel(pr->dec, 2);
     677      125083 :     if (f[1] == N)
     678             :     {
     679       38836 :       if (p == C2) break;
     680       36813 :       continue; /* p inert */
     681             :     }
     682       86247 :     l = (long)(L/pr->logp); /* p^f <= C2  <=> f <= l */
     683      149121 :     for (k=0, m=1; m < lg(f) && f[m]<=l; m++) k += nb[m];
     684       86247 :     if (!k)
     685             :     { /* too inert to appear in FB */
     686       24262 :       if (p == C2) break;
     687       24150 :       continue;
     688             :     }
     689       61985 :     prim[2] = p; LP = idealprimedec_limit_f(nf,prim, l);
     690             :     /* keep non-inert ideals with Norm <= C2 */
     691       61985 :     if (m == lg(f)) setisclone(LP); /* flag it: all prime divisors in FB */
     692       61985 :     F->FB[++i]= p;
     693       61985 :     F->LV[p]  = LP;
     694       61985 :     F->iLP[p] = ip; ip += k;
     695       61985 :     if (p == C2) break;
     696             :   }
     697       12705 :   if (!F->KC) { F->KCZ = i; F->KC = ip; }
     698             :   /* Note F->KC > 0 otherwise GRHchk is false */
     699       12705 :   setlg(F->FB, F->KCZ+1); F->KCZ2 = i;
     700       12705 :   if (DEBUGLEVEL>1)
     701             :   {
     702           0 :     err_printf("\n");
     703           0 :     if (DEBUGLEVEL>6)
     704             :     {
     705           0 :       err_printf("########## FACTORBASE ##########\n\n");
     706           0 :       err_printf("KC2=%ld, KC=%ld, KCZ=%ld, KCZ2=%ld\n",
     707             :                   ip, F->KC, F->KCZ, F->KCZ2);
     708           0 :       for (i=1; i<=F->KCZ; i++) err_printf("++ LV[%ld] = %Ps",i,F->LV[F->FB[i]]);
     709             :     }
     710             :   }
     711       12705 :   F->perm = NULL; F->L_jid = NULL;
     712       12705 : }
     713             : 
     714             : static int
     715       92519 : GRHchk(GEN nf, GRHcheck_t *S, ulong LIMC)
     716             : {
     717       92519 :   double logC = log((double)LIMC), SA = 0, SB = 0;
     718       92519 :   GRHprime_t *pr = S->primes;
     719             : 
     720       92519 :   cache_prime_dec(S, LIMC, nf);
     721       92519 :   for (pr = S->primes;; pr++)
     722      913906 :   {
     723     1006425 :     ulong p = pr->p;
     724             :     GEN dec, fs, ns;
     725             :     double logCslogp;
     726             :     long j;
     727             : 
     728     1006425 :     if (p > LIMC) break;
     729      934381 :     dec = pr->dec; fs = gel(dec, 1); ns = gel(dec,2);
     730      934381 :     logCslogp = logC/pr->logp;
     731     1413027 :     for (j = 1; j < lg(fs); j++)
     732             :     {
     733     1090824 :       long f = fs[j], M, nb;
     734             :       double logNP, q, A, B;
     735     1090824 :       if (f > logCslogp) break;
     736      478646 :       logNP = f * pr->logp;
     737      478646 :       q = 1/sqrt((double)upowuu(p, f));
     738      478646 :       A = logNP * q; B = logNP * A; M = (long)(logCslogp/f);
     739      478646 :       if (M > 1)
     740             :       {
     741       81256 :         double inv1_q = 1 / (1-q);
     742       81256 :         A *= (1 - pow(q, (double)M)) * inv1_q;
     743       81256 :         B *= (1 - pow(q, (double)M)*(M+1 - M*q)) * inv1_q * inv1_q;
     744             :       }
     745      478646 :       nb = ns[j];
     746      478646 :       SA += nb * A;
     747      478646 :       SB += nb * B;
     748             :     }
     749      934381 :     if (p == LIMC) break;
     750             :   }
     751       92519 :   return GRHok(S, logC, SA, SB);
     752             : }
     753             : 
     754             : /*  SMOOTH IDEALS */
     755             : static void
     756     6737295 : store(long i, long e, FACT *fact)
     757             : {
     758     6737295 :   ++fact[0].pr;
     759     6737295 :   fact[fact[0].pr].pr = i; /* index */
     760     6737295 :   fact[fact[0].pr].ex = e; /* exponent */
     761     6737295 : }
     762             : 
     763             : /* divide out x by all P|p, where x as in can_factor().  k = v_p(Nx) */
     764             : static int
     765     3353452 : divide_p_elt(GEN LP, long ip, long k, GEN m, FACT *fact)
     766             : {
     767     3353452 :   long j, l = lg(LP);
     768    17955311 :   for (j=1; j<l; j++)
     769             :   {
     770    17946491 :     GEN P = gel(LP,j);
     771    17946491 :     long v = ZC_nfval(m, P);
     772    17946491 :     if (!v) continue;
     773     6259002 :     store(ip + j, v, fact); /* v = v_P(m) > 0 */
     774     6259002 :     k -= v * pr_get_f(P);
     775     6259002 :     if (!k) return 1;
     776             :   }
     777        8820 :   return 0;
     778             : }
     779             : static int
     780      105318 : divide_p_id(GEN LP, long ip, long k, GEN nf, GEN I, FACT *fact)
     781             : {
     782      105318 :   long j, l = lg(LP);
     783      160865 :   for (j=1; j<l; j++)
     784             :   {
     785      154250 :     GEN P = gel(LP,j);
     786      154250 :     long v = idealval(nf,I, P);
     787      154250 :     if (!v) continue;
     788       99624 :     store(ip + j, v, fact); /* v = v_P(I) > 0 */
     789       99624 :     k -= v * pr_get_f(P);
     790       99624 :     if (!k) return 1;
     791             :   }
     792        6615 :   return 0;
     793             : }
     794             : static int
     795      341144 : divide_p_quo(GEN LP, long ip, long k, GEN nf, GEN I, GEN m, FACT *fact)
     796             : {
     797      341144 :   long j, l = lg(LP);
     798      524741 :   for (j=1; j<l; j++)
     799             :   {
     800      524446 :     GEN P = gel(LP,j);
     801      524446 :     long v = ZC_nfval(m, P);
     802      524446 :     if (!v) continue;
     803      370413 :     v -= idealval(nf,I, P);
     804      370413 :     if (!v) continue;
     805      363998 :     store(ip + j, v, fact); /* v = v_P(m / I) > 0 */
     806      363998 :     k -= v * pr_get_f(P);
     807      363998 :     if (!k) return 1;
     808             :   }
     809         295 :   return 0;
     810             : }
     811             : 
     812             : /* |*N| != 0 is the norm of a primitive ideal, in particular not divisible by
     813             :  * any inert prime. Is |*N| a smooth rational integer wrt F ? (put the
     814             :  * exponents in *ex) */
     815             : static int
     816     4360288 : smooth_norm(FB_t *F, GEN *N, GEN *ex)
     817             : {
     818     4360288 :   GEN FB = F->FB;
     819     4360288 :   const long KCZ = F->KCZ;
     820     4360288 :   const ulong limp = uel(FB,KCZ); /* last p in FB */
     821             :   long i;
     822             : 
     823     4360288 :   *ex = new_chunk(KCZ+1);
     824     4360288 :   for (i=1; ; i++)
     825   240901366 :   {
     826             :     int stop;
     827   245261654 :     ulong p = uel(FB,i);
     828   245261654 :     long v = Z_lvalrem_stop(N, p, &stop);
     829   245261654 :     (*ex)[i] = v;
     830   245261654 :     if (v)
     831             :     {
     832     7985006 :       GEN LP = F->LV[p];
     833     7985006 :       if(!LP) pari_err_BUG("can_factor");
     834    10098520 :       if (lg(LP) == 1) return 0;
     835     7985006 :       if (stop) break;
     836             :     }
     837   243014880 :     if (i == KCZ) return 0;
     838             :   }
     839     2246774 :   (*ex)[0] = i;
     840     2246774 :   return (abscmpiu(*N,limp) <= 0);
     841             : }
     842             : 
     843             : static int
     844     3799914 : divide_p(FB_t *F, long p, long k, GEN nf, GEN I, GEN m, FACT *fact)
     845             : {
     846     3799914 :   GEN LP = F->LV[p];
     847     3799914 :   long ip = F->iLP[p];
     848     3799914 :   if (!m) return divide_p_id (LP,ip,k,nf,I,fact);
     849     3694596 :   if (!I) return divide_p_elt(LP,ip,k,m,fact);
     850      341144 :   return divide_p_quo(LP,ip,k,nf,I,m,fact);
     851             : }
     852             : 
     853             : /* Let x = m if I == NULL,
     854             :  *         I if m == NULL,
     855             :  *         m/I otherwise.
     856             :  * Can we factor the integral primitive ideal x ? |N| = Norm x > 0 */
     857             : static long
     858     4471844 : can_factor(FB_t *F, GEN nf, GEN I, GEN m, GEN N, FACT *fact)
     859             : {
     860             :   GEN ex;
     861     4471844 :   long i, res = 0;
     862     4471844 :   fact[0].pr = 0;
     863     4471844 :   if (is_pm1(N)) return 1;
     864     4360288 :   if (!smooth_norm(F, &N, &ex)) goto END;
     865    32520877 :   for (i=1; i<=ex[0]; i++)
     866    30600502 :     if (ex[i] && !divide_p(F, F->FB[i], ex[i], nf, I, m, fact)) goto END;
     867     1920375 :   res = is_pm1(N) || divide_p(F, itou(N), 1, nf, I, m, fact);
     868     4360288 : END:
     869     4360288 :   if (!res && DEBUGLEVEL > 1) err_printf(".");
     870     4360288 :   return res;
     871             : }
     872             : 
     873             : /* can we factor m/I ? [m in I from idealpseudomin_nonscalar], NI = norm I */
     874             : static long
     875     1850163 : factorgen(FB_t *F, GEN nf, GEN I, GEN NI, GEN m, FACT *fact)
     876             : {
     877     1850163 :   long e, r1 = nf_get_r1(nf);
     878     1850163 :   GEN M = nf_get_M(nf);
     879     1850163 :   GEN N = divri(embed_norm(RgM_RgC_mul(M,m), r1), NI); /* ~ N(m/I) */
     880     1850163 :   N = grndtoi(N, &e);
     881     1850163 :   if (e > -1)
     882             :   {
     883           0 :     if (DEBUGLEVEL > 1) err_printf("+");
     884           0 :     return 0;
     885             :   }
     886     1850163 :   return can_factor(F, nf, I, m, N, fact);
     887             : }
     888             : 
     889             : /*  FUNDAMENTAL UNITS */
     890             : 
     891             : /* a, m real. Return  (Re(x) + a) + I * (Im(x) % m) */
     892             : static GEN
     893     1436632 : addRe_modIm(GEN x, GEN a, GEN m)
     894             : {
     895             :   GEN re, im, z;
     896     1436632 :   if (typ(x) == t_COMPLEX)
     897             :   {
     898     1124229 :     im = modRr_safe(gel(x,2), m);
     899     1124229 :     if (!im) return NULL;
     900     1124228 :     re = gadd(gel(x,1), a);
     901     1124228 :     z = gequal0(im)? re: mkcomplex(re, im);
     902             :   }
     903             :   else
     904      312403 :     z = gadd(x, a);
     905     1436631 :   return z;
     906             : }
     907             : 
     908             : /* clean archimedean components */
     909             : static GEN
     910      595708 : cleanarch(GEN x, long N, long prec)
     911             : {
     912             :   long i, l, R1, RU;
     913      595708 :   GEN s, pi2, y = cgetg_copy(x, &l);
     914             : 
     915      595708 :   if (typ(x) == t_MAT)
     916             :   {
     917      144937 :     for (i = 1; i < l; i++)
     918      120702 :       if (!(gel(y,i) = cleanarch(gel(x,i), N, prec))) return NULL;
     919       24235 :     return y;
     920             :   }
     921      571472 :   RU = l-1; R1 = (RU<<1) - N; pi2 = Pi2n(1, prec);
     922      571472 :   s = gdivgs(RgV_sum(real_i(x)), -N); /* -log |norm(x)| / N */
     923     1598574 :   for (i = 1; i <= R1; i++)
     924     1027102 :     if (!(gel(y,i) = addRe_modIm(gel(x,i), s, pi2))) return NULL;
     925      571472 :   if (i <= RU)
     926             :   {
     927      229106 :     GEN pi4 = Pi2n(2, prec), s2 = gmul2n(s, 1);
     928      638635 :     for (   ; i <= RU; i++)
     929      409530 :       if (!(gel(y,i) = addRe_modIm(gel(x,i), s2, pi4))) return NULL;
     930             :   }
     931      571471 :   return y;
     932             : }
     933             : GEN
     934           0 : nf_cxlog_normalize(GEN nf, GEN x, long prec)
     935             : {
     936           0 :   long N = nf_get_degree(checknf(nf));
     937           0 :   return cleanarch(x, N, prec);
     938             : }
     939             : 
     940             : static GEN
     941         160 : not_given(long reason)
     942             : {
     943         160 :   if (DEBUGLEVEL)
     944           0 :     switch(reason)
     945             :     {
     946           0 :       case fupb_LARGE:
     947           0 :         pari_warn(warner,"fundamental units too large, not given");
     948           0 :         break;
     949           0 :       case fupb_PRECI:
     950           0 :         pari_warn(warner,"insufficient precision for fundamental units, not given");
     951           0 :         break;
     952             :     }
     953         160 :   return NULL;
     954             : }
     955             : 
     956             : /* check whether exp(x) will 1) get too big (real(x) large), 2) require
     957             :  * large accuracy for argument reduction (imag(x) large) */
     958             : static long
     959      597180 : expbitprec(GEN x, long *e)
     960             : {
     961             :   GEN re, im;
     962      597180 :   if (typ(x) != t_COMPLEX) re = x;
     963             :   else
     964             :   {
     965      412227 :     im = gel(x,2); *e = maxss(*e, expo(im) + 5 - bit_prec(im));
     966      412227 :     re = gel(x,1);
     967             :   }
     968      597180 :   return (expo(re) <= 20);
     969             : 
     970             : }
     971             : static long
     972      242384 : RgC_expbitprec(GEN x)
     973             : {
     974      242384 :   long l = lg(x), i, e = - (long)HIGHEXPOBIT;
     975      813342 :   for (i = 1; i < l; i++)
     976      570993 :     if (!expbitprec(gel(x,i), &e)) return LONG_MAX;
     977      242349 :   return e;
     978             : }
     979             : static long
     980        4277 : RgM_expbitprec(GEN x)
     981             : {
     982        4277 :   long i, j, I, J, e = - (long)HIGHEXPOBIT;
     983        4277 :   RgM_dimensions(x, &I,&J);
     984       10822 :   for (j = 1; j <= J; j++)
     985       32732 :     for (i = 1; i <= I; i++)
     986       26187 :       if (!expbitprec(gcoeff(x,i,j), &e)) return LONG_MAX;
     987        4256 :   return e;
     988             : }
     989             : 
     990             : static GEN
     991        1275 : FlxqX_chinese_unit(GEN X, GEN U, GEN invzk, GEN D, GEN T, ulong p)
     992             : {
     993        1275 :   long i, lU = lg(U), lX = lg(X), d = lg(invzk)-1;
     994        1275 :   GEN M = cgetg(lU, t_MAT);
     995        1275 :   if (D)
     996             :   {
     997        1131 :     D = Flv_inv(D, p);
     998       55581 :     for (i = 1; i < lX; i++)
     999       54450 :       if (uel(D, i) != 1)
    1000       45155 :         gel(X,i) = Flx_Fl_mul(gel(X,i), uel(D,i), p);
    1001             :   }
    1002        3589 :   for (i = 1; i < lU; i++)
    1003             :   {
    1004        2314 :     GEN H = FlxqV_factorback(X, gel(U, i), T, p);
    1005        2314 :     gel(M, i) = Flm_Flc_mul(invzk, Flx_to_Flv(H, d), p);
    1006             :   }
    1007        1275 :   return M;
    1008             : }
    1009             : 
    1010             : static GEN
    1011         263 : chinese_unit_slice(GEN A, GEN U, GEN B, GEN D, GEN C, GEN P, GEN *mod)
    1012             : {
    1013         263 :   pari_sp av = avma;
    1014         263 :   long i, n = lg(P)-1, v = varn(C);
    1015             :   GEN H, T;
    1016         263 :   if (n == 1)
    1017             :   {
    1018           0 :     ulong p = uel(P,1);
    1019           0 :     GEN a = ZXV_to_FlxV(A, p), b = ZM_to_Flm(B, p), c = ZX_to_Flx(C, p);
    1020           0 :     GEN d = D ? ZV_to_Flv(D, p): NULL;
    1021           0 :     GEN Hp = FlxqX_chinese_unit(a, U, b, d, c, p);
    1022           0 :     H = gerepileupto(av, Flm_to_ZM(Hp));
    1023           0 :     *mod = utoi(p);
    1024           0 :     return H;
    1025             :   }
    1026         263 :   T = ZV_producttree(P);
    1027         263 :   A = ZXC_nv_mod_tree(A, P, T, v);
    1028         263 :   B = ZM_nv_mod_tree(B, P, T);
    1029         263 :   D = D ? ZV_nv_mod_tree(D, P, T): NULL;
    1030         263 :   C = ZX_nv_mod_tree(C, P, T);
    1031             : 
    1032         263 :   H = cgetg(n+1, t_VEC);
    1033        1538 :   for(i=1; i <= n; i++)
    1034             :   {
    1035        1275 :     ulong p = P[i];
    1036        1275 :     GEN a = gel(A,i), b = gel(B,i), c = gel(C,i), d = D ? gel(D,i): NULL;
    1037        1275 :     gel(H,i) = FlxqX_chinese_unit(a, U, b, d, c, p);
    1038             :   }
    1039         263 :   H = nmV_chinese_center_tree_seq(H, P, T, ZV_chinesetree(P, T));
    1040         263 :   *mod = gmael(T, lg(T)-1, 1);
    1041         263 :   gerepileall(av, 2, &H, mod);
    1042         263 :   return H;
    1043             : }
    1044             : 
    1045             : GEN
    1046         263 : chinese_unit_worker(GEN P, GEN A, GEN U, GEN B, GEN D, GEN C)
    1047             : {
    1048         263 :   GEN V = cgetg(3, t_VEC);
    1049         263 :   gel(V,1) = chinese_unit_slice(A, U, B, isintzero(D) ? NULL: D, C, P, &gel(V,2));
    1050         263 :   return V;
    1051             : }
    1052             : 
    1053             : /* Let x = \prod X[i]^E[i] = u, return u.
    1054             :  * If dX != NULL, X[i] = nX[i] / dX[i] where nX[i] is a ZX, dX[i] in Z */
    1055             : static GEN
    1056         106 : chinese_unit(GEN nf, GEN nX, GEN dX, GEN U, ulong bnd)
    1057             : {
    1058         106 :   pari_sp av = avma;
    1059         106 :   GEN f = nf_get_index(nf), T = nf_get_pol(nf), invzk = nf_get_invzk(nf);
    1060             :   GEN H, mod;
    1061             :   forprime_t S;
    1062         106 :   GEN worker = snm_closure(is_entry("_chinese_unit_worker"),
    1063             :                mkcol5(nX, U, invzk, dX? dX: gen_0, T));
    1064         106 :   init_modular_big(&S);
    1065         106 :   H = gen_crt("chinese_units", worker, &S, f, bnd, 0, &mod, nmV_chinese_center, FpM_center);
    1066         106 :   settyp(H, t_VEC); return gerepilecopy(av, H);
    1067             : }
    1068             : 
    1069             : /* *pE a ZM */
    1070             : static void
    1071         148 : ZM_remove_unused(GEN *pE, GEN *pX)
    1072             : {
    1073         148 :   long j, k, l = lg(*pX);
    1074         148 :   GEN E = *pE, v = cgetg(l, t_VECSMALL);
    1075       14318 :   for (j = k = 1; j < l; j++)
    1076       14170 :     if (!ZMrow_equal0(E, j)) v[k++] = j;
    1077         148 :   if (k < l)
    1078             :   {
    1079         148 :     setlg(v, k);
    1080         148 :     *pX = vecpermute(*pX,v);
    1081         148 :     *pE = rowpermute(E,v);
    1082             :   }
    1083         148 : }
    1084             : 
    1085             : /* s = -log|norm(x)|/N */
    1086             : static GEN
    1087      248950 : fixarch(GEN x, GEN s, long R1)
    1088             : {
    1089             :   long i, l;
    1090      248950 :   GEN y = cgetg_copy(x, &l);
    1091      719269 :   for (i = 1; i <= R1; i++) gel(y,i) = gadd(s, gel(x,i));
    1092      375867 :   for (     ; i <   l; i++) gel(y,i) = gadd(s, gmul2n(gel(x,i),-1));
    1093      248950 :   return y;
    1094             : }
    1095             : 
    1096             : static GEN
    1097       12041 : getfu(GEN nf, GEN *ptA, GEN *ptU, long prec)
    1098             : {
    1099       12041 :   GEN U, y, matep, A, T = nf_get_pol(nf), M = nf_get_M(nf);
    1100       12041 :   long e, j, R1, RU, N = degpol(T);
    1101             : 
    1102       12041 :   R1 = nf_get_r1(nf); RU = (N+R1) >> 1;
    1103       12041 :   if (RU == 1) return cgetg(1,t_VEC);
    1104             : 
    1105        4277 :   A = *ptA;
    1106        4277 :   matep = cgetg(RU,t_MAT);
    1107       10843 :   for (j = 1; j < RU; j++)
    1108             :   {
    1109        6566 :     GEN Aj = gel(A,j), s = gdivgs(RgV_sum(real_i(Aj)), -N);
    1110        6566 :     gel(matep,j) = fixarch(Aj, s, R1);
    1111             :   }
    1112        4277 :   U = lll(real_i(matep));
    1113        4277 :   if (lg(U) < RU) return not_given(fupb_PRECI);
    1114        4277 :   if (ptU) { *ptU = U; *ptA = A = RgM_ZM_mul(A,U); }
    1115        4277 :   y = RgM_ZM_mul(matep,U);
    1116        4277 :   e = RgM_expbitprec(y);
    1117        4277 :   if (e >= 0) return not_given(e == LONG_MAX? fupb_LARGE: fupb_PRECI);
    1118        4256 :   if (prec <= 0) prec = gprecision(A);
    1119        4256 :   y = RgM_solve_realimag(M, gexp(y,prec));
    1120        4256 :   if (!y) return not_given(fupb_PRECI);
    1121        4256 :   y = grndtoi(y, &e); if (e >= 0) return not_given(fupb_PRECI);
    1122        4136 :   settyp(y, t_VEC);
    1123             : 
    1124        4136 :   if (!ptU) *ptA = A = RgM_ZM_mul(A, U);
    1125       10395 :   for (j = 1; j < RU; j++)
    1126             :   { /* y[i] are hopefully unit generators. Normalize: smallest T2 norm */
    1127        6278 :     GEN u = gel(y,j), v = zk_inv(nf, u);
    1128        6278 :     if (!v || !is_pm1(Q_denom(v)) || ZV_isscalar(u))
    1129          19 :       return not_given(fupb_PRECI);
    1130        6259 :     if (gcmp(RgC_fpnorml2(v,DEFAULTPREC), RgC_fpnorml2(u,DEFAULTPREC)) < 0)
    1131             :     {
    1132        2184 :       gel(A,j) = RgC_neg(gel(A,j));
    1133        2184 :       if (ptU) gel(U,j) = ZC_neg(gel(U,j));
    1134        2184 :       u = v;
    1135             :     }
    1136        6259 :     gel(y,j) = nf_to_scalar_or_alg(nf, u);
    1137             :   }
    1138        4117 :   return y;
    1139             : }
    1140             : 
    1141             : static void
    1142           0 : err_units() { pari_err_PREC("makeunits [cannot get units, use bnfinit(,1)]"); }
    1143             : 
    1144             : /* bound for log2 |sigma(u)|, sigma complex embedding, u fundamental unit
    1145             :  * attached to bnf_get_logfu */
    1146             : static double
    1147         106 : log2fubound(GEN bnf)
    1148             : {
    1149         106 :   GEN LU = bnf_get_logfu(bnf);
    1150         106 :   long i, j, l = lg(LU), r1 = nf_get_r1(bnf_get_nf(bnf));
    1151         106 :   double e = 0.0;
    1152         360 :   for (j = 1; j < l; j++)
    1153             :   {
    1154         254 :     GEN u = gel(LU,j);
    1155         669 :     for (i = 1; i <= r1; i++)
    1156             :     {
    1157         415 :       GEN E = real_i(gel(u,i));
    1158         415 :       e = maxdd(e, gtodouble(E));
    1159             :     }
    1160         895 :     for (     ; i <= l; i++)
    1161             :     {
    1162         641 :       GEN E = real_i(gel(u,i));
    1163         641 :       e = maxdd(e, gtodouble(E) / 2);
    1164             :     }
    1165             :   }
    1166         106 :   return e / M_LN2;
    1167             : }
    1168             : /* bound for log2(|split_real_imag(M, y)|_oo / |y|_oo)*/
    1169             : static double
    1170         106 : log2Mbound(GEN nf)
    1171             : {
    1172         106 :   GEN G = nf_get_G(nf), D = nf_get_disc(nf);
    1173         106 :   long r2 = nf_get_r2(nf), l = lg(G), i;
    1174         106 :   double e, d = dbllog2(D)/2 - r2 * M_LN2; /* log2 |det(split_real_imag(M))| */
    1175         106 :   e = log2(nf_get_degree(nf));
    1176         565 :   for (i = 2; i < l; i++) e += dbllog2(gnorml2(gel(G,i))); /* Hadamard bound */
    1177         106 :   return e / 2 - d;
    1178             : }
    1179             : 
    1180             : static GEN
    1181         106 : vec_chinese_unit(GEN bnf)
    1182             : {
    1183         106 :   GEN nf = bnf_get_nf(bnf), SUnits = bnf_get_sunits(bnf);
    1184         106 :   ulong bnd = (ulong)ceil(log2Mbound(nf) + log2fubound(bnf));
    1185         106 :   GEN X, dX, Y, U, f = nf_get_index(nf);
    1186         106 :   long j, l, v = nf_get_varn(nf);
    1187         106 :   if (!SUnits) err_units(); /* no compact units */
    1188         106 :   Y = gel(SUnits,1);
    1189         106 :   U = gel(SUnits,2);
    1190         106 :   ZM_remove_unused(&U, &Y); l = lg(Y); X = cgetg(l, t_VEC);
    1191         106 :   if (is_pm1(f)) f = dX = NULL; else dX = cgetg(l, t_VEC);
    1192        5522 :   for (j = 1; j < l; j++)
    1193             :   {
    1194        5416 :     GEN t = nf_to_scalar_or_alg(nf, gel(Y,j));
    1195        5416 :     if (f)
    1196             :     {
    1197             :       GEN den;
    1198        4366 :       t = Q_remove_denom(t, &den);
    1199        4366 :       gel(dX,j) = den ? den: gen_1;
    1200             :     }
    1201        5416 :     gel(X,j) = typ(t) == t_INT? scalarpol_shallow(t,v): t;
    1202             :   }
    1203         106 :   return chinese_unit(nf, X, dX, U, bnd);
    1204             : }
    1205             : 
    1206             : static GEN
    1207         855 : makeunits(GEN bnf)
    1208             : {
    1209         855 :   GEN nf = bnf_get_nf(bnf), fu = bnf_get_fu_nocheck(bnf);
    1210         855 :   GEN tu = nf_to_scalar_or_basis(nf, bnf_get_tuU(bnf));
    1211         855 :   fu = (typ(fu) == t_MAT)? vec_chinese_unit(bnf): matalgtobasis(nf, fu);
    1212         855 :   return vec_prepend(fu, tu);
    1213             : }
    1214             : 
    1215             : /*******************************************************************/
    1216             : /*                                                                 */
    1217             : /*           PRINCIPAL IDEAL ALGORITHM (DISCRETE LOG)              */
    1218             : /*                                                                 */
    1219             : /*******************************************************************/
    1220             : 
    1221             : /* G: prime ideals, E: vector of non-negative exponents.
    1222             :  * C = possible extra prime (^1) or NULL
    1223             :  * Return Norm (product) */
    1224             : static GEN
    1225          21 : get_norm_fact_primes(GEN G, GEN E, GEN C)
    1226             : {
    1227          21 :   pari_sp av=avma;
    1228          21 :   GEN N = gen_1, P, p;
    1229          21 :   long i, c = lg(E);
    1230          21 :   for (i=1; i<c; i++)
    1231             :   {
    1232           0 :     GEN ex = gel(E,i);
    1233           0 :     long s = signe(ex);
    1234           0 :     if (!s) continue;
    1235             : 
    1236           0 :     P = gel(G,i); p = pr_get_p(P);
    1237           0 :     N = mulii(N, powii(p, mului(pr_get_f(P), ex)));
    1238             :   }
    1239          21 :   if (C) N = mulii(N, pr_norm(C));
    1240          21 :   return gerepileuptoint(av, N);
    1241             : }
    1242             : 
    1243             : /* gen: HNF ideals */
    1244             : static GEN
    1245      240228 : get_norm_fact(GEN gen, GEN ex, GEN *pd)
    1246             : {
    1247      240228 :   long i, c = lg(ex);
    1248             :   GEN d,N,I,e,n,ne,de;
    1249      240228 :   d = N = gen_1;
    1250      393748 :   for (i=1; i<c; i++)
    1251      153520 :     if (signe(gel(ex,i)))
    1252             :     {
    1253      103580 :       I = gel(gen,i); e = gel(ex,i); n = ZM_det_triangular(I);
    1254      103580 :       ne = powii(n,e);
    1255      103580 :       de = equalii(n, gcoeff(I,1,1))? ne: powii(gcoeff(I,1,1), e);
    1256      103580 :       N = mulii(N, ne);
    1257      103580 :       d = mulii(d, de);
    1258             :     }
    1259      240228 :   *pd = d; return N;
    1260             : }
    1261             : 
    1262             : static GEN
    1263      392846 : get_pr_lists(GEN FB, long N, int list_pr)
    1264             : {
    1265             :   GEN pr, L;
    1266      392846 :   long i, l = lg(FB), p, pmax;
    1267             : 
    1268      392846 :   pmax = 0;
    1269     3488533 :   for (i=1; i<l; i++)
    1270             :   {
    1271     3095687 :     pr = gel(FB,i); p = pr_get_smallp(pr);
    1272     3095687 :     if (p > pmax) pmax = p;
    1273             :   }
    1274      392846 :   L = const_vec(pmax, NULL);
    1275      392846 :   if (list_pr)
    1276             :   {
    1277           0 :     for (i=1; i<l; i++)
    1278             :     {
    1279           0 :       pr = gel(FB,i); p = pr_get_smallp(pr);
    1280           0 :       if (!L[p]) gel(L,p) = vectrunc_init(N+1);
    1281           0 :       vectrunc_append(gel(L,p), pr);
    1282             :     }
    1283           0 :     for (p=1; p<=pmax; p++)
    1284           0 :       if (L[p]) gen_sort_inplace(gel(L,p), (void*)&cmp_prime_over_p,
    1285             :                                  &cmp_nodata, NULL);
    1286             :   }
    1287             :   else
    1288             :   {
    1289     3488533 :     for (i=1; i<l; i++)
    1290             :     {
    1291     3095687 :       pr = gel(FB,i); p = pr_get_smallp(pr);
    1292     3095687 :       if (!L[p]) gel(L,p) = vecsmalltrunc_init(N+1);
    1293     3095687 :       vecsmalltrunc_append(gel(L,p), i);
    1294             :     }
    1295             :   }
    1296      392846 :   return L;
    1297             : }
    1298             : 
    1299             : /* recover FB, LV, iLP, KCZ from Vbase */
    1300             : static GEN
    1301      392846 : recover_partFB(FB_t *F, GEN Vbase, long N)
    1302             : {
    1303      392846 :   GEN FB, LV, iLP, L = get_pr_lists(Vbase, N, 0);
    1304      392846 :   long l = lg(L), p, ip, i;
    1305             : 
    1306      392846 :   i = ip = 0;
    1307      392846 :   FB = cgetg(l, t_VECSMALL);
    1308      392846 :   iLP= cgetg(l, t_VECSMALL);
    1309      392846 :   LV = cgetg(l, t_VEC);
    1310     7084845 :   for (p = 2; p < l; p++)
    1311             :   {
    1312     6691999 :     if (!L[p]) continue;
    1313     1684109 :     FB[++i] = p;
    1314     1684109 :     gel(LV,p) = vecpermute(Vbase, gel(L,p));
    1315     1684109 :     iLP[p]= ip; ip += lg(gel(L,p))-1;
    1316             :   }
    1317      392846 :   F->KCZ = i;
    1318      392846 :   F->KC = ip;
    1319      392846 :   F->FB = FB; setlg(FB, i+1);
    1320      392846 :   F->LV = (GEN*)LV;
    1321      392846 :   F->iLP= iLP; return L;
    1322             : }
    1323             : 
    1324             : /* add v^e to factorization */
    1325             : static void
    1326       15336 : add_to_fact(long v, long e, FACT *fact)
    1327             : {
    1328       15336 :   long i, l = fact[0].pr;
    1329       55826 :   for (i=1; i<=l && fact[i].pr < v; i++)/*empty*/;
    1330       15336 :   if (i <= l && fact[i].pr == v) fact[i].ex += e; else store(v, e, fact);
    1331       15336 : }
    1332             : static void
    1333           0 : inv_fact(FACT *fact)
    1334             : {
    1335           0 :   long i, l = fact[0].pr;
    1336           0 :   for (i=1; i<=l; i++) fact[i].ex = -fact[i].ex;
    1337           0 : }
    1338             : 
    1339             : /* L (small) list of primes above the same p including pr. Return pr index */
    1340             : static int
    1341        3024 : pr_index(GEN L, GEN pr)
    1342             : {
    1343        3024 :   long j, l = lg(L);
    1344        3024 :   GEN al = pr_get_gen(pr);
    1345        3024 :   for (j=1; j<l; j++)
    1346        3024 :     if (ZV_equal(al, pr_get_gen(gel(L,j)))) return j;
    1347           0 :   pari_err_BUG("codeprime");
    1348             :   return 0; /* LCOV_EXCL_LINE */
    1349             : }
    1350             : 
    1351             : static long
    1352        3024 : Vbase_to_FB(FB_t *F, GEN pr)
    1353             : {
    1354        3024 :   long p = pr_get_smallp(pr);
    1355        3024 :   return F->iLP[p] + pr_index(F->LV[p], pr);
    1356             : }
    1357             : 
    1358             : /* x, y 2 extended ideals whose first component is an integral HNF and second
    1359             :  * a famat */
    1360             : static GEN
    1361        1887 : idealHNF_mulred(GEN nf, GEN x, GEN y)
    1362             : {
    1363        1887 :   GEN A = idealHNF_mul(nf, gel(x,1), gel(y,1));
    1364        1887 :   GEN F = famat_mul_shallow(gel(x,2), gel(y,2));
    1365        1887 :   return idealred(nf, mkvec2(A, F));
    1366             : }
    1367             : /* idealred(x * pr^n), n > 0 is small, x extended ideal. Reduction in order to
    1368             :  * avoid prec pb: don't let id become too large as lgsub increases */
    1369             : static GEN
    1370        4477 : idealmulpowprime2(GEN nf, GEN x, GEN pr, ulong n)
    1371             : {
    1372        4477 :   GEN A = idealmulpowprime(nf, gel(x,1), pr, utoipos(n));
    1373        4477 :   return mkvec2(A, gel(x,2));
    1374             : }
    1375             : static GEN
    1376       13615 : init_famat(GEN x) { return mkvec2(x, trivial_fact()); }
    1377             : /* optimized idealfactorback + reduction; z = init_famat() */
    1378             : static GEN
    1379       11750 : powred(GEN z, GEN nf, GEN p, GEN e)
    1380       11750 : { gel(z,1) = p; return idealpowred(nf, z, e); }
    1381             : static GEN
    1382        9863 : genback(GEN z, GEN nf, GEN P, GEN E)
    1383             : {
    1384        9863 :   long i, l = lg(E);
    1385        9863 :   GEN I = NULL;
    1386       26524 :   for (i = 1; i < l; i++)
    1387       16661 :     if (signe(gel(E,i)))
    1388             :     {
    1389       11750 :       GEN J = powred(z, nf, gel(P,i), gel(E,i));
    1390       11750 :       I = I? idealHNF_mulred(nf, I, J): J;
    1391             :     }
    1392        9863 :   return I; /* != NULL since a generator */
    1393             : }
    1394             : 
    1395             : /* return famat y (principal ideal) such that y / x is smooth [wrt Vbase] */
    1396             : static GEN
    1397      409128 : SPLIT(FB_t *F, GEN nf, GEN x, GEN Vbase, FACT *fact)
    1398             : {
    1399      409128 :   GEN vecG, ex, Ly, y, x0, Nx = ZM_det_triangular(x);
    1400             :   long nbtest_lim, nbtest, i, j, k, ru, lgsub;
    1401             :   pari_sp av;
    1402             : 
    1403             :   /* try without reduction if x is small */
    1404      818144 :   if (gexpo(gcoeff(x,1,1)) < 100 &&
    1405      502891 :       can_factor(F, nf, x, NULL, Nx, fact)) return NULL;
    1406             : 
    1407      315253 :   av = avma;
    1408      315253 :   Ly = idealpseudominvec(x, nf_get_roundG(nf));
    1409      359258 :   for(k=1; k<lg(Ly); k++)
    1410             :   {
    1411      350648 :     y = gel(Ly,k);
    1412      350648 :     if (factorgen(F, nf, x, Nx, y, fact)) return y;
    1413             :   }
    1414        8610 :   set_avma(av);
    1415             : 
    1416             :   /* reduce in various directions */
    1417        8610 :   ru = lg(nf_get_roots(nf));
    1418        8610 :   vecG = cgetg(ru, t_VEC);
    1419       13573 :   for (j=1; j<ru; j++)
    1420             :   {
    1421       11998 :     gel(vecG,j) = nf_get_Gtwist1(nf, j);
    1422       11998 :     av = avma;
    1423       11998 :     Ly = idealpseudominvec(x, gel(vecG,j));
    1424       33495 :     for(k=1; k<lg(Ly); k++)
    1425             :     {
    1426       28532 :       y = gel(Ly,k);
    1427       28532 :       if (factorgen(F, nf, x, Nx, y, fact)) return y;
    1428             :     }
    1429        4963 :     set_avma(av);
    1430             :   }
    1431             : 
    1432             :   /* tough case, multiply by random products */
    1433        1575 :   lgsub = 3;
    1434        1575 :   ex = cgetg(lgsub, t_VECSMALL);
    1435        1575 :   x0 = init_famat(x);
    1436        1575 :   nbtest = 1; nbtest_lim = 4;
    1437             :   for(;;)
    1438         637 :   {
    1439        2212 :     GEN Ired, I, NI, id = x0;
    1440        2212 :     av = avma;
    1441        2212 :     if (DEBUGLEVEL>2) err_printf("# ideals tried = %ld\n",nbtest);
    1442        6945 :     for (i=1; i<lgsub; i++)
    1443             :     {
    1444        4733 :       ex[i] = random_bits(RANDOM_BITS);
    1445        4733 :       if (ex[i]) id = idealmulpowprime2(nf, id, gel(Vbase,i), ex[i]);
    1446             :     }
    1447        2212 :     if (id == x0) continue;
    1448             :     /* I^(-1) * \prod Vbase[i]^ex[i] = (id[2]) / x */
    1449             : 
    1450        2205 :     I = gel(id,1); NI = ZM_det_triangular(I);
    1451        2205 :     if (can_factor(F, nf, I, NULL, NI, fact))
    1452             :     {
    1453           0 :       inv_fact(fact); /* I^(-1) */
    1454           0 :       for (i=1; i<lgsub; i++)
    1455           0 :         if (ex[i]) add_to_fact(Vbase_to_FB(F,gel(Vbase,i)), ex[i], fact);
    1456           0 :       return gel(id,2);
    1457             :     }
    1458        2205 :     Ired = ru == 2? I: ZM_lll(I, 0.99, LLL_INPLACE);
    1459        3669 :     for (j=1; j<ru; j++)
    1460             :     {
    1461        3039 :       pari_sp av2 = avma;
    1462        3039 :       Ly = idealpseudominvec(Ired, gel(vecG,j));
    1463        8302 :       for (k=1; k < lg(Ly); k++)
    1464             :       {
    1465        6838 :         y = gel(Ly,k);
    1466        6838 :         if (factorgen(F, nf, I, NI, y, fact))
    1467             :         {
    1468        4759 :           for (i=1; i<lgsub; i++)
    1469        3184 :             if (ex[i]) add_to_fact(Vbase_to_FB(F,gel(Vbase,i)), ex[i], fact);
    1470        1575 :           return famat_mul_shallow(gel(id,2), y);
    1471             :         }
    1472             :       }
    1473        1464 :       set_avma(av2);
    1474             :     }
    1475         630 :     set_avma(av);
    1476         630 :     if (++nbtest > nbtest_lim)
    1477             :     {
    1478          34 :       nbtest = 0;
    1479          34 :       if (++lgsub < minss(8, lg(Vbase)-1))
    1480             :       {
    1481          34 :         nbtest_lim <<= 1;
    1482          34 :         ex = cgetg(lgsub, t_VECSMALL);
    1483             :       }
    1484           0 :       else nbtest_lim = LONG_MAX; /* don't increase further */
    1485          34 :       if (DEBUGLEVEL>2) err_printf("SPLIT: increasing factor base [%ld]\n",lgsub);
    1486             :     }
    1487             :   }
    1488             : }
    1489             : 
    1490             : INLINE GEN
    1491      392804 : bnf_get_W(GEN bnf) { return gel(bnf,1); }
    1492             : INLINE GEN
    1493      785594 : bnf_get_B(GEN bnf) { return gel(bnf,2); }
    1494             : INLINE GEN
    1495      801015 : bnf_get_C(GEN bnf) { return gel(bnf,4); }
    1496             : INLINE GEN
    1497      392860 : bnf_get_vbase(GEN bnf) { return gel(bnf,5); }
    1498             : INLINE GEN
    1499      392790 : bnf_get_Ur(GEN bnf) { return gmael(bnf,9,1); }
    1500             : INLINE GEN
    1501      145307 : bnf_get_ga(GEN bnf) { return gmael(bnf,9,2); }
    1502             : INLINE GEN
    1503      147778 : bnf_get_GD(GEN bnf) { return gmael(bnf,9,3); }
    1504             : 
    1505             : /* Return y (as an elt of K or a t_MAT representing an elt in Z[K])
    1506             :  * such that x / (y) is smooth and store the exponents of  its factorization
    1507             :  * on g_W and g_B in Wex / Bex; return NULL for y = 1 */
    1508             : static GEN
    1509      392790 : split_ideal(GEN bnf, GEN x, GEN *pWex, GEN *pBex)
    1510             : {
    1511      392790 :   GEN L, y, Vbase = bnf_get_vbase(bnf);
    1512      392790 :   GEN Wex, W  = bnf_get_W(bnf);
    1513      392790 :   GEN Bex, B  = bnf_get_B(bnf);
    1514             :   long p, j, i, l, nW, nB;
    1515             :   FACT *fact;
    1516             :   FB_t F;
    1517             : 
    1518      392790 :   L = recover_partFB(&F, Vbase, lg(x)-1);
    1519      392790 :   fact = (FACT*)stack_malloc((F.KC+1)*sizeof(FACT));
    1520      392790 :   y = SPLIT(&F, bnf_get_nf(bnf), x, Vbase, fact);
    1521      392790 :   nW = lg(W)-1; *pWex = Wex = zero_zv(nW);
    1522      392790 :   nB = lg(B)-1; *pBex = Bex = zero_zv(nB); l = lg(F.FB);
    1523      392790 :   p = j = 0; /* -Wall */
    1524      765409 :   for (i = 1; i <= fact[0].pr; i++)
    1525             :   { /* decode index C = ip+j --> (p,j) */
    1526      372619 :     long a, b, t, C = fact[i].pr;
    1527     1160823 :     for (t = 1; t < l; t++)
    1528             :     {
    1529     1113679 :       long q = F.FB[t], k = C - F.iLP[q];
    1530     1113679 :       if (k <= 0) break;
    1531      788204 :       p = q;
    1532      788204 :       j = k;
    1533             :     }
    1534      372619 :     a = gel(L, p)[j];
    1535      372619 :     b = a - nW;
    1536      372619 :     if (b <= 0) Wex[a] = y? -fact[i].ex: fact[i].ex;
    1537      278302 :     else        Bex[b] = y? -fact[i].ex: fact[i].ex;
    1538             :   }
    1539      392790 :   return y;
    1540             : }
    1541             : 
    1542             : GEN
    1543      208365 : init_red_mod_units(GEN bnf, long prec)
    1544             : {
    1545      208365 :   GEN s = gen_0, p1,s1,mat, logfu = bnf_get_logfu(bnf);
    1546      208365 :   long i,j, RU = lg(logfu);
    1547             : 
    1548      208365 :   if (RU == 1) return NULL;
    1549      208365 :   mat = cgetg(RU,t_MAT);
    1550      537009 :   for (j=1; j<RU; j++)
    1551             :   {
    1552      328644 :     p1 = cgetg(RU+1,t_COL); gel(mat,j) = p1;
    1553      328644 :     s1 = gen_0;
    1554      942716 :     for (i=1; i<RU; i++)
    1555             :     {
    1556      614072 :       gel(p1,i) = real_i(gcoeff(logfu,i,j));
    1557      614072 :       s1 = mpadd(s1, mpsqr(gel(p1,i)));
    1558             :     }
    1559      328644 :     gel(p1,RU) = gen_0; if (mpcmp(s1,s) > 0) s = s1;
    1560             :   }
    1561      208365 :   s = gsqrt(gmul2n(s,RU),prec);
    1562      208365 :   if (expo(s) < 27) s = utoipos(1UL << 27);
    1563      208365 :   return mkvec2(mat, s);
    1564             : }
    1565             : 
    1566             : /* z computed above. Return unit exponents that would reduce col (arch) */
    1567             : GEN
    1568      208365 : red_mod_units(GEN col, GEN z)
    1569             : {
    1570             :   long i,RU;
    1571             :   GEN x,mat,N2;
    1572             : 
    1573      208365 :   if (!z) return NULL;
    1574      208365 :   mat= gel(z,1);
    1575      208365 :   N2 = gel(z,2);
    1576      208365 :   RU = lg(mat); x = cgetg(RU+1,t_COL);
    1577      537009 :   for (i=1; i<RU; i++) gel(x,i) = real_i(gel(col,i));
    1578      208365 :   gel(x,RU) = N2;
    1579      208365 :   x = lll(shallowconcat(mat,x));
    1580      208365 :   if (typ(x) != t_MAT || lg(x) <= RU) return NULL;
    1581      208365 :   x = gel(x,RU);
    1582      208365 :   if (signe(gel(x,RU)) < 0) x = gneg_i(x);
    1583      208365 :   if (!gequal1(gel(x,RU))) pari_err_BUG("red_mod_units");
    1584      208365 :   setlg(x,RU); return x;
    1585             : }
    1586             : 
    1587             : static GEN
    1588      717212 : add(GEN a, GEN t) { return a = a? RgC_add(a,t): t; }
    1589             : 
    1590             : /* [x] archimedian components, A column vector. return [x] A */
    1591             : static GEN
    1592      598842 : act_arch(GEN A, GEN x)
    1593             : {
    1594             :   GEN a;
    1595      598842 :   long i,l = lg(A), tA = typ(A);
    1596      598842 :   if (tA == t_MAT)
    1597             :   { /* assume lg(x) >= l */
    1598       36393 :     a = cgetg(l, t_MAT);
    1599       68000 :     for (i=1; i<l; i++) gel(a,i) = act_arch(gel(A,i), x);
    1600       36393 :     return a;
    1601             :   }
    1602      562449 :   if (l==1) return cgetg(1, t_COL);
    1603      562449 :   a = NULL;
    1604      562449 :   if (tA == t_VECSMALL)
    1605             :   {
    1606     1771469 :     for (i=1; i<l; i++)
    1607             :     {
    1608     1531241 :       long c = A[i];
    1609     1531241 :       if (c) a = add(a, gmulsg(c, gel(x,i)));
    1610             :     }
    1611             :   }
    1612             :   else
    1613             :   { /* A a t_COL of t_INT. Assume lg(A)==lg(x) */
    1614      687861 :     for (i=1; i<l; i++)
    1615             :     {
    1616      365640 :       GEN c = gel(A,i);
    1617      365640 :       if (signe(c)) a = add(a, gmul(c, gel(x,i)));
    1618             :     }
    1619             :   }
    1620      562449 :   return a? a: zerocol(lgcols(x)-1);
    1621             : }
    1622             : /* act_arch(matdiagonal(v), x) */
    1623             : static GEN
    1624       12131 : diagact_arch(GEN v, GEN x)
    1625             : {
    1626       12131 :   long i, l = lg(v);
    1627       12131 :   GEN a = cgetg(l, t_MAT);
    1628       22064 :   for (i = 1; i < l; i++) gel(a,i) = gmul(gel(x,i), gel(v,i));
    1629       12131 :   return a;
    1630             : }
    1631             : 
    1632             : static long
    1633      400756 : prec_arch(GEN bnf)
    1634             : {
    1635      400756 :   GEN a = bnf_get_C(bnf);
    1636      400756 :   long i, l = lg(a), prec;
    1637             : 
    1638      400756 :   for (i=1; i<l; i++)
    1639      400672 :     if ( (prec = gprecision(gel(a,i))) ) return prec;
    1640          84 :   return DEFAULTPREC;
    1641             : }
    1642             : 
    1643             : static long
    1644         869 : needed_bitprec(GEN x)
    1645             : {
    1646         869 :   long i, e = 0, l = lg(x);
    1647        7546 :   for (i = 1; i < l; i++)
    1648             :   {
    1649        6677 :     GEN c = gel(x,i);
    1650        6677 :     long f = gexpo(c) - prec2nbits(gprecision(c));
    1651        6677 :     if (f > e) e = f;
    1652             :   }
    1653         869 :   return e;
    1654             : }
    1655             : 
    1656             : /* col = archimedian components of x, Nx its norm, dx a multiple of its
    1657             :  * denominator. Return x or NULL (fail) */
    1658             : GEN
    1659      242384 : isprincipalarch(GEN bnf, GEN col, GEN kNx, GEN e, GEN dx, long *pe)
    1660             : {
    1661             :   GEN nf, x, y, logfu, s, M;
    1662      242384 :   long N, prec = gprecision(col);
    1663      242384 :   bnf = checkbnf(bnf); nf = bnf_get_nf(bnf); M = nf_get_M(nf);
    1664      242384 :   if (!prec) prec = prec_arch(bnf);
    1665      242384 :   *pe = 128;
    1666      242384 :   logfu = bnf_get_logfu(bnf);
    1667      242384 :   N = nf_get_degree(nf);
    1668      242384 :   if (!(col = cleanarch(col,N,prec))) return NULL;
    1669      242384 :   if (lg(col) > 2)
    1670             :   { /* reduce mod units */
    1671      208365 :     GEN u, z = init_red_mod_units(bnf,prec);
    1672      208365 :     if (!(u = red_mod_units(col,z))) return NULL;
    1673      208365 :     col = RgC_add(col, RgM_RgC_mul(logfu, u));
    1674      208365 :     if (!(col = cleanarch(col,N,prec))) return NULL;
    1675             :   }
    1676      242384 :   s = divru(mulir(e, glog(kNx,prec)), N);
    1677      242384 :   col = fixarch(col, s, nf_get_r1(nf));
    1678      242384 :   if (RgC_expbitprec(col) >= 0) return NULL;
    1679      242349 :   col = gexp(col, prec);
    1680             :   /* d.alpha such that x = alpha \prod gj^ej */
    1681      242349 :   x = RgM_solve_realimag(M,col); if (!x) return NULL;
    1682      242349 :   x = RgC_Rg_mul(x, dx);
    1683      242349 :   y = grndtoi(x, pe);
    1684      242349 :   if (*pe > -5) { *pe = needed_bitprec(x); return NULL; }
    1685      241480 :   return RgC_Rg_div(y, dx);
    1686             : }
    1687             : 
    1688             : /* y = C \prod g[i]^e[i] ? */
    1689             : static int
    1690      241480 : fact_ok(GEN nf, GEN y, GEN C, GEN g, GEN e)
    1691             : {
    1692      241480 :   pari_sp av = avma;
    1693      241480 :   long i, c = lg(e);
    1694      241480 :   GEN z = C? C: gen_1;
    1695      394945 :   for (i=1; i<c; i++)
    1696      153465 :     if (signe(gel(e,i))) z = idealmul(nf, z, idealpow(nf, gel(g,i), gel(e,i)));
    1697      241480 :   if (typ(z) != t_MAT) z = idealhnf_shallow(nf,z);
    1698      241480 :   if (typ(y) != t_MAT) y = idealhnf_shallow(nf,y);
    1699      241480 :   return gc_bool(av, ZM_equal(y,z));
    1700             : }
    1701             : static GEN
    1702      392790 : ZV_divrem(GEN A, GEN B, GEN *pR)
    1703             : {
    1704      392790 :   long i, l = lg(A);
    1705      392790 :   GEN Q = cgetg(l, t_COL), R = cgetg(l, t_COL);
    1706      734523 :   for (i = 1; i < l; i++) gel(Q,i) = truedvmdii(gel(A,i), gel(B,i), &gel(R,i));
    1707      392790 :   *pR = R; return Q;
    1708             : }
    1709             : 
    1710             : static GEN
    1711      392790 : Ur_ZC_mul(GEN bnf, GEN v)
    1712             : {
    1713      392790 :   GEN w, U = bnf_get_Ur(bnf);
    1714      392790 :   long i, l = lg(bnf_get_cyc(bnf)); /* may be < lgcols(U) */
    1715             : 
    1716      392790 :   w = cgetg(l, t_COL);
    1717      734523 :   for (i = 1; i < l; i++) gel(w,i) = ZMrow_ZC_mul(U, v, i);
    1718      392790 :   return w;
    1719             : }
    1720             : 
    1721             : static GEN
    1722         697 : ZV_mul(GEN x, GEN y)
    1723             : {
    1724         697 :   long i, l = lg(x);
    1725         697 :   GEN z = cgetg(l, t_COL);
    1726        2763 :   for (i = 1; i < l; i++) gel(z,i) = mulii(gel(x,i), gel(y,i));
    1727         697 :   return z;
    1728             : }
    1729             : 
    1730             : /* assume x in HNF; cf class_group_gen for notations. Return NULL iff
    1731             :  * flag & nf_FORCE and computation of principal ideal generator fails */
    1732             : static GEN
    1733      400245 : isprincipalall(GEN bnf, GEN x, long *pprec, long flag)
    1734             : {
    1735             :   GEN xar, Wex, Bex, gen, xc, col, A, Q, R, UA;
    1736      400245 :   GEN C = bnf_get_C(bnf), nf = bnf_get_nf(bnf), cyc = bnf_get_cyc(bnf);
    1737             :   long nB, nW, e;
    1738             : 
    1739      400245 :   if (lg(cyc) == 1 && !(flag & (nf_GEN|nf_GENMAT|nf_GEN_IF_PRINCIPAL)))
    1740         931 :     return cgetg(1,t_COL);
    1741      399314 :   if (lg(x) == 2)
    1742             :   { /* nf = Q */
    1743          84 :     col = gel(x,1);
    1744          84 :     if (flag & nf_GENMAT) col = to_famat_shallow(col, gen_1);
    1745          84 :     return (flag & nf_GEN_IF_PRINCIPAL)? col: mkvec2(cgetg(1,t_COL), col);
    1746             :   }
    1747             : 
    1748      399230 :   x = Q_primitive_part(x, &xc);
    1749      399230 :   if (equali1(gcoeff(x,1,1))) /* trivial ideal */
    1750             :   {
    1751        6440 :     R = zerocol(lg(cyc)-1);
    1752        6440 :     if (!(flag & (nf_GEN|nf_GENMAT|nf_GEN_IF_PRINCIPAL))) return R;
    1753        6419 :     if (flag & nf_GEN_IF_PRINCIPAL)
    1754        4641 :       return scalarcol_shallow(xc? xc: gen_1, nf_get_degree(nf));
    1755        1778 :     if (flag & nf_GENMAT)
    1756        1141 :       col = xc? to_famat_shallow(xc, gen_1): trivial_fact();
    1757             :     else
    1758         637 :       col = scalarcol_shallow(xc? xc: gen_1, nf_get_degree(nf));
    1759        1778 :     return mkvec2(R, col);
    1760             :   }
    1761      392790 :   xar = split_ideal(bnf, x, &Wex, &Bex);
    1762             :   /* x = g_W Wex + g_B Bex + [xar] = g_W (Wex - B*Bex) + [xar] + [C_B]Bex */
    1763      392790 :   A = zc_to_ZC(Wex); nB = lg(Bex)-1;
    1764      392790 :   if (nB) A = ZC_sub(A, ZM_zc_mul(bnf_get_B(bnf), Bex));
    1765      392790 :   UA = Ur_ZC_mul(bnf, A);
    1766      392790 :   Q = ZV_divrem(UA, cyc, &R);
    1767             :   /* g_W (Wex - B*Bex) = G Ur A - [ga]A = G R + [GD]Q - [ga]A
    1768             :    * Finally: x = G R + [xar] + [C_B]Bex + [GD]Q - [ga]A */
    1769      392790 :   if (!(flag & (nf_GEN|nf_GENMAT|nf_GEN_IF_PRINCIPAL))) return R;
    1770      240742 :   if ((flag & nf_GEN_IF_PRINCIPAL) && !ZV_equal0(R)) return gen_0;
    1771             : 
    1772      240735 :   nW = lg(Wex)-1;
    1773      240735 :   gen = bnf_get_gen(bnf);
    1774      240735 :   col = NULL;
    1775      240735 :   if (lg(R) == 1
    1776      145814 :       || abscmpiu(gel(R,vecindexmax(R)), 4 * bit_accuracy(*pprec)) < 0)
    1777             :   { /* q = N (x / prod gj^ej) = N(alpha), denom(alpha) | d */
    1778      240228 :     GEN d, q = gdiv(ZM_det_triangular(x), get_norm_fact(gen, R, &d));
    1779      240228 :     col = xar? nf_cxlog(nf, xar, *pprec): NULL;
    1780      240228 :     if (nB) col = add(col, act_arch(Bex, nW? vecslice(C,nW+1,lg(C)-1): C));
    1781      240228 :     if (nW) col = add(col, RgC_sub(act_arch(Q, bnf_get_GD(bnf)),
    1782             :                                    act_arch(A, bnf_get_ga(bnf))));
    1783      240228 :     col = isprincipalarch(bnf, col, q, gen_1, d, &e);
    1784      240228 :     if (col && !fact_ok(nf,x, col,gen,R)) col = NULL;
    1785             :   }
    1786      240735 :   if (!col && (flag & nf_GENMAT))
    1787             :   {
    1788        1272 :     GEN SUnits = bnf_get_sunits(bnf);
    1789        1272 :     if (SUnits)
    1790             :     {
    1791         740 :       GEN X = gel(SUnits,1), U = gel(SUnits,2), C = gel(SUnits,3);
    1792         740 :       GEN v = gel(bnf,9), Ge = gel(v,4), M1 = gel(v,5), M2 = gel(v,6);
    1793         740 :       GEN z = NULL, F = NULL;
    1794         740 :       if (nB)
    1795             :       {
    1796         740 :         GEN C2 = nW? vecslice(C, nW+1, lg(C)-1): C;
    1797         740 :         z = ZM_zc_mul(C2, Bex);
    1798             :       }
    1799         740 :       if (nW)
    1800             :       { /* [GD]Q - [ga]A = ([X]M1 - [Ge]D) Q - ([X]M2 - [Ge]Ur) A */
    1801         697 :         GEN C1 = vecslice(C, 1, nW);
    1802         697 :         GEN v = ZC_sub(ZM_ZC_mul(M1,Q), ZM_ZC_mul(M2,A));
    1803         697 :         z = add(z, ZM_ZC_mul(C1, v));
    1804         697 :         F = famat_reduce(famatV_factorback(Ge, ZC_sub(UA, ZV_mul(cyc,Q))));
    1805         697 :         if (lgcols(F) == 1) F = NULL;
    1806             :       }
    1807             :       /* reduce modulo units and Q^* */
    1808         740 :       if (lg(U) != 1) z = ZC_sub(z, ZM_ZC_mul(U, RgM_Babai(U,z)));
    1809         740 :       col = mkmat2(X, z);
    1810         740 :       if (F) col = famat_mul_shallow(col, F);
    1811         740 :       col = famat_remove_trivial(col);
    1812         740 :       if (xar) col = famat_mul_shallow(col, xar);
    1813             :     }
    1814         532 :     else if (!ZV_equal0(R))
    1815             :     { /* in case isprincipalfact calls bnfinit() due to prec trouble...*/
    1816         532 :       GEN y = isprincipalfact(bnf, x, gen, ZC_neg(R), flag);
    1817         532 :       if (typ(y) != t_VEC) return y;
    1818         532 :       col = gel(y,2);
    1819             :     }
    1820             :   }
    1821      240735 :   if (col)
    1822             :   { /* add back missing content */
    1823      241099 :     if (xc) col = (typ(col)==t_MAT)? famat_mul_shallow(col,xc)
    1824         448 :                                    : RgC_Rg_mul(col,xc);
    1825      240651 :     if (typ(col) != t_MAT && lg(col) != 1 && (flag & nf_GENMAT))
    1826      228333 :       col = to_famat_shallow(col, gen_1);
    1827             :   }
    1828             :   else
    1829             :   {
    1830          84 :     if (e < 0) e = 0;
    1831          84 :     *pprec += nbits2extraprec(e + 128);
    1832          84 :     if (flag & nf_FORCE)
    1833             :     {
    1834          70 :       if (DEBUGLEVEL)
    1835           0 :         pari_warn(warner,"precision too low for generators, e = %ld",e);
    1836          70 :       return NULL;
    1837             :     }
    1838          14 :     pari_warn(warner,"precision too low for generators, not given");
    1839          14 :     col = cgetg(1, t_COL);
    1840             :   }
    1841      240665 :   return (flag & nf_GEN_IF_PRINCIPAL)? col: mkvec2(R, col);
    1842             : }
    1843             : 
    1844             : static GEN
    1845       62377 : triv_gen(GEN bnf, GEN x, long flag)
    1846             : {
    1847       62377 :   pari_sp av = avma;
    1848       62377 :   GEN nf = bnf_get_nf(bnf);
    1849             :   long c;
    1850       62377 :   if (flag & nf_GEN_IF_PRINCIPAL)
    1851             :   {
    1852           7 :     if (!(flag & nf_GENMAT)) return algtobasis(nf,x);
    1853           7 :     x = nf_to_scalar_or_basis(nf,x);
    1854           7 :     if (typ(x) == t_INT && is_pm1(x)) return trivial_fact();
    1855           0 :     return gerepilecopy(av, to_famat_shallow(x, gen_1));
    1856             :   }
    1857       62370 :   c = lg(bnf_get_cyc(bnf)) - 1;
    1858       62370 :   if (flag & nf_GENMAT)
    1859       52850 :     retmkvec2(zerocol(c), to_famat_shallow(algtobasis(nf,x), gen_1));
    1860        9520 :   if (flag & nf_GEN)
    1861          21 :     retmkvec2(zerocol(c), algtobasis(nf,x));
    1862        9499 :   return zerocol(c);
    1863             : }
    1864             : 
    1865             : GEN
    1866      440740 : bnfisprincipal0(GEN bnf,GEN x,long flag)
    1867             : {
    1868      440740 :   pari_sp av = avma;
    1869             :   GEN arch, c, nf;
    1870             :   long pr;
    1871             : 
    1872      440740 :   bnf = checkbnf(bnf);
    1873      440740 :   nf = bnf_get_nf(bnf);
    1874      440740 :   switch( idealtyp(&x, &arch) )
    1875             :   {
    1876       50778 :     case id_PRINCIPAL:
    1877       50778 :       if (gequal0(x)) pari_err_DOMAIN("bnfisprincipal","ideal","=",gen_0,x);
    1878       50778 :       return triv_gen(bnf, x, flag);
    1879      378237 :     case id_PRIME:
    1880      378237 :       if (pr_is_inert(x)) return triv_gen(bnf, pr_get_p(x), flag);
    1881      366645 :       x = pr_hnf(nf, x);
    1882      366645 :       break;
    1883       11725 :     case id_MAT:
    1884       11725 :       if (lg(x)==1) pari_err_DOMAIN("bnfisprincipal","ideal","=",gen_0,x);
    1885       11725 :       if (nf_get_degree(nf) != lg(x)-1)
    1886           0 :         pari_err_TYPE("idealtyp [dimension != degree]", x);
    1887             :   }
    1888      378370 :   pr = prec_arch(bnf); /* precision of unit matrix */
    1889      378370 :   c = getrand();
    1890             :   for (;;)
    1891           0 :   {
    1892      378370 :     pari_sp av1 = avma;
    1893      378370 :     GEN y = isprincipalall(bnf,x,&pr,flag);
    1894      378370 :     if (y) return gerepilecopy(av, y);
    1895             : 
    1896           0 :     if (DEBUGLEVEL) pari_warn(warnprec,"isprincipal",pr);
    1897           0 :     set_avma(av1); bnf = bnfnewprec_shallow(bnf,pr); setrand(c);
    1898             :   }
    1899             : }
    1900             : GEN
    1901      162457 : isprincipal(GEN bnf,GEN x) { return bnfisprincipal0(bnf,x,0); }
    1902             : 
    1903             : /* FIXME: OBSOLETE */
    1904             : GEN
    1905           0 : isprincipalgen(GEN bnf,GEN x)
    1906           0 : { return bnfisprincipal0(bnf,x,nf_GEN); }
    1907             : GEN
    1908           0 : isprincipalforce(GEN bnf,GEN x)
    1909           0 : { return bnfisprincipal0(bnf,x,nf_FORCE); }
    1910             : GEN
    1911           0 : isprincipalgenforce(GEN bnf,GEN x)
    1912           0 : { return bnfisprincipal0(bnf,x,nf_GEN | nf_FORCE); }
    1913             : 
    1914             : /* lg(u) > 1 */
    1915             : static int
    1916          77 : RgV_is1(GEN u) { return isint1(gel(u,1)) && RgV_isscalar(u); }
    1917             : static GEN
    1918       21805 : add_principal_part(GEN nf, GEN u, GEN v, long flag)
    1919             : {
    1920       21805 :   if (flag & nf_GENMAT)
    1921        8232 :     return (typ(u) == t_COL && RgV_is1(u))? v: famat_mul_shallow(v,u);
    1922             :   else
    1923       13573 :     return nfmul(nf, v, u);
    1924             : }
    1925             : 
    1926             : #if 0
    1927             : /* compute C prod P[i]^e[i],  e[i] >=0 for all i. C may be NULL (omitted)
    1928             :  * e destroyed ! */
    1929             : static GEN
    1930             : expand(GEN nf, GEN C, GEN P, GEN e)
    1931             : {
    1932             :   long i, l = lg(e), done = 1;
    1933             :   GEN id = C;
    1934             :   for (i=1; i<l; i++)
    1935             :   {
    1936             :     GEN ei = gel(e,i);
    1937             :     if (signe(ei))
    1938             :     {
    1939             :       if (mod2(ei)) id = id? idealmul(nf, id, gel(P,i)): gel(P,i);
    1940             :       ei = shifti(ei,-1);
    1941             :       if (signe(ei)) done = 0;
    1942             :       gel(e,i) = ei;
    1943             :     }
    1944             :   }
    1945             :   if (id != C) id = idealred(nf, id);
    1946             :   if (done) return id;
    1947             :   return idealmulred(nf, id, idealsqr(nf, expand(nf,id,P,e)));
    1948             : }
    1949             : /* C is an extended ideal, possibly with C[1] = NULL */
    1950             : static GEN
    1951             : expandext(GEN nf, GEN C, GEN P, GEN e)
    1952             : {
    1953             :   long i, l = lg(e), done = 1;
    1954             :   GEN A = gel(C,1);
    1955             :   for (i=1; i<l; i++)
    1956             :   {
    1957             :     GEN ei = gel(e,i);
    1958             :     if (signe(ei))
    1959             :     {
    1960             :       if (mod2(ei)) A = A? idealmul(nf, A, gel(P,i)): gel(P,i);
    1961             :       ei = shifti(ei,-1);
    1962             :       if (signe(ei)) done = 0;
    1963             :       gel(e,i) = ei;
    1964             :     }
    1965             :   }
    1966             :   if (A == gel(C,1))
    1967             :     A = C;
    1968             :   else
    1969             :     A = idealred(nf, mkvec2(A, gel(C,2)));
    1970             :   if (done) return A;
    1971             :   return idealmulred(nf, A, idealsqr(nf, expand(nf,A,P,e)));
    1972             : }
    1973             : #endif
    1974             : 
    1975             : static GEN
    1976           0 : expand(GEN nf, GEN C, GEN P, GEN e)
    1977             : {
    1978           0 :   long i, l = lg(e);
    1979           0 :   GEN B, A = C;
    1980           0 :   for (i=1; i<l; i++) /* compute prod P[i]^e[i] */
    1981           0 :     if (signe(gel(e,i)))
    1982             :     {
    1983           0 :       B = idealpowred(nf, gel(P,i), gel(e,i));
    1984           0 :       A = A? idealmulred(nf,A,B): B;
    1985             :     }
    1986           0 :   return A;
    1987             : }
    1988             : static GEN
    1989       21819 : expandext(GEN nf, GEN C, GEN P, GEN e)
    1990             : {
    1991       21819 :   long i, l = lg(e);
    1992       21819 :   GEN B, A = gel(C,1), C1 = A;
    1993       72251 :   for (i=1; i<l; i++) /* compute prod P[i]^e[i] */
    1994       50432 :     if (signe(gel(e,i)))
    1995             :     {
    1996       28361 :       gel(C,1) = gel(P,i);
    1997       28361 :       B = idealpowred(nf, C, gel(e,i));
    1998       28361 :       A = A? idealmulred(nf,A,B): B;
    1999             :     }
    2000       21819 :   return A == C1? C: A;
    2001             : }
    2002             : 
    2003             : /* isprincipal for C * \prod P[i]^e[i] (C omitted if NULL) */
    2004             : GEN
    2005       21819 : isprincipalfact(GEN bnf, GEN C, GEN P, GEN e, long flag)
    2006             : {
    2007       21819 :   const long gen = flag & (nf_GEN|nf_GENMAT|nf_GEN_IF_PRINCIPAL);
    2008             :   long prec;
    2009       21819 :   pari_sp av = avma;
    2010       21819 :   GEN C0, Cext, c, id, nf = checknf(bnf);
    2011             : 
    2012       21819 :   if (gen)
    2013             :   {
    2014        8239 :     Cext = (flag & nf_GENMAT)? trivial_fact()
    2015       21819 :                              : mkpolmod(gen_1,nf_get_pol(nf));
    2016       21819 :     C0 = mkvec2(C, Cext);
    2017       21819 :     id = expandext(nf, C0, P, e);
    2018             :   } else {
    2019           0 :     Cext = NULL;
    2020           0 :     C0 = C;
    2021           0 :     id = expand(nf, C, P, e);
    2022             :   }
    2023       21819 :   if (id == C0) /* e = 0 */
    2024             :   {
    2025        8151 :     if (!C) return bnfisprincipal0(bnf, gen_1, flag);
    2026        8144 :     switch(typ(C))
    2027             :     {
    2028           7 :       case t_INT: case t_FRAC: case t_POL: case t_POLMOD: case t_COL:
    2029           7 :         return triv_gen(bnf, C, flag);
    2030             :     }
    2031        8137 :     C = idealhnf_shallow(nf,C);
    2032             :   }
    2033             :   else
    2034             :   {
    2035       13668 :     if (gen) { C = gel(id,1); Cext = gel(id,2); } else C = id;
    2036             :   }
    2037       21805 :   prec = prec_arch(bnf);
    2038       21805 :   c = getrand();
    2039             :   for (;;)
    2040          70 :   {
    2041       21875 :     pari_sp av1 = avma;
    2042       21875 :     GEN y = isprincipalall(bnf, C, &prec, flag);
    2043       21875 :     if (y)
    2044             :     {
    2045       21805 :       if (flag & nf_GEN_IF_PRINCIPAL)
    2046             :       {
    2047       18235 :         if (typ(y) == t_INT) return gc_NULL(av);
    2048       18235 :         y = add_principal_part(nf, y, Cext, flag);
    2049             :       }
    2050             :       else
    2051             :       {
    2052        3570 :         GEN u = gel(y,2);
    2053        3570 :         if (!gen || typ(y) != t_VEC) return gerepileupto(av,y);
    2054        3570 :         if (lg(u) != 1) gel(y,2) = add_principal_part(nf, u, Cext, flag);
    2055             :       }
    2056       21805 :       return gerepilecopy(av, y);
    2057             :     }
    2058          70 :     if (DEBUGLEVEL) pari_warn(warnprec,"isprincipal",prec);
    2059          70 :     set_avma(av1); bnf = bnfnewprec_shallow(bnf,prec); setrand(c);
    2060             :   }
    2061             : }
    2062             : GEN
    2063           0 : isprincipalfact_or_fail(GEN bnf, GEN C, GEN P, GEN e)
    2064             : {
    2065           0 :   const long flag = nf_GENMAT|nf_FORCE;
    2066             :   long prec;
    2067           0 :   pari_sp av = avma;
    2068           0 :   GEN u, y, id, C0, Cext, nf = bnf_get_nf(bnf);
    2069             : 
    2070           0 :   Cext = trivial_fact();
    2071           0 :   C0 = mkvec2(C, Cext);
    2072           0 :   id = expandext(nf, C0, P, e);
    2073           0 :   if (id == C0) /* e = 0 */
    2074           0 :     C = idealhnf_shallow(nf,C);
    2075             :   else {
    2076           0 :     C = gel(id,1); Cext = gel(id,2);
    2077             :   }
    2078           0 :   prec = prec_arch(bnf);
    2079           0 :   y = isprincipalall(bnf, C, &prec, flag);
    2080           0 :   if (!y) { set_avma(av); return utoipos(prec); }
    2081           0 :   u = gel(y,2);
    2082           0 :   if (lg(u) != 1) gel(y,2) = add_principal_part(nf, u, Cext, flag);
    2083           0 :   return gerepilecopy(av, y);
    2084             : }
    2085             : 
    2086             : GEN
    2087       24577 : nfsign_from_logarch(GEN LA, GEN invpi, GEN archp)
    2088             : {
    2089       24577 :   long l = lg(archp), i;
    2090       24577 :   GEN y = cgetg(l, t_VECSMALL);
    2091       24577 :   pari_sp av = avma;
    2092             : 
    2093       55552 :   for (i=1; i<l; i++)
    2094             :   {
    2095       30975 :     GEN c = ground( gmul(imag_i(gel(LA,archp[i])), invpi) );
    2096       30975 :     y[i] = mpodd(c)? 1: 0;
    2097             :   }
    2098       24577 :   set_avma(av); return y;
    2099             : }
    2100             : 
    2101             : GEN
    2102       37317 : nfsign_tu(GEN bnf, GEN archp)
    2103             : {
    2104             :   long n;
    2105       37317 :   if (bnf_get_tuN(bnf) != 2) return cgetg(1, t_VECSMALL);
    2106       33152 :   n = archp? lg(archp) - 1: nf_get_r1(bnf_get_nf(bnf));
    2107       33152 :   return const_vecsmall(n, 1);
    2108             : }
    2109             : GEN
    2110       38556 : nfsign_fu(GEN bnf, GEN archp)
    2111             : {
    2112       38556 :   GEN invpi, y, A = bnf_get_logfu(bnf), nf = bnf_get_nf(bnf);
    2113       38556 :   long j = 1, RU = lg(A);
    2114             : 
    2115       38556 :   if (!archp) archp = identity_perm( nf_get_r1(nf) );
    2116       38556 :   invpi = invr( mppi(nf_get_prec(nf)) );
    2117       38556 :   y = cgetg(RU,t_MAT);
    2118       63035 :   for (j = 1; j < RU; j++)
    2119       24479 :     gel(y,j) = nfsign_from_logarch(gel(A,j), invpi, archp);
    2120       38556 :   return y;
    2121             : }
    2122             : GEN
    2123          35 : nfsign_units(GEN bnf, GEN archp, int add_zu)
    2124             : {
    2125          35 :   GEN sfu = nfsign_fu(bnf, archp);
    2126          35 :   return add_zu? vec_prepend(sfu, nfsign_tu(bnf, archp)): sfu;
    2127             : }
    2128             : 
    2129             : /* obsolete */
    2130             : GEN
    2131           7 : signunits(GEN bnf)
    2132             : {
    2133             :   pari_sp av;
    2134             :   GEN S, y, nf;
    2135             :   long i, j, r1, r2;
    2136             : 
    2137           7 :   bnf = checkbnf(bnf); nf = bnf_get_nf(bnf);
    2138           7 :   nf_get_sign(nf, &r1,&r2);
    2139           7 :   S = zeromatcopy(r1, r1+r2-1); av = avma;
    2140           7 :   y = nfsign_fu(bnf, NULL);
    2141          14 :   for (j = 1; j < lg(y); j++)
    2142             :   {
    2143           7 :     GEN Sj = gel(S,j), yj = gel(y,j);
    2144          21 :     for (i = 1; i <= r1; i++) gel(Sj,i) = yj[i]? gen_m1: gen_1;
    2145             :   }
    2146           7 :   set_avma(av); return S;
    2147             : }
    2148             : 
    2149             : static GEN
    2150      162679 : get_log_embed(REL_t *rel, GEN M, long RU, long R1, long prec)
    2151             : {
    2152      162679 :   GEN arch, C, z = rel->m;
    2153             :   long i;
    2154      162679 :   arch = typ(z) == t_COL? RgM_RgC_mul(M, z): const_col(nbrows(M), z);
    2155      162679 :   C = cgetg(RU+1, t_COL); arch = glog(arch, prec);
    2156      357492 :   for (i=1; i<=R1; i++) gel(C,i) = gel(arch,i);
    2157      355832 :   for (   ; i<=RU; i++) gel(C,i) = gmul2n(gel(arch,i), 1);
    2158      162679 :   return C;
    2159             : }
    2160             : static GEN
    2161      243129 : rel_embed(REL_t *rel, FB_t *F, GEN embs, long ind, GEN M, long RU, long R1,
    2162             :           long prec)
    2163             : {
    2164             :   GEN C, D, perm;
    2165             :   long i, n;
    2166      243129 :   if (!rel->relaut) return get_log_embed(rel, M, RU, R1, prec);
    2167             :   /* image of another relation by automorphism */
    2168       80450 :   C = gel(embs, ind - rel->relorig);
    2169       80450 :   perm = gel(F->embperm, rel->relaut);
    2170       80450 :   D = cgetg_copy(C, &n);
    2171      355493 :   for (i = 1; i < n; i++)
    2172             :   {
    2173      275043 :     long v = perm[i];
    2174      275043 :     gel(D,i) = (v > 0)? gel(C,v): conj_i(gel(C,-v));
    2175             :   }
    2176       80450 :   return D;
    2177             : }
    2178             : static GEN
    2179       31266 : get_embs(FB_t *F, RELCACHE_t *cache, GEN nf, GEN embs, long PREC)
    2180             : {
    2181       31266 :   long ru, j, k, l = cache->last - cache->chk + 1, r1 = nf_get_r1(nf);
    2182       31266 :   GEN M = nf_get_M(nf), nembs = cgetg(cache->last - cache->base+1, t_MAT);
    2183             :   REL_t *rel;
    2184             : 
    2185     3216974 :   for (k = 1; k <= cache->chk - cache->base; k++) gel(nembs,k) = gel(embs,k);
    2186       31266 :   embs = nembs; ru = nbrows(M);
    2187      264397 :   for (j=1,rel = cache->chk + 1; j < l; rel++,j++,k++)
    2188      233131 :     gel(embs,k) = rel_embed(rel, F, embs, k, M, ru, r1, PREC);
    2189       31266 :   return embs;
    2190             : }
    2191             : static void
    2192       78800 : set_rel_alpha(REL_t *rel, GEN auts, GEN vA, long ind)
    2193             : {
    2194             :   GEN u;
    2195       78800 :   if (!rel->relaut)
    2196       45744 :     u = rel->m;
    2197             :   else
    2198       33056 :     u = ZM_ZC_mul(gel(auts, rel->relaut), gel(vA, ind - rel->relorig));
    2199       78800 :   gel(vA, ind) = u;
    2200       78800 : }
    2201             : static GEN
    2202     1622788 : set_fact(FB_t *F, FACT *fact, GEN e, long *pnz)
    2203             : {
    2204     1622788 :   long n = fact[0].pr;
    2205     1622788 :   GEN c = zero_Flv(F->KC);
    2206     1622788 :   if (!n) /* trivial factorization */
    2207           0 :     *pnz = F->KC+1;
    2208             :   else
    2209             :   {
    2210     1622788 :     long i, nz = minss(fact[1].pr, fact[n].pr);
    2211     7913602 :     for (i = 1; i <= n; i++) c[fact[i].pr] = fact[i].ex;
    2212     1622788 :     if (e)
    2213             :     {
    2214       12312 :       long l = lg(e);
    2215       57723 :       for (i = 1; i < l; i++)
    2216       45411 :         if (e[i]) { long v = F->subFB[i]; c[v] += e[i]; if (v < nz) nz = v; }
    2217             :     }
    2218     1622788 :     *pnz = nz;
    2219             :   }
    2220     1622788 :   return c;
    2221             : }
    2222             : 
    2223             : /* Is cols already in the cache ? bs = index of first non zero coeff in cols
    2224             :  * General check for colinearity useless since exceedingly rare */
    2225             : static int
    2226     1815044 : already_known(RELCACHE_t *cache, long bs, GEN cols)
    2227             : {
    2228             :   REL_t *r;
    2229     1815044 :   long l = lg(cols);
    2230   173474730 :   for (r = cache->last; r > cache->base; r--)
    2231   171833911 :     if (bs == r->nz)
    2232             :     {
    2233     9759538 :       GEN coll = r->R;
    2234     9759538 :       long b = bs;
    2235   107687084 :       while (b < l && cols[b] == coll[b]) b++;
    2236     9759538 :       if (b == l) return 1;
    2237             :     }
    2238     1640819 :   return 0;
    2239             : }
    2240             : 
    2241             : /* Add relation R to cache, nz = index of first non zero coeff in R.
    2242             :  * If relation is a linear combination of the previous ones, return 0.
    2243             :  * Otherwise, update basis and return > 0. Compute mod p (much faster)
    2244             :  * so some kernel vector might not be genuine. */
    2245             : static int
    2246     1815261 : add_rel_i(RELCACHE_t *cache, GEN R, long nz, GEN m, long orig, long aut, REL_t **relp, long in_rnd_rel)
    2247             : {
    2248     1815261 :   long i, k, n = lg(R)-1;
    2249             : 
    2250     1815261 :   if (nz == n+1) { k = 0; goto ADD_REL; }
    2251     1815044 :   if (already_known(cache, nz, R)) return -1;
    2252     1640819 :   if (cache->last >= cache->base + cache->len) return 0;
    2253     1640819 :   if (DEBUGLEVEL>6)
    2254             :   {
    2255           0 :     err_printf("adding vector = %Ps\n",R);
    2256           0 :     err_printf("generators =\n%Ps\n", cache->basis);
    2257             :   }
    2258     1640819 :   if (cache->missing)
    2259             :   {
    2260     1571333 :     GEN a = leafcopy(R), basis = cache->basis;
    2261     1571333 :     k = lg(a);
    2262    97212811 :     do --k; while (!a[k]);
    2263     5029207 :     while (k)
    2264             :     {
    2265     3579224 :       GEN c = gel(basis, k);
    2266     3579224 :       if (c[k])
    2267             :       {
    2268     3457874 :         long ak = a[k];
    2269   300038412 :         for (i=1; i < k; i++) if (c[i]) a[i] = (a[i] + ak*(mod_p-c[i])) % mod_p;
    2270     3457874 :         a[k] = 0;
    2271   171212030 :         do --k; while (!a[k]); /* k cannot go below 0: codeword is a sentinel */
    2272             :       }
    2273             :       else
    2274             :       {
    2275      121350 :         ulong invak = Fl_inv(uel(a,k), mod_p);
    2276             :         /* Cleanup a */
    2277     7343467 :         for (i = k; i-- > 1; )
    2278             :         {
    2279     7222117 :           long j, ai = a[i];
    2280     7222117 :           c = gel(basis, i);
    2281     7222117 :           if (!ai || !c[i]) continue;
    2282       99187 :           ai = mod_p-ai;
    2283     3447626 :           for (j = 1; j < i; j++) if (c[j]) a[j] = (a[j] + ai*c[j]) % mod_p;
    2284       99187 :           a[i] = 0;
    2285             :         }
    2286             :         /* Insert a/a[k] as k-th column */
    2287      121350 :         c = gel(basis, k);
    2288     7343467 :         for (i = 1; i<k; i++) if (a[i]) c[i] = (a[i] * invak) % mod_p;
    2289      121350 :         c[k] = 1; a = c;
    2290             :         /* Cleanup above k */
    2291     7263843 :         for (i = k+1; i<n; i++)
    2292             :         {
    2293             :           long j, ck;
    2294     7142493 :           c = gel(basis, i);
    2295     7142493 :           ck = c[k];
    2296     7142493 :           if (!ck) continue;
    2297     1060285 :           ck = mod_p-ck;
    2298    40136048 :           for (j = 1; j < k; j++) if (a[j]) c[j] = (c[j] + ck*a[j]) % mod_p;
    2299     1060285 :           c[k] = 0;
    2300             :         }
    2301      121350 :         cache->missing--;
    2302      121350 :         break;
    2303             :       }
    2304             :     }
    2305             :   }
    2306             :   else
    2307       69486 :     k = (cache->last - cache->base) + 1;
    2308     1640819 :   if (k || cache->relsup > 0 || (m && in_rnd_rel))
    2309             :   {
    2310             :     REL_t *rel;
    2311             : 
    2312      213248 : ADD_REL:
    2313      213465 :     rel = ++cache->last;
    2314      213465 :     if (!k && cache->relsup && nz < n+1)
    2315             :     {
    2316       21889 :       cache->relsup--;
    2317       21889 :       k = (rel - cache->base) + cache->missing;
    2318             :     }
    2319      213465 :     rel->R  = gclone(R);
    2320      213465 :     rel->m  =  m ? gclone(m) : NULL;
    2321      213465 :     rel->nz = nz;
    2322      213465 :     if (aut)
    2323             :     {
    2324       78489 :       rel->relorig = (rel - cache->base) - orig;
    2325       78489 :       rel->relaut = aut;
    2326             :     }
    2327             :     else
    2328      134976 :       rel->relaut = 0;
    2329      213465 :     if (relp) *relp = rel;
    2330      213465 :     if (DEBUGLEVEL) dbg_newrel(cache);
    2331             :   }
    2332     1641036 :   return k;
    2333             : }
    2334             : 
    2335             : static int
    2336     1668967 : add_rel(RELCACHE_t *cache, FB_t *F, GEN R, long nz, GEN m, long in_rnd_rel)
    2337             : {
    2338             :   REL_t *rel;
    2339             :   long k, l, reln;
    2340     1668967 :   const long lauts = lg(F->idealperm), KC = F->KC;
    2341             : 
    2342     1668967 :   k = add_rel_i(cache, R, nz, m, 0, 0, &rel, in_rnd_rel);
    2343     1668967 :   if (k > 0 && typ(m) != t_INT)
    2344             :   {
    2345       88274 :     GEN Rl = cgetg(KC+1, t_VECSMALL);
    2346       88274 :     reln = rel - cache->base;
    2347      234568 :     for (l = 1; l < lauts; l++)
    2348             :     {
    2349      146294 :       GEN perml = gel(F->idealperm, l);
    2350      146294 :       long i, nzl = perml[nz];
    2351             : 
    2352    12399334 :       for (i = 1; i <= KC; i++) Rl[i] = 0;
    2353    11161278 :       for (i = nz; i <= KC; i++)
    2354    11014984 :         if (R[i])
    2355             :         {
    2356      532713 :           long v = perml[i];
    2357             : 
    2358      532713 :           if (v < nzl) nzl = v;
    2359      532713 :           Rl[v] = R[i];
    2360             :         }
    2361      146294 :       (void)add_rel_i(cache, Rl, nzl, NULL, reln, l, NULL, in_rnd_rel);
    2362             :     }
    2363             :   }
    2364     1668967 :   return k;
    2365             : }
    2366             : 
    2367             : INLINE void
    2368    16711975 : step(GEN x, double *y, GEN inc, long k)
    2369             : {
    2370    16711975 :   if (!y[k])
    2371     8639060 :     x[k]++; /* leading coeff > 0 */
    2372             :   else
    2373             :   {
    2374     8072915 :     long i = inc[k];
    2375     8072915 :     x[k] += i;
    2376     8072915 :     inc[k] = (i > 0)? -1-i: 1-i;
    2377             :   }
    2378    16711975 : }
    2379             : 
    2380             : INLINE long
    2381     1918828 : Fincke_Pohst_ideal(RELCACHE_t *cache, FB_t *F, GEN nf, GEN M, GEN I,
    2382             :     GEN NI, FACT *fact, long Nrelid, FP_t *fp, RNDREL_t *rr, long prec,
    2383             :     long *Nsmall, long *Nfact)
    2384             : {
    2385             :   pari_sp av;
    2386     1918828 :   const long N = nf_get_degree(nf), R1 = nf_get_r1(nf);
    2387     1918828 :   GEN G = nf_get_G(nf), G0 = nf_get_roundG(nf), r, u, gx, inc, ideal;
    2388             :   double BOUND, B1, B2;
    2389     1918828 :   long j, k, skipfirst, relid=0, dependent=0, try_elt=0, try_factor=0;
    2390             : 
    2391     1918828 :   inc = const_vecsmall(N, 1);
    2392     1918828 :   u = ZM_lll(ZM_mul(G0, I), 0.99, LLL_IM);
    2393     1918828 :   ideal = ZM_mul(I,u); /* approximate T2-LLL reduction */
    2394     1918828 :   r = gaussred_from_QR(RgM_mul(G, ideal), prec); /* Cholesky for T2 | ideal */
    2395     1918828 :   if (!r) pari_err_BUG("small_norm (precision too low)");
    2396             : 
    2397     6721089 :   for (k=1; k<=N; k++)
    2398             :   {
    2399     4802261 :     fp->v[k] = gtodouble(gcoeff(r,k,k));
    2400    10759239 :     for (j=1; j<k; j++) fp->q[j][k] = gtodouble(gcoeff(r,j,k));
    2401     4802261 :     if (DEBUGLEVEL>3) err_printf("v[%ld]=%.4g ",k,fp->v[k]);
    2402             :   }
    2403     1918828 :   B1 = fp->v[1]; /* T2(ideal[1]) */
    2404     1918828 :   B2 = fp->v[2] + B1 * fp->q[1][2] * fp->q[1][2]; /* T2(ideal[2]) */
    2405     1918828 :   if (ZV_isscalar(gel(ideal,1))) /* probable */
    2406             :   {
    2407     1250233 :     skipfirst = 1;
    2408     1250233 :     BOUND = mindd(BMULT * B1, 2 * B2);
    2409             :   }
    2410             :   else
    2411             :   {
    2412      668595 :     BOUND = mindd(BMULT * B1, 2 * B2);
    2413      668595 :     skipfirst = 0;
    2414             :   }
    2415             :   /* BOUND at most BMULT fp->x smallest known vector */
    2416     1918828 :   if (DEBUGLEVEL>1)
    2417             :   {
    2418           0 :     if (DEBUGLEVEL>3) err_printf("\n");
    2419           0 :     err_printf("BOUND = %.4g\n",BOUND);
    2420             :   }
    2421     1918828 :   BOUND *= 1 + 1e-6;
    2422     1918828 :   k = N; fp->y[N] = fp->z[N] = 0; fp->x[N] = 0;
    2423     6590174 :   for (av = avma;; set_avma(av), step(fp->x,fp->y,inc,k))
    2424     4671346 :   {
    2425             :     GEN R;
    2426             :     long nz;
    2427             :     do
    2428             :     { /* look for primitive element of small norm, cf minim00 */
    2429    10114206 :       int fl = 0;
    2430             :       double p;
    2431    10114206 :       if (k > 1)
    2432             :       {
    2433     5442860 :         long l = k-1;
    2434     5442860 :         fp->z[l] = 0;
    2435    27122360 :         for (j=k; j<=N; j++) fp->z[l] += fp->q[l][j]*fp->x[j];
    2436     5442860 :         p = (double)fp->x[k] + fp->z[k];
    2437     5442860 :         fp->y[l] = fp->y[k] + p*p*fp->v[k];
    2438     5442860 :         if (l <= skipfirst && !fp->y[1]) fl = 1;
    2439     5442860 :         fp->x[l] = (long)floor(-fp->z[l] + 0.5);
    2440     5442860 :         k = l;
    2441             :       }
    2442     5257832 :       for(;; step(fp->x,fp->y,inc,k))
    2443             :       {
    2444    17274493 :         if (++try_elt > maxtry_ELEMENT) return 0;
    2445    15372038 :         if (!fl)
    2446             :         {
    2447    14121805 :           p = (double)fp->x[k] + fp->z[k];
    2448    14121805 :           if (fp->y[k] + p*p*fp->v[k] <= BOUND) break;
    2449             : 
    2450     6782797 :           step(fp->x,fp->y,inc,k);
    2451             : 
    2452     6782797 :           p = (double)fp->x[k] + fp->z[k];
    2453     6782797 :           if (fp->y[k] + p*p*fp->v[k] <= BOUND) break;
    2454             :         }
    2455     7141339 :         fl = 0; inc[k] = 1;
    2456     7141339 :         if (++k > N) return 0;
    2457             :       }
    2458     8230699 :     } while (k > 1);
    2459             : 
    2460             :     /* element complete */
    2461     8300414 :     if (zv_content(fp->x) !=1) continue; /* not primitive */
    2462     3688462 :     gx = ZM_zc_mul(ideal,fp->x);
    2463     3688462 :     if (ZV_isscalar(gx)) continue;
    2464     3667498 :     if (++try_factor > maxtry_FACT) return 0;
    2465             : 
    2466     3667463 :     if (!Nrelid)
    2467             :     {
    2468        1446 :       if (!factorgen(F,nf,I,NI,gx,fact)) continue;
    2469          15 :       return 1;
    2470             :     }
    2471     3666017 :     else if (rr)
    2472             :     {
    2473     1462699 :       if (!factorgen(F,nf,I,NI,gx,fact)) continue;
    2474       12312 :       add_to_fact(rr->jid, 1, fact);
    2475             :     }
    2476             :     else
    2477             :     {
    2478     2203318 :       GEN Nx, xembed = RgM_RgC_mul(M, gx);
    2479             :       long e;
    2480     2203318 :       if (Nsmall) (*Nsmall)++;
    2481     2203318 :       Nx = grndtoi(embed_norm(xembed, R1), &e);
    2482     2203318 :       if (e >= 0) {
    2483           0 :         if (DEBUGLEVEL > 1) err_printf("+");
    2484      599984 :         continue;
    2485             :       }
    2486     2203318 :       if (!can_factor(F, nf, NULL, gx, Nx, fact)) continue;
    2487             :     }
    2488             : 
    2489             :     /* smooth element */
    2490     1615646 :     R = set_fact(F, fact, rr ? rr->ex : NULL, &nz);
    2491             :     /* make sure we get maximal rank first, then allow all relations */
    2492     1615646 :     if (add_rel(cache, F, R, nz, gx, rr ? 1 : 0) <= 0)
    2493             :     { /* probably Q-dependent from previous ones: forget it */
    2494     1530977 :       if (DEBUGLEVEL>1) err_printf("*");
    2495     1530977 :       if (++dependent > maxtry_DEP) break;
    2496     1520981 :       continue;
    2497             :     }
    2498       84669 :     dependent = 0;
    2499       84669 :     if (DEBUGLEVEL && Nfact) (*Nfact)++;
    2500       84669 :     if (cache->last >= cache->end) return 1; /* we have enough */
    2501       65771 :     if (++relid == Nrelid) break;
    2502             :   }
    2503       16373 :   return 0;
    2504             : }
    2505             : 
    2506             : static void
    2507       52397 : small_norm(RELCACHE_t *cache, FB_t *F, GEN nf, long Nrelid, GEN M,
    2508             :            FACT *fact, GEN p0)
    2509             : {
    2510       52397 :   const long prec = nf_get_prec(nf);
    2511             :   FP_t fp;
    2512             :   pari_sp av;
    2513       52397 :   GEN L_jid = F->L_jid, Np0;
    2514       52397 :   long Nsmall, Nfact, n = lg(L_jid);
    2515             :   pari_timer T;
    2516             : 
    2517       52397 :   if (DEBUGLEVEL)
    2518             :   {
    2519           0 :     timer_start(&T);
    2520           0 :     err_printf("#### Look for %ld relations in %ld ideals (small_norm)\n",
    2521           0 :                cache->end - cache->last, lg(L_jid)-1);
    2522           0 :     if (p0) err_printf("Look in p0 = %Ps\n", vecslice(p0,1,4));
    2523             :   }
    2524       52397 :   Nsmall = Nfact = 0;
    2525       52397 :   minim_alloc(lg(M), &fp.q, &fp.x, &fp.y, &fp.z, &fp.v);
    2526       52397 :   Np0 = p0? pr_norm(p0): NULL;
    2527     1468638 :   for (av = avma; --n; set_avma(av))
    2528             :   {
    2529     1430103 :     long j = L_jid[n];
    2530     1430103 :     GEN id = gel(F->LP, j), Nid;
    2531     1430103 :     if (DEBUGLEVEL>1)
    2532           0 :       err_printf("\n*** Ideal no %ld: %Ps\n", j, vecslice(id,1,4));
    2533     1430103 :     if (p0)
    2534     1365717 :     { Nid = mulii(Np0, pr_norm(id)); id = idealmul(nf, p0, id); }
    2535             :     else
    2536       64386 :     { Nid = pr_norm(id); id = pr_hnf(nf, id);}
    2537     1430103 :     if (Fincke_Pohst_ideal(cache, F, nf, M, id, Nid, fact, Nrelid, &fp,
    2538       13862 :                            NULL, prec, &Nsmall, &Nfact)) break;
    2539             :   }
    2540       52397 :   if (DEBUGLEVEL && Nsmall)
    2541             :   {
    2542           0 :     if (DEBUGLEVEL == 1)
    2543           0 :     { if (Nfact) err_printf("\n"); }
    2544             :     else
    2545           0 :       err_printf("  \nnb. fact./nb. small norm = %ld/%ld = %.3f\n",
    2546           0 :                   Nfact,Nsmall,((double)Nfact)/Nsmall);
    2547           0 :     if (timer_get(&T)>1) timer_printf(&T,"small_norm");
    2548             :   }
    2549       52397 : }
    2550             : 
    2551             : static GEN
    2552       25465 : get_random_ideal(FB_t *F, GEN nf, GEN ex)
    2553             : {
    2554       25465 :   long i, l = lg(ex);
    2555             :   for (;;)
    2556           1 :   {
    2557       25466 :     GEN I = NULL;
    2558      141067 :     for (i = 1; i < l; i++)
    2559      115601 :       if ((ex[i] = random_bits(RANDOM_BITS)))
    2560             :       {
    2561      108585 :         GEN pr = gel(F->LP, F->subFB[i]), e = utoipos(ex[i]);
    2562      108585 :         I = I? idealmulpowprime(nf, I, pr, e): idealpow(nf, pr, e);
    2563             :       }
    2564       25466 :     if (I && !ZM_isscalar(I,NULL)) return I; /* != (n)Z_K */
    2565             :   }
    2566             : }
    2567             : 
    2568             : static void
    2569       25465 : rnd_rel(RELCACHE_t *cache, FB_t *F, GEN nf, FACT *fact)
    2570             : {
    2571             :   pari_timer T;
    2572       25465 :   GEN L_jid = F->L_jid, M = nf_get_M(nf), R, NR;
    2573       25465 :   long i, l = lg(L_jid), prec = nf_get_prec(nf);
    2574             :   RNDREL_t rr;
    2575             :   FP_t fp;
    2576             :   pari_sp av;
    2577             : 
    2578       25465 :   if (DEBUGLEVEL) {
    2579           0 :     timer_start(&T);
    2580           0 :     err_printf("\n#### Look for %ld relations in %ld ideals (rnd_rel)\n",
    2581           0 :                cache->end - cache->last, l-1);
    2582             :   }
    2583       25465 :   rr.ex = cgetg(lg(F->subFB), t_VECSMALL);
    2584       25465 :   R = get_random_ideal(F, nf, rr.ex); /* random product from subFB */
    2585       25465 :   NR = ZM_det_triangular(R);
    2586       25465 :   minim_alloc(lg(M), &fp.q, &fp.x, &fp.y, &fp.z, &fp.v);
    2587      508032 :   for (av = avma, i = 1; i < l; i++, set_avma(av))
    2588             :   { /* try P[j] * base */
    2589      487603 :     long j = L_jid[i];
    2590      487603 :     GEN P = gel(F->LP, j), Nid = mulii(NR, pr_norm(P));
    2591      487603 :     if (DEBUGLEVEL>1) err_printf("\n*** Ideal %ld: %Ps\n", j, vecslice(P,1,4));
    2592      487603 :     rr.jid = j;
    2593      487603 :     if (Fincke_Pohst_ideal(cache, F, nf, M, idealHNF_mul(nf, R, P), Nid, fact,
    2594        5036 :                            RND_REL_RELPID, &fp, &rr, prec, NULL, NULL)) break;
    2595             :   }
    2596       25465 :   if (DEBUGLEVEL)
    2597             :   {
    2598           0 :     err_printf("\n");
    2599           0 :     if (timer_get(&T) > 1) timer_printf(&T,"for remaining ideals");
    2600             :   }
    2601       25465 : }
    2602             : 
    2603             : static GEN
    2604       12040 : automorphism_perms(GEN M, GEN auts, GEN cyclic, long r1, long r2, long N)
    2605             : {
    2606       12040 :   long L = lgcols(M), lauts = lg(auts), lcyc = lg(cyclic), i, j, l, m;
    2607       12040 :   GEN Mt, perms = cgetg(lauts, t_VEC);
    2608             :   pari_sp av;
    2609             : 
    2610       27006 :   for (l = 1; l < lauts; l++) gel(perms, l) = cgetg(L, t_VECSMALL);
    2611       12040 :   av = avma;
    2612       12040 :   Mt = shallowtrans(gprec_w(M, LOWDEFAULTPREC));
    2613       12040 :   Mt = shallowconcat(Mt, conj_i(vecslice(Mt, r1+1, r1+r2)));
    2614       25536 :   for (l = 1; l < lcyc; l++)
    2615             :   {
    2616       13496 :     GEN thiscyc = gel(cyclic, l), thisperm, perm, prev, Nt;
    2617       13496 :     long k = thiscyc[1];
    2618             : 
    2619       13496 :     Nt = RgM_mul(shallowtrans(gel(auts, k)), Mt);
    2620       13496 :     perm = gel(perms, k);
    2621       38325 :     for (i = 1; i < L; i++)
    2622             :     {
    2623       24829 :       GEN v = gel(Nt, i), minD;
    2624       24829 :       minD = gnorml2(gsub(v, gel(Mt, 1)));
    2625       24829 :       perm[i] = 1;
    2626      125804 :       for (j = 2; j <= N; j++)
    2627             :       {
    2628      100975 :         GEN D = gnorml2(gsub(v, gel(Mt, j)));
    2629      100975 :         if (gcmp(D, minD) < 0) { minD = D; perm[i] = j >= L ? r2-j : j; }
    2630             :       }
    2631             :     }
    2632       15148 :     for (prev = perm, m = 2; m < lg(thiscyc); m++, prev = thisperm)
    2633             :     {
    2634        1652 :       thisperm = gel(perms, thiscyc[m]);
    2635       11004 :       for (i = 1; i < L; i++)
    2636             :       {
    2637        9352 :         long pp = labs(prev[i]);
    2638        9352 :         thisperm[i] = prev[i] < 0 ? -perm[pp] : perm[pp];
    2639             :       }
    2640             :     }
    2641             :   }
    2642       12040 :   set_avma(av); return perms;
    2643             : }
    2644             : 
    2645             : /* Determine the field automorphisms as matrices on the integral basis */
    2646             : static GEN
    2647       12096 : automorphism_matrices(GEN nf, GEN *cycp)
    2648             : {
    2649       12096 :   pari_sp av = avma;
    2650       12096 :   GEN auts = galoisconj(nf, NULL), mats, cyclic, cyclicidx;
    2651       12096 :   long nauts = lg(auts)-1, i, j, k, l;
    2652             : 
    2653       12096 :   cyclic = cgetg(nauts+1, t_VEC);
    2654       12096 :   cyclicidx = zero_Flv(nauts);
    2655       22862 :   for (l = 1; l <= nauts; l++)
    2656             :   {
    2657       22862 :     GEN aut = gel(auts, l);
    2658       22862 :     if (gequalX(aut)) { swap(gel(auts, l), gel(auts, nauts)); break; }
    2659             :   }
    2660             :   /* trivial automorphism is last */
    2661       39186 :   for (l = 1; l <= nauts; l++) gel(auts, l) = algtobasis(nf, gel(auts, l));
    2662             :   /* Compute maximal cyclic subgroups */
    2663       27090 :   for (l = nauts; --l > 0; ) if (!cyclicidx[l])
    2664             :   {
    2665       13713 :     GEN elt = gel(auts, l), aut = elt, cyc = cgetg(nauts+1, t_VECSMALL);
    2666       13713 :     cyc[1] = cyclicidx[l] = l; j = 1;
    2667             :     do
    2668             :     {
    2669       15414 :       elt = galoisapply(nf, elt, aut);
    2670       48489 :       for (k = 1; k <= nauts; k++) if (gequal(elt, gel(auts, k))) break;
    2671       15414 :       cyclicidx[k] = l; cyc[++j] = k;
    2672             :     }
    2673       15414 :     while (k != nauts);
    2674       13713 :     setlg(cyc, j);
    2675       13713 :     gel(cyclic, l) = cyc;
    2676             :   }
    2677       27090 :   for (i = j = 1; i < nauts; i++)
    2678       14994 :     if (cyclicidx[i] == i) cyclic[j++] = cyclic[i];
    2679       12096 :   setlg(cyclic, j);
    2680       12096 :   mats = cgetg(nauts, t_VEC);
    2681       25620 :   while (--j > 0)
    2682             :   {
    2683       13524 :     GEN cyc = gel(cyclic, j);
    2684       13524 :     long id = cyc[1];
    2685       13524 :     GEN M, Mi, aut = gel(auts, id);
    2686             : 
    2687       13524 :     gel(mats, id) = Mi = M = nfgaloismatrix(nf, aut);
    2688       15176 :     for (i = 2; i < lg(cyc); i++) gel(mats, cyc[i]) = Mi = ZM_mul(Mi, M);
    2689             :   }
    2690       12096 :   gerepileall(av, 2, &mats, &cyclic);
    2691       12096 :   if (cycp) *cycp = cyclic;
    2692       12096 :   return mats;
    2693             : }
    2694             : 
    2695             : /* vP a list of maximal ideals above the same p from idealprimedec: f(P/p) is
    2696             :  * increasing; 1 <= j <= #vP; orbit a zc of length <= #vP; auts a vector of
    2697             :  * automorphisms in ZM form.
    2698             :  * Set orbit[i] = 1 for all vP[i], i >= j, in the orbit of pr = vP[j] wrt auts.
    2699             :  * N.B.1 orbit need not be initialized to 0: useful to incrementally run
    2700             :  * through successive orbits
    2701             :  * N.B.2 i >= j, so primes with index < j will be missed; run incrementally
    2702             :  * starting from j = 1 ! */
    2703             : static void
    2704       11887 : pr_orbit_fill(GEN orbit, GEN auts, GEN vP, long j)
    2705             : {
    2706       11887 :   GEN pr = gel(vP,j), gen = pr_get_gen(pr);
    2707       11887 :   long i, l = lg(auts), J = lg(orbit), f = pr_get_f(pr);
    2708       11887 :   orbit[j] = 1;
    2709       23774 :   for (i = 1; i < l; i++)
    2710             :   {
    2711       11887 :     GEN g = ZM_ZC_mul(gel(auts,i), gen);
    2712             :     long k;
    2713       11894 :     for (k = j+1; k < J; k++)
    2714             :     {
    2715          21 :       GEN prk = gel(vP,k);
    2716          21 :       if (pr_get_f(prk) > f) break; /* f(P[k]) increases with k */
    2717             :       /* don't check that e matches: (almost) always 1 ! */
    2718          21 :       if (!orbit[k] && ZC_prdvd(g, prk)) { orbit[k] = 1; break; }
    2719             :     }
    2720             :   }
    2721       11887 : }
    2722             : /* remark: F->KCZ changes if be_honest() fails */
    2723             : static int
    2724          28 : be_honest(FB_t *F, GEN nf, GEN auts, FACT *fact)
    2725             : {
    2726             :   long i, iz, nbtest;
    2727          28 :   long lgsub = lg(F->subFB), KCZ0 = F->KCZ;
    2728          28 :   long N = nf_get_degree(nf), prec = nf_get_prec(nf);
    2729          28 :   GEN M = nf_get_M(nf);
    2730             :   FP_t fp;
    2731             :   pari_sp av;
    2732             : 
    2733          28 :   if (DEBUGLEVEL) {
    2734           0 :     err_printf("Be honest for %ld primes from %ld to %ld\n", F->KCZ2 - F->KCZ,
    2735           0 :                F->FB[ F->KCZ+1 ], F->FB[ F->KCZ2 ]);
    2736             :   }
    2737          28 :   minim_alloc(N+1, &fp.q, &fp.x, &fp.y, &fp.z, &fp.v);
    2738          28 :   if (lg(auts) == 1) auts = NULL;
    2739          28 :   av = avma;
    2740          36 :   for (iz=F->KCZ+1; iz<=F->KCZ2; iz++, set_avma(av))
    2741             :   {
    2742          29 :     long p = F->FB[iz];
    2743          29 :     GEN pr_orbit, P = F->LV[p];
    2744          29 :     long j, J = lg(P); /* > 1 */
    2745             :     /* the P|p, NP > C2 are assumed in subgroup generated by FB + last P
    2746             :      * with NP <= C2 is unramified --> check all but last */
    2747          29 :     if (pr_get_e(gel(P,J-1)) == 1) J--;
    2748          29 :     if (J == 1) continue;
    2749          29 :     if (DEBUGLEVEL>1) err_printf("%ld ", p);
    2750          29 :     pr_orbit = auts? zero_zv(J-1): NULL;
    2751          51 :     for (j = 1; j < J; j++)
    2752             :     {
    2753             :       GEN Nid, id, id0;
    2754          43 :       if (pr_orbit)
    2755             :       {
    2756          43 :         if (pr_orbit[j]) continue;
    2757             :         /* discard all primes in automorphism orbit simultaneously */
    2758          36 :         pr_orbit_fill(pr_orbit, auts, P, j);
    2759             :       }
    2760          36 :       id = id0 = pr_hnf(nf,gel(P,j));
    2761          36 :       Nid = pr_norm(gel(P,j));
    2762          36 :       for (nbtest=0;;)
    2763             :       {
    2764        1122 :         if (Fincke_Pohst_ideal(NULL, F, nf, M, id, Nid, fact, 0, &fp,
    2765          15 :                                NULL, prec, NULL, NULL)) break;
    2766        1107 :         if (++nbtest > maxtry_HONEST)
    2767             :         {
    2768          21 :           if (DEBUGLEVEL)
    2769           0 :             pari_warn(warner,"be_honest() failure on prime %Ps\n", gel(P,j));
    2770          21 :           return 0;
    2771             :         }
    2772             :         /* occurs at most once in the whole function */
    2773        6830 :         for (i = 1, id = id0; i < lgsub; i++)
    2774             :         {
    2775        5744 :           long ex = random_bits(RANDOM_BITS);
    2776        5744 :           if (ex)
    2777             :           {
    2778        5417 :             GEN pr = gel(F->LP, F->subFB[i]);
    2779        5417 :             id = idealmulpowprime(nf, id, pr, utoipos(ex));
    2780             :           }
    2781             :         }
    2782        1086 :         if (!equali1(gcoeff(id,N,N))) id = Q_primpart(id);
    2783        1086 :         if (expi(gcoeff(id,1,1)) > 100) id = idealred(nf, id);
    2784        1086 :         Nid = ZM_det_triangular(id);
    2785             :       }
    2786             :     }
    2787           8 :     F->KCZ++; /* SUCCESS, "enlarge" factorbase */
    2788             :   }
    2789           7 :   F->KCZ = KCZ0; return gc_bool(av,1);
    2790             : }
    2791             : 
    2792             : /* all primes with N(P) <= BOUND factor on factorbase ? */
    2793             : void
    2794          56 : bnftestprimes(GEN bnf, GEN BOUND)
    2795             : {
    2796          56 :   pari_sp av0 = avma, av;
    2797          56 :   ulong count = 0;
    2798          56 :   GEN auts, p, nf = bnf_get_nf(bnf), Vbase = bnf_get_vbase(bnf);
    2799          56 :   GEN fb = gen_sort_shallow(Vbase, (void*)&cmp_prime_ideal, cmp_nodata);
    2800          56 :   ulong pmax = pr_get_smallp(gel(fb, lg(fb)-1)); /*largest p in factorbase*/
    2801             :   forprime_t S;
    2802             :   FACT *fact;
    2803             :   FB_t F;
    2804             : 
    2805          56 :   (void)recover_partFB(&F, Vbase, nf_get_degree(nf));
    2806          56 :   fact = (FACT*)stack_malloc((F.KC+1)*sizeof(FACT));
    2807          56 :   forprime_init(&S, gen_2, BOUND);
    2808          56 :   auts = automorphism_matrices(nf, NULL);
    2809          56 :   if (lg(auts) == 1) auts = NULL;
    2810          56 :   av = avma;
    2811       37226 :   while (( p = forprime_next(&S) ))
    2812             :   {
    2813             :     GEN pr_orbit, vP;
    2814             :     long j, J;
    2815             : 
    2816       37170 :     if (DEBUGLEVEL == 1 && ++count > 1000)
    2817             :     {
    2818           0 :       err_printf("passing p = %Ps / %Ps\n", p, BOUND);
    2819           0 :       count = 0;
    2820             :     }
    2821       37170 :     set_avma(av);
    2822       37170 :     vP = idealprimedec_limit_norm(bnf, p, BOUND);
    2823       37170 :     J = lg(vP);
    2824             :     /* if last is unramified, all P|p in subgroup generated by FB: skip last */
    2825       37170 :     if (J > 1 && pr_get_e(gel(vP,J-1)) == 1) J--;
    2826       37170 :     if (J == 1) continue;
    2827       14448 :     if (DEBUGLEVEL>1) err_printf("*** p = %Ps\n",p);
    2828       14448 :     pr_orbit = auts? zero_zv(J-1): NULL;
    2829       31353 :     for (j = 1; j < J; j++)
    2830             :     {
    2831       16905 :       GEN P = gel(vP,j);
    2832       16905 :       long k = 0;
    2833       16905 :       if (pr_orbit)
    2834             :       {
    2835       11858 :         if (pr_orbit[j]) continue;
    2836             :         /* discard all primes in automorphism orbit simultaneously */
    2837       11851 :         pr_orbit_fill(pr_orbit, auts, vP, j);
    2838             :       }
    2839       16898 :       if (abscmpiu(p, pmax) > 0 || !(k = tablesearch(fb, P, &cmp_prime_ideal)))
    2840       16338 :         (void)SPLIT(&F, nf, pr_hnf(nf,P), Vbase, fact);
    2841       16898 :       if (DEBUGLEVEL>1)
    2842             :       {
    2843           0 :         err_printf("  Testing P = %Ps\n",P);
    2844           0 :         if (k) err_printf("    #%ld in factor base\n",k);
    2845           0 :         else err_printf("    is %Ps\n", isprincipal(bnf,P));
    2846             :       }
    2847             :     }
    2848             :   }
    2849          56 :   set_avma(av0);
    2850          56 : }
    2851             : 
    2852             : /* A t_MAT of complex floats, in fact reals. Extract a submatrix B
    2853             :  * whose columns are definitely non-0, i.e. gexpo(A[j]) >= -2
    2854             :  *
    2855             :  * If possible precision problem (t_REAL 0 with large exponent), set
    2856             :  * *precpb to 1 */
    2857             : static GEN
    2858       17920 : clean_cols(GEN A, int *precpb)
    2859             : {
    2860       17920 :   long l = lg(A), h, i, j, k;
    2861             :   GEN B;
    2862       17920 :   *precpb = 0;
    2863       17920 :   if (l == 1) return A;
    2864       17920 :   h = lgcols(A);;
    2865       17920 :   B = cgetg(l, t_MAT);
    2866     1786404 :   for (i = k = 1; i < l; i++)
    2867             :   {
    2868     1768484 :     GEN Ai = gel(A,i);
    2869     1768484 :     int non0 = 0;
    2870     8938552 :     for (j = 1; j < h; j++)
    2871             :     {
    2872     7170068 :       GEN c = gel(Ai,j);
    2873     7170068 :       if (gexpo(c) >= -2)
    2874             :       {
    2875     6558250 :         if (gequal0(c)) *precpb = 1; else non0 = 1;
    2876             :       }
    2877             :     }
    2878     1768484 :     if (non0) gel(B, k++) = Ai;
    2879             :   }
    2880       17920 :   setlg(B, k); return B;
    2881             : }
    2882             : 
    2883             : static long
    2884     1596866 : compute_multiple_of_R_pivot(GEN X, GEN x0/*unused*/, long ix, GEN c)
    2885             : {
    2886     1596866 :   GEN x = gel(X,ix);
    2887     1596866 :   long i, k = 0, ex = - (long)HIGHEXPOBIT, lx = lg(x);
    2888             :   (void)x0;
    2889     8232105 :   for (i=1; i<lx; i++)
    2890     6635239 :     if (!c[i] && !gequal0(gel(x,i)))
    2891             :     {
    2892     1203494 :       long e = gexpo(gel(x,i));
    2893     1203494 :       if (e > ex) { ex = e; k = i; }
    2894             :     }
    2895     1596866 :   return (k && ex > -32)? k: lx;
    2896             : }
    2897             : 
    2898             : /* Ar = (log |sigma_i(u_j)|) for units (u_j) found so far,
    2899             :  * RU = R1+R2 = unit rank, N = field degree
    2900             :  * need = unit rank defect
    2901             :  * L = NULL (prec problem) or B^(-1) * A with approximate rational entries
    2902             :  * (as t_REAL), B a submatrix of A, with (probably) maximal rank RU */
    2903             : static GEN
    2904       25943 : compute_multiple_of_R(GEN Ar, long RU, long N, long *pneed, long *bit, GEN *ptL)
    2905             : {
    2906             :   GEN T, d, mdet, Im_mdet, kR, L;
    2907       25943 :   long i, j, r, R1 = 2*RU - N;
    2908             :   int precpb;
    2909       25943 :   pari_sp av = avma;
    2910             : 
    2911       25943 :   if (RU == 1) { *ptL = zeromat(0, lg(Ar)-1); return gen_1; }
    2912             : 
    2913       17920 :   if (DEBUGLEVEL) err_printf("\n#### Computing regulator multiple\n");
    2914       17920 :   mdet = clean_cols(Ar, &precpb);
    2915             :   /* will cause precision to increase on later failure, but we may succeed! */
    2916       17920 :   *ptL = precpb? NULL: gen_1;
    2917       17920 :   T = cgetg(RU+1,t_COL);
    2918       55298 :   for (i=1; i<=R1; i++) gel(T,i) = gen_1;
    2919       37627 :   for (   ; i<=RU; i++) gel(T,i) = gen_2;
    2920       17920 :   mdet = shallowconcat(T, mdet); /* det(Span(mdet)) = N * R */
    2921             : 
    2922             :   /* could be using indexrank(), but need custom "get_pivot" function */
    2923       17920 :   d = RgM_pivots(mdet, NULL, &r, &compute_multiple_of_R_pivot);
    2924             :   /* # of independent columns == unit rank ? */
    2925       17920 :   if (lg(mdet)-1 - r != RU)
    2926             :   {
    2927       10308 :     if (DEBUGLEVEL)
    2928           0 :       err_printf("Unit group rank = %ld < %ld\n",lg(mdet)-1 - r, RU);
    2929       10308 :     *pneed = RU - (lg(mdet)-1-r); return gc_NULL(av);
    2930             :   }
    2931             : 
    2932        7612 :   Im_mdet = cgetg(RU+1, t_MAT); /* extract independent columns */
    2933             :   /* N.B: d[1] = 1, corresponding to T above */
    2934        7612 :   gel(Im_mdet, 1) = T;
    2935       49636 :   for (i = j = 2; i <= RU; j++)
    2936       42024 :     if (d[j]) gel(Im_mdet, i++) = gel(mdet,j);
    2937             : 
    2938             :   /* integral multiple of R: the cols we picked form a Q-basis, they have an
    2939             :    * index in the full lattice. First column is T */
    2940        7612 :   kR = divru(det2(Im_mdet), N);
    2941             :   /* R > 0.2 uniformly */
    2942        7612 :   if (!signe(kR) || expo(kR) < -3)
    2943             :   {
    2944           0 :     if (DEBUGLEVEL) err_printf("Regulator is zero.\n");
    2945           0 :     *pneed = 0; return gc_NULL(av);
    2946             :   }
    2947        7612 :   d = det2(rowslice(vecslice(Im_mdet, 2, RU), 2, RU));
    2948        7612 :   setabssign(d); setabssign(kR);
    2949        7612 :   if (gexpo(gsub(d,kR)) - gexpo(d) > -20) { *ptL = NULL; return gc_NULL(av); }
    2950        7578 :   L = RgM_inv(Im_mdet);
    2951             :   /* estimate # of correct bits in result */
    2952        7578 :   if (!L || (*bit = - gexpo(RgM_Rg_sub(RgM_mul(L,Im_mdet), gen_1))) < 16)
    2953           0 :   { *ptL = NULL; return gc_NULL(av); }
    2954             : 
    2955        7578 :   L = RgM_mul(rowslice(L,2,RU), Ar); /* approximate rational entries */
    2956        7578 :   gerepileall(av,2, &L, &kR);
    2957        7578 :   *ptL = L; return kR;
    2958             : }
    2959             : 
    2960             : /* leave small integer n as is, convert huge n to t_REAL (for readability) */
    2961             : static GEN
    2962           0 : i2print(GEN n)
    2963           0 : { return lgefint(n) <= DEFAULTPREC? n: itor(n,LOWDEFAULTPREC); }
    2964             : 
    2965             : static long
    2966       15584 : bad_check(GEN c)
    2967             : {
    2968       15584 :   long ec = gexpo(c);
    2969       15584 :   if (DEBUGLEVEL) err_printf("\n ***** check = %.28Pg\n",c);
    2970             :   /* safe check for c < 0.75 : avoid underflow in gtodouble() */
    2971       15584 :   if (ec < -1 || (ec == -1 && gtodouble(c) < 0.75)) return fupb_PRECI;
    2972             :   /* safe check for c > 1.3 : avoid overflow */
    2973       15584 :   if (ec > 0 || (ec == 0 && gtodouble(c) > 1.3)) return fupb_RELAT;
    2974       12062 :   return fupb_NONE;
    2975             : }
    2976             : /* Input:
    2977             :  * lambda = approximate rational entries: coords of units found so far on a
    2978             :  * sublattice of maximal rank (sublambda)
    2979             :  * *ptkR = regulator of sublambda = multiple of regulator of lambda
    2980             :  * Compute R = true regulator of lambda.
    2981             :  *
    2982             :  * If c := Rz ~ 1, by Dirichlet's formula, then lambda is the full group of
    2983             :  * units AND the full set of relations for the class group has been computed.
    2984             :  *
    2985             :  * In fact z is a very rough approximation and we only expect 0.75 < Rz < 1.3
    2986             :  * bit is an estimate for the actual accuracy of lambda
    2987             :  *
    2988             :  * Output: *ptkR = R, *ptU = basis of fundamental units (in terms lambda) */
    2989             : static long
    2990       15601 : compute_R(GEN lambda, GEN z, long bit, GEN *ptL, GEN *ptkR)
    2991             : {
    2992       15601 :   pari_sp av = avma;
    2993       15601 :   long r, reason, RU = lg(lambda) == 1? 1: lgcols(lambda);
    2994             :   GEN L, H, D, den, R, c;
    2995             : 
    2996       15601 :   *ptL = NULL;
    2997       15601 :   if (DEBUGLEVEL) err_printf("\n#### Computing check\n");
    2998       15601 :   if (RU == 1) { *ptkR = gen_1; *ptL = lambda; return bad_check(z); }
    2999        7578 :   D = gmul2n(mpmul(*ptkR,z), 1); /* bound for denom(lambda) */
    3000        7578 :   if (expo(D) < 0 && rtodbl(D) < 0.95) return fupb_PRECI;
    3001        7578 :   lambda = bestappr(lambda,D);
    3002        7578 :   if (lg(lambda) == 1)
    3003             :   {
    3004           0 :     if (DEBUGLEVEL) err_printf("truncation error in bestappr\n");
    3005           0 :     return fupb_PRECI;
    3006             :   }
    3007        7578 :   den = Q_denom(lambda);
    3008        7578 :   if (mpcmp(den,D) > 0)
    3009             :   {
    3010           7 :     if (DEBUGLEVEL) err_printf("D = %Ps\nden = %Ps\n",D, i2print(den));
    3011           7 :     return fupb_PRECI;
    3012             :   }
    3013        7571 :   L = Q_muli_to_int(lambda, den);
    3014        7571 :   if (bit > 0)
    3015             :   {
    3016        4260 :     if (lg(L) > 1)
    3017             :     {
    3018        4260 :       if (RU > 5) bit -= 64;
    3019        4103 :       else if (RU > 3) bit -= 32;
    3020             :     }
    3021        4260 :     if (gexpo(L) + expi(den) > bit - 32)
    3022             :     {
    3023          10 :       if (DEBUGLEVEL) err_printf("dubious bestappr; den = %Ps\n", i2print(den));
    3024          10 :       return fupb_PRECI;
    3025             :     }
    3026             :   }
    3027        7561 :   H = ZM_hnf(L); r = lg(H)-1;
    3028        7561 :   if (!r || r != nbrows(H))
    3029           0 :     R = gen_0; /* wrong rank */
    3030             :   else
    3031        7561 :     R = gmul(*ptkR, gdiv(ZM_det_triangular(H), powiu(den, r)));
    3032             :   /* R = tentative regulator; regulator > 0.2 uniformly */
    3033        7561 :   if (gexpo(R) < -3) {
    3034           0 :     if (DEBUGLEVEL) err_printf("\n#### Tentative regulator: %.28Pg\n", R);
    3035           0 :     return gc_long(av, fupb_PRECI);
    3036             :   }
    3037        7561 :   c = gmul(R,z); /* should be n (= 1 if we are done) */
    3038        7561 :   if (DEBUGLEVEL) err_printf("\n#### Tentative regulator: %.28Pg\n", R);
    3039        7561 :   if ((reason = bad_check(c))) return gc_long(av, reason);
    3040        4298 :   *ptkR = R; *ptL = L; return fupb_NONE;
    3041             : }
    3042             : static GEN
    3043       12131 : get_clg2(GEN cyc, GEN Ga, GEN C, GEN Ur, GEN Ge, GEN M1, GEN M2)
    3044             : {
    3045       12131 :   GEN GD = gsub(act_arch(M1, C), diagact_arch(cyc, Ga));
    3046       12131 :   GEN ga = gsub(act_arch(M2, C), act_arch(Ur, Ga));
    3047       12131 :   return mkvecn(6, Ur, ga, GD, Ge, M1, M2);
    3048             : }
    3049             : /* compute class group (clg1) + data for isprincipal (clg2) */
    3050             : static GEN
    3051       12040 : class_group_gen(GEN nf,GEN W,GEN C,GEN Vbase,long prec, GEN *pclg2)
    3052             : {
    3053             :   GEN M1, M2, z, G, Ga, Ge, cyc, X, Y, D, U, V, Ur, Ui, Uir;
    3054             :   long j, l;
    3055             : 
    3056       12040 :   D = ZM_snfall(W,&U,&V); /* UWV=D, D diagonal, G = g Ui (G=new gens, g=old) */
    3057       12040 :   Ui = ZM_inv(U, NULL);
    3058       12040 :   l = lg(D); cyc = cgetg(l, t_VEC); /* elementary divisors */
    3059       21903 :   for (j = 1; j < l; j++)
    3060             :   {
    3061       10725 :     gel(cyc,j) = gcoeff(D,j,j); /* strip useless components */
    3062       10725 :     if (is_pm1(gel(cyc,j))) break;
    3063             :   }
    3064       12040 :   l = j;
    3065       12040 :   Ur  = ZM_hnfdivrem(U, D, &Y);
    3066       12040 :   Uir = ZM_hnfdivrem(Ui,W, &X);
    3067             :  /* {x} = logarithmic embedding of x (arch. component)
    3068             :   * NB: [J,z] = idealred(I) --> I = y J, with {y} = - z
    3069             :   * G = g Uir - {Ga},  Uir = Ui + WX
    3070             :   * g = G Ur  - {ga},  Ur  = U + DY */
    3071       12040 :   G = cgetg(l,t_VEC);
    3072       12040 :   Ga= cgetg(l,t_MAT);
    3073       12040 :   Ge= cgetg(l,t_COL);
    3074       12040 :   z = init_famat(NULL);
    3075       21903 :   for (j = 1; j < l; j++)
    3076             :   {
    3077        9863 :     GEN I = genback(z, nf, Vbase, gel(Uir,j));
    3078        9863 :     gel(G,j) = gel(I,1); /* generator, order cyc[j] */
    3079        9863 :     gel(Ge,j)= gel(I,2);
    3080        9863 :     gel(Ga,j)= nf_cxlog(nf, gel(I,2), prec);
    3081        9863 :     if (!gel(Ga,j)) pari_err_PREC("class_group_gen");
    3082             :   }
    3083             :   /* {ga} = {GD}Y + G U - g = {GD}Y - {Ga} U + gW X U
    3084             :                             = gW (X Ur + V Y) - {Ga}Ur */
    3085       12040 :   M2 = ZM_add(ZM_mul(X,Ur), ZM_mul(V,Y));
    3086       12040 :   setlg(cyc,l); setlg(V,l); setlg(D,l);
    3087             :   /* G D =: {GD} = g (Ui + W X) D - {Ga}D = g W (V + X D) - {Ga}D
    3088             :    * NB: Ui D = W V. gW is given by (first l-1 cols of) C */
    3089       12040 :   M1 = ZM_add(V, ZM_mul(X,D));
    3090       12040 :   *pclg2 = get_clg2(cyc, Ga, C, Ur, Ge, M1, M2);
    3091       12040 :   return mkvec3(ZV_prod(cyc), cyc, G);
    3092             : }
    3093             : 
    3094             : /* compute principal ideals corresponding to (gen[i]^cyc[i]) */
    3095             : static GEN
    3096        2471 : makecycgen(GEN bnf)
    3097             : {
    3098             :   GEN cyc, gen, h, nf, y, GD;
    3099             :   long e,i,l;
    3100             : 
    3101        2471 :   if (DEBUGLEVEL) pari_warn(warner,"completing bnf (building cycgen)");
    3102        2471 :   nf = bnf_get_nf(bnf);
    3103        2471 :   cyc = bnf_get_cyc(bnf);
    3104        2471 :   gen = bnf_get_gen(bnf);
    3105        2471 :   GD = bnf_get_GD(bnf);
    3106        2471 :   h = cgetg_copy(gen, &l);
    3107        5383 :   for (i = 1; i < l; i++)
    3108             :   {
    3109        2912 :     GEN gi = gel(gen,i), ci = gel(cyc,i);
    3110        2912 :     if (abscmpiu(ci, 5) < 0)
    3111             :     {
    3112        2135 :       GEN N = ZM_det_triangular(gi);
    3113        2135 :       y = isprincipalarch(bnf,gel(GD,i), N, ci, gen_1, &e);
    3114        2135 :       if (y && fact_ok(nf,y,NULL,mkvec(gi),mkvec(ci)))
    3115             :       {
    3116        2079 :         gel(h,i) = to_famat_shallow(y,gen_1);
    3117        2079 :         continue;
    3118             :       }
    3119             :     }
    3120         833 :     y = isprincipalfact(bnf, NULL, mkvec(gi), mkvec(ci), nf_GENMAT|nf_FORCE);
    3121         833 :     gel(h,i) = gel(y,2);
    3122             :   }
    3123        2471 :   return h;
    3124             : }
    3125             : 
    3126             : static GEN
    3127          21 : get_y(GEN bnf, GEN W, GEN B, GEN C, GEN pFB, long j)
    3128             : {
    3129          21 :   GEN y, nf  = bnf_get_nf(bnf);
    3130          21 :   long e, lW = lg(W)-1;
    3131          21 :   GEN ex = (j<=lW)? gel(W,j): gel(B,j-lW);
    3132          21 :   GEN P = (j<=lW)? NULL: gel(pFB,j);
    3133          21 :   if (C)
    3134             :   { /* archimedean embeddings known: cheap trial */
    3135          21 :     GEN Nx = get_norm_fact_primes(pFB, ex, P);
    3136          21 :     y = isprincipalarch(bnf,gel(C,j), Nx,gen_1, gen_1, &e);
    3137          21 :     if (y && fact_ok(nf,y,P,pFB,ex)) return y;
    3138             :   }
    3139           0 :   y = isprincipalfact_or_fail(bnf, P, pFB, ex);
    3140           0 :   return typ(y) == t_INT? y: gel(y,2);
    3141             : }
    3142             : /* compute principal ideals corresponding to bnf relations */
    3143             : static GEN
    3144          14 : makematal(GEN bnf)
    3145             : {
    3146          14 :   GEN W = bnf_get_W(bnf), B = bnf_get_B(bnf), C = bnf_get_C(bnf);
    3147             :   GEN pFB, ma, retry;
    3148          14 :   long lma, j, prec = 0;
    3149             : 
    3150          14 :   if (DEBUGLEVEL) pari_warn(warner,"completing bnf (building matal)");
    3151          14 :   lma=lg(W)+lg(B)-1;
    3152          14 :   pFB = bnf_get_vbase(bnf);
    3153          14 :   ma = cgetg(lma,t_VEC);
    3154          14 :   retry = vecsmalltrunc_init(lma);
    3155          35 :   for (j=lma-1; j>0; j--)
    3156             :   {
    3157          21 :     pari_sp av = avma;
    3158          21 :     GEN y = get_y(bnf, W, B, C, pFB, j);
    3159          21 :     if (typ(y) == t_INT)
    3160             :     {
    3161           0 :       long E = itos(y);
    3162           0 :       if (DEBUGLEVEL>1) err_printf("\n%ld done later at prec %ld\n",j,E);
    3163           0 :       set_avma(av);
    3164           0 :       vecsmalltrunc_append(retry, j);
    3165           0 :       if (E > prec) prec = E;
    3166             :     }
    3167             :     else
    3168             :     {
    3169          21 :       if (DEBUGLEVEL>1) err_printf("%ld ",j);
    3170          21 :       gel(ma,j) = gerepileupto(av,y);
    3171             :     }
    3172             :   }
    3173          14 :   if (prec)
    3174             :   {
    3175           0 :     long k, l = lg(retry);
    3176           0 :     GEN y, nf = bnf_get_nf(bnf);
    3177           0 :     if (DEBUGLEVEL) pari_warn(warnprec,"makematal",prec);
    3178           0 :     nf = nfnewprec_shallow(nf,prec);
    3179           0 :     bnf = Buchall(nf, nf_FORCE, prec);
    3180           0 :     if (DEBUGLEVEL) err_printf("makematal, adding missing entries:");
    3181           0 :     for (k=1; k<l; k++)
    3182             :     {
    3183           0 :       pari_sp av = avma;
    3184           0 :       long j = retry[k];
    3185           0 :       y = get_y(bnf,W,B,NULL, pFB, j);
    3186           0 :       if (typ(y) == t_INT) pari_err_PREC("makematal");
    3187           0 :       if (DEBUGLEVEL>1) err_printf("%ld ",j);
    3188           0 :       gel(ma,j) = gerepileupto(av,y);
    3189             :     }
    3190             :   }
    3191          14 :   if (DEBUGLEVEL>1) err_printf("\n");
    3192          14 :   return ma;
    3193             : }
    3194             : 
    3195             : enum { MATAL = 1, CYCGEN, UNITS };
    3196             : GEN
    3197       13321 : bnf_build_cycgen(GEN bnf)
    3198       13321 : { return obj_checkbuild(bnf, CYCGEN, &makecycgen); }
    3199             : GEN
    3200          14 : bnf_build_matalpha(GEN bnf)
    3201          14 : { return obj_checkbuild(bnf, MATAL, &makematal); }
    3202             : GEN
    3203        7729 : bnf_build_units(GEN bnf)
    3204        7729 : { return obj_checkbuild(bnf, UNITS, &makeunits); }
    3205             : 
    3206             : /* return fu in compact form if available; in terms of a fixed basis
    3207             :  * of S-units */
    3208             : GEN
    3209          42 : bnf_compactfu_mat(GEN bnf)
    3210             : {
    3211          42 :   GEN X, U, SUnits = bnf_get_sunits(bnf);
    3212          42 :   if (!SUnits) return NULL;
    3213          42 :   X = gel(SUnits,1);
    3214          42 :   U = gel(SUnits,2); ZM_remove_unused(&U, &X);
    3215          42 :   return mkvec2(X, U);
    3216             : }
    3217             : /* return fu in compact form if available; individually as famat */
    3218             : GEN
    3219         819 : bnf_compactfu(GEN bnf)
    3220             : {
    3221         819 :   GEN fu, X, U, SUnits = bnf_get_sunits(bnf);
    3222             :   long i, l;
    3223         819 :   if (!SUnits) return NULL;
    3224         805 :   X = gel(SUnits,1);
    3225         805 :   U = gel(SUnits,2); l = lg(U); fu = cgetg(l, t_VEC);
    3226        2492 :   for (i = 1; i < l; i++)
    3227        1687 :     gel(fu,i) = famat_remove_trivial(mkmat2(X, gel(U,i)));
    3228         805 :   return fu;
    3229             : }
    3230             : /* return expanded fu if available */
    3231             : GEN
    3232       38955 : bnf_has_fu(GEN bnf)
    3233             : {
    3234       38955 :   GEN fu = obj_check(bnf, UNITS);
    3235       38955 :   if (fu) return vecsplice(fu, 1);
    3236       38114 :   fu = bnf_get_fu_nocheck(bnf);
    3237       38114 :   return (typ(fu) == t_MAT)? NULL: fu;
    3238             : }
    3239             : /* return expanded fu if available; build if cheap */
    3240             : GEN
    3241       38892 : bnf_build_cheapfu(GEN bnf)
    3242             : {
    3243             :   GEN fu, SUnits;
    3244       38892 :   if ((fu = bnf_has_fu(bnf))) return fu;
    3245         120 :   if ((SUnits = bnf_get_sunits(bnf)))
    3246             :   {
    3247         120 :     pari_sp av = avma;
    3248         120 :     long e = gexpo(real_i(bnf_get_logfu(bnf)));
    3249         120 :     set_avma(av); if (e < 13) return vecsplice(bnf_build_units(bnf), 1);
    3250             :   }
    3251          28 :   return NULL;
    3252             : }
    3253             : 
    3254             : static GEN
    3255          91 : get_regulator(GEN mun)
    3256             : {
    3257          91 :   pari_sp av = avma;
    3258             :   GEN R;
    3259             : 
    3260          91 :   if (lg(mun) == 1) return gen_1;
    3261          84 :   R = det( rowslice(real_i(mun), 1, lgcols(mun)-2) );
    3262          84 :   setabssign(R); return gerepileuptoleaf(av, R);
    3263             : }
    3264             : 
    3265             : /* return corrected archimedian components for elts of x (vector)
    3266             :  * (= log(sigma_i(x)) - log(|Nx|) / [K:Q]) */
    3267             : static GEN
    3268          28 : get_archclean(GEN nf, GEN x, long prec, int units)
    3269             : {
    3270          28 :   long k, N, l = lg(x);
    3271          28 :   GEN M = cgetg(l, t_MAT);
    3272             : 
    3273          28 :   if (l == 1) return M;
    3274          14 :   N = nf_get_degree(nf);
    3275          42 :   for (k = 1; k < l; k++)
    3276             :   {
    3277          28 :     pari_sp av = avma;
    3278          28 :     GEN c = nf_cxlog(nf, gel(x,k), prec);
    3279          28 :     if (!c || (!units && !(c = cleanarch(c, N, prec)))) return NULL;
    3280          28 :     gel(M,k) = gerepilecopy(av, c);
    3281             :   }
    3282          14 :   return M;
    3283             : }
    3284             : static void
    3285          77 : Sunits_archclean(GEN nf, GEN Sunits, GEN *pmun, GEN *pC, long prec)
    3286             : {
    3287          77 :   GEN M, X = gel(Sunits,1), U = gel(Sunits,2), G = gel(Sunits,3);
    3288          77 :   long k, N = nf_get_degree(nf), l = lg(X);
    3289             : 
    3290          77 :   M = cgetg(l, t_MAT);
    3291        3290 :   for (k = 1; k < l; k++)
    3292        3213 :     if (!(gel(M,k) = nf_cxlog(nf, gel(X,k), prec))) return;
    3293          77 :   *pmun = cleanarch(RgM_ZM_mul(M, U), N, prec);
    3294          77 :   if (*pmun) *pC = cleanarch(RgM_ZM_mul(M, G), N, prec);
    3295             : }
    3296             : 
    3297             : GEN
    3298          91 : bnfnewprec_shallow(GEN bnf, long prec)
    3299             : {
    3300          91 :   GEN nf0 = bnf_get_nf(bnf), nf, v, fu, matal, y, mun, C;
    3301          91 :   GEN Sunits = bnf_get_sunits(bnf), Ur, Ga, Ge, M1, M2;
    3302          91 :   long r1, r2, prec0 = prec;
    3303             : 
    3304          91 :   nf_get_sign(nf0, &r1, &r2);
    3305          91 :   if (Sunits)
    3306             :   {
    3307          77 :     fu = matal = NULL;
    3308          77 :     prec += nbits2extraprec(gexpo(Sunits));
    3309             :   }
    3310             :   else
    3311             :   {
    3312          14 :     fu = bnf_build_units(bnf);
    3313          14 :     fu = vecslice(fu, 2, lg(fu)-1);
    3314          14 :     if (r1 + r2 > 1) {
    3315           7 :       long e = gexpo(bnf_get_logfu(bnf)) + 1 - TWOPOTBITS_IN_LONG;
    3316           7 :       if (e >= 0) prec += nbits2extraprec(e);
    3317             :     }
    3318          14 :     matal = bnf_build_matalpha(bnf);
    3319             :   }
    3320             : 
    3321          91 :   if (DEBUGLEVEL && prec0 != prec) pari_warn(warnprec,"bnfnewprec",prec);
    3322          91 :   for(C = NULL;;)
    3323           0 :   {
    3324          91 :     pari_sp av = avma;
    3325          91 :     nf = nfnewprec_shallow(nf0,prec);
    3326          91 :     if (Sunits)
    3327          77 :       Sunits_archclean(nf, Sunits, &mun, &C, prec);
    3328             :     else
    3329             :     {
    3330          14 :       mun = get_archclean(nf, fu, prec, 1);
    3331          14 :       if (mun) C = get_archclean(nf, matal, prec, 0);
    3332             :     }
    3333          91 :     if (C) break;
    3334           0 :     set_avma(av); prec = precdbl(prec);
    3335           0 :     if (DEBUGLEVEL) pari_warn(warnprec,"bnfnewprec(extra)",prec);
    3336             :   }
    3337          91 :   y = leafcopy(bnf);
    3338          91 :   gel(y,3) = mun;
    3339          91 :   gel(y,4) = C;
    3340          91 :   gel(y,7) = nf;
    3341          91 :   gel(y,8) = v = leafcopy(gel(bnf,8));
    3342          91 :   gel(v,2) = get_regulator(mun);
    3343          91 :   v = gel(bnf,9);
    3344          91 :   if (lg(v) < 7) pari_err_TYPE("bnfnewprec [obsolete bnf format]", bnf);
    3345          91 :   Ur = gel(v,1);
    3346          91 :   Ge = gel(v,4);
    3347          91 :   Ga = nfV_cxlog(nf, Ge, prec);
    3348          91 :   M1 = gel(v,5);
    3349          91 :   M2 = gel(v,6);
    3350          91 :   gel(y,9) = get_clg2(bnf_get_cyc(bnf), Ga, C, Ur, Ge, M1, M2);
    3351          91 :   return y;
    3352             : }
    3353             : GEN
    3354          21 : bnfnewprec(GEN bnf, long prec)
    3355             : {
    3356          21 :   pari_sp av = avma;
    3357          21 :   return gerepilecopy(av, bnfnewprec_shallow(checkbnf(bnf), prec));
    3358             : }
    3359             : 
    3360             : GEN
    3361           0 : bnrnewprec_shallow(GEN bnr, long prec)
    3362             : {
    3363           0 :   GEN y = cgetg(7,t_VEC);
    3364             :   long i;
    3365           0 :   gel(y,1) = bnfnewprec_shallow(bnr_get_bnf(bnr), prec);
    3366           0 :   for (i=2; i<7; i++) gel(y,i) = gel(bnr,i);
    3367           0 :   return y;
    3368             : }
    3369             : GEN
    3370           7 : bnrnewprec(GEN bnr, long prec)
    3371             : {
    3372           7 :   GEN y = cgetg(7,t_VEC);
    3373             :   long i;
    3374           7 :   checkbnr(bnr);
    3375           7 :   gel(y,1) = bnfnewprec(bnr_get_bnf(bnr), prec);
    3376          42 :   for (i=2; i<7; i++) gel(y,i) = gcopy(gel(bnr,i));
    3377           7 :   return y;
    3378             : }
    3379             : 
    3380             : static GEN
    3381       12817 : buchall_end(GEN nf,GEN res, GEN clg2, GEN W, GEN B, GEN A, GEN C,GEN Vbase)
    3382             : {
    3383       12817 :   GEN z = obj_init(9, 3);
    3384       12817 :   gel(z,1) = W;
    3385       12817 :   gel(z,2) = B;
    3386       12817 :   gel(z,3) = A;
    3387       12817 :   gel(z,4) = C;
    3388       12817 :   gel(z,5) = Vbase;
    3389       12817 :   gel(z,6) = gen_0;
    3390       12817 :   gel(z,7) = nf;
    3391       12817 :   gel(z,8) = res;
    3392       12817 :   gel(z,9) = clg2;
    3393       12817 :   return z;
    3394             : }
    3395             : 
    3396             : GEN
    3397        1792 : bnfinit0(GEN P, long flag, GEN data, long prec)
    3398             : {
    3399        1792 :   double c1 = 0., c2 = 0.;
    3400        1792 :   long fl, relpid = BNF_RELPID;
    3401             : 
    3402        1792 :   if (data)
    3403             :   {
    3404          35 :     long lx = lg(data);
    3405          35 :     if (typ(data) != t_VEC || lx > 5) pari_err_TYPE("bnfinit",data);
    3406          35 :     switch(lx)
    3407             :     {
    3408           0 :       case 4: relpid = itos(gel(data,3));
    3409          21 :       case 3: c2 = gtodouble(gel(data,2));
    3410          28 :       case 2: c1 = gtodouble(gel(data,1));
    3411             :     }
    3412             :   }
    3413        1792 :   switch(flag)
    3414             :   {
    3415        1470 :     case 2:
    3416        1470 :     case 0: fl = 0; break;
    3417         322 :     case 1: fl = nf_FORCE; break;
    3418           0 :     default: pari_err_FLAG("bnfinit");
    3419             :       return NULL; /* LCOV_EXCL_LINE */
    3420             :   }
    3421        1792 :   return Buchall_param(P, c1, c2, relpid, fl, prec);
    3422             : }
    3423             : GEN
    3424       11032 : Buchall(GEN P, long flag, long prec)
    3425       11032 : { return Buchall_param(P, 0., 0., BNF_RELPID, flag & nf_FORCE, prec); }
    3426             : 
    3427             : static GEN
    3428         777 : Buchall_deg1(GEN nf)
    3429             : {
    3430         777 :   GEN v = cgetg(1,t_VEC), m = cgetg(1,t_MAT);
    3431         777 :   GEN res, W, A, B, C, Vbase = cgetg(1,t_COL);
    3432         777 :   GEN fu = v, R = gen_1, zu = mkvec2(gen_2, gen_m1);
    3433         777 :   GEN clg1 = mkvec3(gen_1,v,v), clg2 = mkvecn(6, m,m,m,v,m,m);
    3434             : 
    3435         777 :   W = A = B = C = m; res = mkvec5(clg1, R, gen_1, zu, fu);
    3436         777 :   return buchall_end(nf,res,clg2,W,B,A,C,Vbase);
    3437             : }
    3438             : 
    3439             : /* return (small set of) indices of columns generating the same lattice as x.
    3440             :  * Assume HNF(x) is inexpensive (few rows, many columns).
    3441             :  * Dichotomy approach since interesting columns may be at the very end */
    3442             : GEN
    3443       12041 : extract_full_lattice(GEN x)
    3444             : {
    3445       12041 :   long dj, j, k, l = lg(x);
    3446             :   GEN h, h2, H, v;
    3447             : 
    3448       12041 :   if (l < 200) return NULL; /* not worth it */
    3449             : 
    3450           7 :   v = vecsmalltrunc_init(l);
    3451           7 :   H = ZM_hnf(x);
    3452           7 :   h = cgetg(1, t_MAT);
    3453           7 :   dj = 1;
    3454         315 :   for (j = 1; j < l; )
    3455             :   {
    3456         315 :     pari_sp av = avma;
    3457         315 :     long lv = lg(v);
    3458             : 
    3459        3192 :     for (k = 0; k < dj; k++) v[lv+k] = j+k;
    3460         315 :     setlg(v, lv + dj);
    3461         315 :     h2 = ZM_hnf(vecpermute(x, v));
    3462         315 :     if (ZM_equal(h, h2))
    3463             :     { /* these dj columns can be eliminated */
    3464         161 :       set_avma(av); setlg(v, lv);
    3465         161 :       j += dj;
    3466         161 :       if (j >= l) break;
    3467         161 :       dj <<= 1;
    3468         161 :       if (j + dj >= l) { dj = (l - j) >> 1; if (!dj) dj = 1; }
    3469             :     }
    3470         154 :     else if (dj > 1)
    3471             :     { /* at least one interesting column, try with first half of this set */
    3472         112 :       set_avma(av); setlg(v, lv);
    3473         112 :       dj >>= 1; /* > 0 */
    3474             :     }
    3475             :     else
    3476             :     { /* this column should be kept */
    3477          42 :       if (ZM_equal(h2, H)) break;
    3478          35 :       h = h2; j++;
    3479             :     }
    3480             :   }
    3481           7 :   return v;
    3482             : }
    3483             : 
    3484             : static void
    3485       12082 : init_rel(RELCACHE_t *cache, FB_t *F, long add_need)
    3486             : {
    3487       12082 :   const long n = F->KC + add_need; /* expected # of needed relations */
    3488             :   long i, j, k, p;
    3489             :   GEN c, P;
    3490             :   GEN R;
    3491             : 
    3492       12082 :   if (DEBUGLEVEL) err_printf("KCZ = %ld, KC = %ld, n = %ld\n", F->KCZ,F->KC,n);
    3493       12082 :   reallocate(cache, 10*n + 50); /* make room for lots of relations */
    3494       12082 :   cache->chk = cache->base;
    3495       12082 :   cache->end = cache->base + n;
    3496       12082 :   cache->relsup = add_need;
    3497       12082 :   cache->last = cache->base;
    3498       12082 :   cache->missing = lg(cache->basis) - 1;
    3499       69958 :   for (i = 1; i <= F->KCZ; i++)
    3500             :   { /* trivial relations (p) = prod P^e */
    3501       57876 :     p = F->FB[i]; P = F->LV[p];
    3502       57876 :     if (!isclone(P)) continue;
    3503             : 
    3504             :     /* all prime divisors in FB */
    3505       45962 :     c = zero_Flv(F->KC); k = F->iLP[p];
    3506       45962 :     R = c; c += k;
    3507      147791 :     for (j = lg(P)-1; j; j--) c[j] = pr_get_e(gel(P,j));
    3508       45962 :     add_rel(cache, F, R, k+1, pr_get_p(gel(P,1)), 0);
    3509             :   }
    3510       12082 : }
    3511             : 
    3512             : /* Let z = \zeta_n in nf. List of not-obviously-dependent generators for
    3513             :  * cyclotomic units modulo torsion in Q(z) [independent when n a prime power]:
    3514             :  * - z^a - 1,  n/(a,n) not a prime power, a \nmid n unless a=1,  1 <= a < n/2
    3515             :  * - (Z^a - 1)/(Z - 1),  p^k || n, Z = z^{n/p^k}, (p,a) = 1, 1 < a <= (p^k-1)/2
    3516             :  */
    3517             : GEN
    3518       12082 : nfcyclotomicunits(GEN nf, GEN zu)
    3519             : {
    3520       12082 :   long n = itos(gel(zu, 1)), n2, lP, i, a;
    3521             :   GEN z, fa, P, E, L, mz, powz;
    3522       12082 :   if (n <= 6) return cgetg(1, t_VEC);
    3523             : 
    3524         182 :   z = algtobasis(nf,gel(zu, 2));
    3525         182 :   if ((n & 3) == 2) { n = n >> 1; z = ZC_neg(z); } /* ensure n != 2 (mod 4) */
    3526         182 :   n2 = n/2;
    3527         182 :   mz = zk_multable(nf, z); /* multiplication by z */
    3528         182 :   powz = cgetg(n2, t_VEC); gel(powz,1) = z;
    3529         322 :   for (i = 2; i < n2; i++) gel(powz,i) = ZM_ZC_mul(mz, gel(powz,i-1));
    3530             :   /* powz[i] = z^i */
    3531             : 
    3532         182 :   L = vectrunc_init(n);
    3533         182 :   fa = factoru(n);
    3534         182 :   P = gel(fa,1); lP = lg(P);
    3535         182 :   E = gel(fa,2);
    3536         385 :   for (i = 1; i < lP; i++)
    3537             :   { /* second kind */
    3538         203 :     long p = P[i], k = E[i], pk = upowuu(p,k), pk2 = (pk-1) / 2;
    3539         203 :     GEN u = gen_1;
    3540         399 :     for (a = 2; a <= pk2; a++)
    3541             :     {
    3542         196 :       u = nfadd(nf, u, gel(powz, (n/pk) * (a-1))); /* = (Z^a-1)/(Z-1) */
    3543         196 :       if (a % p) vectrunc_append(L, u);
    3544             :     }
    3545             :   }
    3546         287 :   if (lP > 2) for (a = 1; a < n2; a++)
    3547             :   { /* first kind, when n not a prime power */
    3548             :     ulong p;
    3549         105 :     if (a > 1 && (n % a == 0 || uisprimepower(n/ugcd(a,n), &p))) continue;
    3550          42 :     vectrunc_append(L, nfadd(nf, gel(powz, a), gen_m1));
    3551             :   }
    3552         182 :   return L;
    3553             : }
    3554             : static void
    3555       12082 : add_cyclotomic_units(GEN nf, GEN zu, RELCACHE_t *cache, FB_t *F)
    3556             : {
    3557       12082 :   pari_sp av = avma;
    3558       12082 :   GEN L = nfcyclotomicunits(nf, zu);
    3559       12082 :   long i, l = lg(L);
    3560       12082 :   if (l > 1)
    3561             :   {
    3562         182 :     GEN R = zero_Flv(F->KC);
    3563         399 :     for(i = 1; i < l; i++) add_rel(cache, F, R, F->KC+1, gel(L,i), 0);
    3564             :   }
    3565       12082 :   set_avma(av);
    3566       12082 : }
    3567             : 
    3568             : static GEN
    3569       78163 : trim_list(FB_t *F)
    3570             : {
    3571       78163 :   pari_sp av = avma;
    3572       78163 :   GEN v, L_jid = F->L_jid, minidx = F->minidx, present = zero_Flv(F->KC);
    3573       78163 :   long i, j, imax = minss(lg(L_jid), F->KC + 1);
    3574             : 
    3575       78163 :   v = cgetg(imax, t_VECSMALL);
    3576     3256924 :   for (i = j = 1; i < imax; i++)
    3577             :   {
    3578     3178761 :     long k = minidx[ L_jid[i] ];
    3579     3178761 :     if (!present[k]) { v[j++] = L_jid[i]; present[k] = 1; }
    3580             :   }
    3581       78163 :   setlg(v, j); return gerepileuptoleaf(av, v);
    3582             : }
    3583             : 
    3584             : static void
    3585       10957 : try_elt(RELCACHE_t *cache, FB_t *F, GEN nf, GEN x, FACT *fact)
    3586             : {
    3587       10957 :   pari_sp av = avma;
    3588             :   GEN R, Nx;
    3589       10957 :   long nz, tx = typ(x);
    3590             : 
    3591       10957 :   if (tx == t_INT || tx == t_FRAC) return;
    3592        7142 :   if (tx != t_COL) x = algtobasis(nf, x);
    3593        7142 :   if (RgV_isscalar(x)) return;
    3594        7142 :   x = Q_primpart(x);
    3595        7142 :   Nx = nfnorm(nf, x);
    3596        7142 :   if (!can_factor(F, nf, NULL, x, Nx, fact)) return;
    3597             : 
    3598             :   /* smooth element */
    3599        7142 :   R = set_fact(F, fact, NULL, &nz);
    3600             :   /* make sure we get maximal rank first, then allow all relations */
    3601        7142 :   (void) add_rel(cache, F, R, nz, x, 0);
    3602        7142 :   set_avma(av);
    3603             : }
    3604             : 
    3605             : static long
    3606     1073705 : scalar_bit(GEN x) { return bit_accuracy(gprecision(x)) - gexpo(x); }
    3607             : static long
    3608       11997 : RgM_bit(GEN x, long bit)
    3609             : {
    3610       11997 :   long i, j, m, b = bit, l = lg(x);
    3611       11997 :   if (l == 1) return b;
    3612       11997 :   m = lgcols(x);
    3613      472886 :   for (j = 1; j < l; j++)
    3614     1534594 :     for (i = 1; i < m; i++ ) b = minss(b, scalar_bit(gcoeff(x,i,j)));
    3615       11997 :   return b;
    3616             : }
    3617             : static void
    3618       10286 : matenlarge(GEN C, long h)
    3619             : {
    3620       10286 :   GEN _0 = zerocol(h);
    3621             :   long i;
    3622     1663775 :   for (i = lg(C); --i; ) gel(C,i) = shallowconcat(gel(C,i), _0);
    3623       10286 : }
    3624             : 
    3625             : /* E = floating point embeddings */
    3626             : static GEN
    3627       10286 : matbotidembs(RELCACHE_t *cache, GEN E)
    3628             : {
    3629       10286 :   long w = cache->last - cache->chk, h = cache->last - cache->base;
    3630       10286 :   long j, d = h - w, hE = nbrows(E);
    3631       10286 :   GEN y = cgetg(w+1,t_MAT), _0 = zerocol(h);
    3632       56396 :   for (j = 1; j <= w; j++)
    3633             :   {
    3634       46110 :     GEN c = shallowconcat(gel(E,j), _0);
    3635       46110 :     if (d + j >= 1) gel(c, d + j + hE) = gen_1;
    3636       46110 :     gel(y,j) = c;
    3637             :   }
    3638       10286 :   return y;
    3639             : }
    3640             : static GEN
    3641        1645 : matbotid(RELCACHE_t *cache)
    3642             : {
    3643        1645 :   long w = cache->last - cache->chk, h = cache->last - cache->base;
    3644        1645 :   long j, d = h - w;
    3645        1645 :   GEN y = cgetg(w+1,t_MAT);
    3646       40593 :   for (j = 1; j <= w; j++)
    3647             :   {
    3648       38948 :     GEN c = zerocol(h);
    3649       38948 :     if (d + j >= 1) gel(c, d + j) = gen_1;
    3650       38948 :     gel(y,j) = c;
    3651             :   }
    3652        1645 :   return y;
    3653             : }
    3654             : 
    3655             : static long
    3656          98 : myprecdbl(long prec, GEN C)
    3657             : {
    3658          98 :   long p = prec2nbits(prec) < 1280? precdbl(prec): (long)(prec * 1.5);
    3659          98 :   if (C) p = maxss(p, minss(3*p, prec + nbits2extraprec(gexpo(C))));
    3660          98 :   return p;
    3661             : }
    3662             : 
    3663             : static GEN
    3664        1366 : _nfnewprec(GEN nf, long prec, long *isclone)
    3665             : {
    3666        1366 :   GEN NF = gclone(nfnewprec_shallow(nf, prec));
    3667        1366 :   if (*isclone) gunclone(nf);
    3668        1366 :   *isclone = 1; return NF;
    3669             : }
    3670             : 
    3671             : /* Nrelid = nb relations per ideal, possibly 0. If flag is set, keep data in
    3672             :  * algebraic form. */
    3673             : GEN
    3674       12824 : Buchall_param(GEN P, double cbach, double cbach2, long Nrelid, long flag, long prec)
    3675             : {
    3676             :   pari_timer T;
    3677       12824 :   pari_sp av0 = avma, av, av2;
    3678             :   long PREC, N, R1, R2, RU, low, high, LIMC0, LIMC, LIMC2, LIMCMAX, zc, i;
    3679       12824 :   long LIMres, bit = 0, flag_nfinit = 0;
    3680       12824 :   long nreldep, sfb_trials, need, old_need, precdouble = 0, TRIES = 0;
    3681       12824 :   long nfisclone = 0;
    3682             :   long done_small, small_fail, fail_limit, squash_index, small_norm_prec;
    3683             :   double LOGD, LOGD2, lim;
    3684       12824 :   GEN computed = NULL, fu = NULL, zu, nf, M_sn, D, A, W, R, h, Ce, PERM;
    3685             :   GEN small_multiplier, auts, cyclic, embs, SUnits;
    3686             :   GEN res, L, invhr, B, C, lambda, dep, clg1, clg2, Vbase;
    3687       12824 :   const char *precpb = NULL;
    3688             :   nfmaxord_t nfT;
    3689             :   RELCACHE_t cache;
    3690             :   FB_t F;
    3691             :   GRHcheck_t GRHcheck;
    3692             :   FACT *fact;
    3693             : 
    3694       12824 :   if (DEBUGLEVEL) timer_start(&T);
    3695       12824 :   P = get_nfpol(P, &nf);
    3696       12817 :   if (nf)
    3697         588 :     D = nf_get_disc(nf);
    3698             :   else
    3699             :   {
    3700       12229 :     nfinit_basic(&nfT, P);
    3701       12229 :     D = nfT.dK;
    3702       12229 :     if (!ZX_is_monic(nfT.T0))
    3703             :     {
    3704          14 :       pari_warn(warner,"non-monic polynomial in bnfinit, using polredbest");
    3705          14 :       flag_nfinit = nf_RED;
    3706             :     }
    3707             :   }
    3708       12817 :   N = degpol(P);
    3709       12817 :   if (N <= 1)
    3710             :   {
    3711         777 :     if (!nf) nf = nfinit_complete(&nfT, flag_nfinit, DEFAULTPREC);
    3712         777 :     return gerepilecopy(av0, Buchall_deg1(nf));
    3713             :   }
    3714       12040 :   D = absi_shallow(D);
    3715       12040 :   LOGD = dbllog2(D) * M_LN2;
    3716       12040 :   LOGD2 = LOGD*LOGD;
    3717       12040 :   LIMCMAX = (long)(12.*LOGD2);
    3718             :   /* In small_norm, LLL reduction produces v0 in I such that
    3719             :    *     T2(v0) <= (4/3)^((n-1)/2) NI^(2/n) disc(K)^(1/n)
    3720             :    * We consider v with T2(v) <= BMULT * T2(v0)
    3721             :    * Hence Nv <= ((4/3)^((n-1)/2) * BMULT / n)^(n/2) NI sqrt(disc(K)).
    3722             :    * NI <= LIMCMAX^2 */
    3723       12040 :   PREC = maxss(DEFAULTPREC, prec);
    3724       12040 :   if (nf) PREC = maxss(PREC, nf_get_prec(nf));
    3725       12040 :   PREC = maxss(PREC, nbits2prec((long)(LOGD2 * 0.02) + N*N));
    3726       12040 :   if (DEBUGLEVEL) err_printf("PREC = %ld\n", PREC);
    3727       12040 :   small_norm_prec = nbits2prec( BITS_IN_LONG +
    3728       12040 :     (N/2. * ((N-1)/2.*log(4./3) + log(BMULT/(double)N))
    3729       12040 :      + 2*log((double) LIMCMAX) + LOGD/2) / M_LN2 ); /*enough to compute norms*/
    3730       12040 :   if (small_norm_prec > PREC) PREC = small_norm_prec;
    3731       12040 :   if (!nf)
    3732       11627 :     nf = nfinit_complete(&nfT, flag_nfinit, PREC);
    3733         413 :   else if (nf_get_prec(nf) < PREC)
    3734          28 :     nf = nfnewprec_shallow(nf, PREC);
    3735       12040 :   M_sn = nf_get_M(nf);
    3736       12040 :   if (PREC > small_norm_prec) M_sn = gprec_w(M_sn, small_norm_prec);
    3737             : 
    3738       12040 :   zu = nfrootsof1(nf);
    3739       12040 :   gel(zu,2) = nf_to_scalar_or_alg(nf, gel(zu,2));
    3740             : 
    3741       12040 :   nf_get_sign(nf, &R1, &R2); RU = R1+R2;
    3742       12040 :   auts = automorphism_matrices(nf, &cyclic);
    3743       12040 :   F.embperm = automorphism_perms(nf_get_M(nf), auts, cyclic, R1, R2, N);
    3744       12040 :   if (DEBUGLEVEL)
    3745             :   {
    3746           0 :     timer_printf(&T, "nfinit & nfrootsof1");
    3747           0 :     err_printf("%s bnf: R1 = %ld, R2 = %ld\nD = %Ps\n",
    3748             :                flag? "Algebraic": "Floating point", R1,R2, D);
    3749             :   }
    3750       12040 :   if (LOGD < 20.)
    3751             :   { /* tiny disc, Minkowski may be smaller than Bach */
    3752       11529 :     lim = exp(-N + R2 * log(4/M_PI) + LOGD/2) * sqrt(2*M_PI*N);
    3753       11529 :     if (lim < 3) lim = 3;
    3754             :   }
    3755             :   else /* to be ignored */
    3756         511 :     lim = -1;
    3757       12040 :   if (cbach > 12.) {
    3758           0 :     if (cbach2 < cbach) cbach2 = cbach;
    3759           0 :     cbach = 12.;
    3760             :   }
    3761       12040 :   if (cbach < 0.)
    3762           0 :     pari_err_DOMAIN("Buchall","Bach constant","<",gen_0,dbltor(cbach));
    3763             : 
    3764       12040 :   cache.base = NULL; F.subFB = NULL; F.LP = NULL; SUnits = Ce = NULL;
    3765       12040 :   init_GRHcheck(&GRHcheck, N, R1, LOGD);
    3766       12040 :   high = low = LIMC0 = maxss((long)(cbach2*LOGD2), 1);
    3767       58247 :   while (!GRHchk(nf, &GRHcheck, high)) { low = high; high *= 2; }
    3768       46312 :   while (high - low > 1)
    3769             :   {
    3770       34272 :     long test = (low+high)/2;
    3771       34272 :     if (GRHchk(nf, &GRHcheck, test)) high = test; else low = test;
    3772             :   }
    3773       12040 :   LIMC2 = (high == LIMC0+1 && GRHchk(nf, &GRHcheck, LIMC0))? LIMC0: high;
    3774       12040 :   if (LIMC2 > LIMCMAX) LIMC2 = LIMCMAX;
    3775             :   /* Assuming GRH, {P, NP <= LIMC2} generate Cl(K) */
    3776       12040 :   if (DEBUGLEVEL) err_printf("LIMC2 = %ld\n", LIMC2);
    3777       12040 :   LIMC0 = (long)(cbach*LOGD2); /* initial value for LIMC */
    3778       12040 :   LIMC = cbach? LIMC0: LIMC2; /* use {P, NP <= LIMC} as a factorbase */
    3779       12040 :   LIMC = maxss(LIMC, nthideal(&GRHcheck, nf, N));
    3780       12040 :   if (DEBUGLEVEL) timer_printf(&T, "computing Bach constant");
    3781       12040 :   LIMres = primeneeded(N, R1, R2, LOGD);
    3782       12040 :   cache_prime_dec(&GRHcheck, LIMres, nf);
    3783             :   /* invhr ~ 2^r1 (2pi)^r2 / sqrt(D) w * Res(zeta_K, s=1) = 1 / hR */
    3784       24080 :   invhr = gmul(gdiv(gmul2n(powru(mppi(DEFAULTPREC), R2), RU),
    3785       12040 :               mulri(gsqrt(D,DEFAULTPREC),gel(zu,1))),
    3786             :               compute_invres(&GRHcheck, LIMres));
    3787       12040 :   if (DEBUGLEVEL) timer_printf(&T, "computing inverse of hR");
    3788       12040 :   av = avma;
    3789             : 
    3790       12705 : START:
    3791       12705 :   if (DEBUGLEVEL) timer_start(&T);
    3792       12705 :   if (TRIES) LIMC = bnf_increase_LIMC(LIMC,LIMCMAX);
    3793       12705 :   if (DEBUGLEVEL && LIMC > LIMC0)
    3794           0 :     err_printf("%s*** Bach constant: %f\n", TRIES?"\n":"", LIMC/LOGD2);
    3795       12705 :   if (cache.base)
    3796             :   {
    3797             :     REL_t *rel;
    3798       20310 :     for (i = 1, rel = cache.base + 1; rel < cache.last; rel++)
    3799       20268 :       if (rel->m) i++;
    3800          42 :     computed = cgetg(i, t_VEC);
    3801       20310 :     for (i = 1, rel = cache.base + 1; rel < cache.last; rel++)
    3802       20268 :       if (rel->m) gel(computed, i++) = rel->m;
    3803          42 :     computed = gclone(computed); delete_cache(&cache);
    3804             :   }
    3805       12705 :   TRIES++; set_avma(av);
    3806       12705 :   if (F.LP) delete_FB(&F);
    3807       12705 :   if (LIMC2 < LIMC) LIMC2 = LIMC;
    3808       12705 :   if (DEBUGLEVEL) { err_printf("LIMC = %ld, LIMC2 = %ld\n",LIMC,LIMC2); }
    3809             : 
    3810       12705 :   FBgen(&F, nf, N, LIMC, LIMC2, &GRHcheck);
    3811       12705 :   if (!F.KC) goto START;
    3812       12705 :   av = avma;
    3813       12705 :   subFBgen(&F,auts,cyclic,lim < 0? LIMC2: mindd(lim,LIMC2),MINSFB);
    3814       12705 :   if (lg(F.subFB) == 1) goto START;
    3815       12082 :   if (DEBUGLEVEL)
    3816           0 :     timer_printf(&T, "factorbase (#subFB = %ld) and ideal permutations",
    3817           0 :                      lg(F.subFB)-1);
    3818             : 
    3819       12082 :   fact = (FACT*)stack_malloc((F.KC+1)*sizeof(FACT));
    3820       12082 :   PERM = leafcopy(F.perm); /* to be restored in case of precision increase */
    3821       12082 :   cache.basis = zero_Flm_copy(F.KC,F.KC);
    3822       12082 :   small_multiplier = zero_Flv(F.KC);
    3823       12082 :   done_small = small_fail = squash_index = zc = sfb_trials = nreldep = 0;
    3824       12082 :   fail_limit = F.KC + 1;
    3825       12082 :   W = A = R = NULL;
    3826       12082 :   av2 = avma;
    3827       12082 :   init_rel(&cache, &F, RELSUP + RU-1);
    3828       12082 :   old_need = need = cache.end - cache.last;
    3829       12082 :   add_cyclotomic_units(nf, zu, &cache, &F);
    3830       12082 :   if (DEBUGLEVEL) err_printf("\n");
    3831       12082 :   cache.end = cache.last + need;
    3832             : 
    3833       12082 :   if (computed)
    3834             :   {
    3835       10999 :     for (i = 1; i < lg(computed); i++)
    3836       10957 :       try_elt(&cache, &F, nf, gel(computed, i), fact);
    3837          42 :     gunclone(computed);
    3838          42 :     if (DEBUGLEVEL && i > 1)
    3839           0 :       timer_printf(&T, "including already computed relations");
    3840          42 :     need = 0;
    3841             :   }
    3842             : 
    3843             :   do
    3844             :   {
    3845             :     GEN Ar, C0;
    3846             :     do
    3847             :     {
    3848       78257 :       pari_sp av4 = avma;
    3849       78257 :       if (need > 0)
    3850             :       {
    3851       78163 :         long oneed = cache.end - cache.last;
    3852             :         /* Test below can be true if small_norm did not find enough linearly
    3853             :          * dependent relations */
    3854       78163 :         if (need < oneed) need = oneed;
    3855       78163 :         pre_allocate(&cache, need+lg(auts)-1+(R ? lg(W)-1 : 0));
    3856       78163 :         cache.end = cache.last + need;
    3857       78163 :         F.L_jid = trim_list(&F);
    3858             :       }
    3859       78257 :       if (need > 0 && Nrelid > 0 && (done_small <= F.KC+1 || A) &&
    3860       56816 :           small_fail <= fail_limit &&
    3861       56816 :           cache.last < cache.base + 2*F.KC+2*RU+RELSUP /* heuristic */)
    3862             :       {
    3863       52677 :         long j, k, LIE = (R && lg(W) > 1 && (done_small % 2));
    3864       52677 :         REL_t *last = cache.last;
    3865       52677 :         pari_sp av3 = avma;
    3866             :         GEN p0;
    3867       52677 :         if (LIE)
    3868             :         { /* We have full rank for class group and unit. The following tries to
    3869             :            * improve the prime group lattice by looking for relations involving
    3870             :            * the primes generating the class group. */
    3871        1391 :           long n = lg(W)-1; /* need n relations to squash the class group */
    3872        1391 :           F.L_jid = vecslice(F.perm, 1, n);
    3873        1391 :           cache.end = cache.last + n;
    3874             :           /* Lie to the add_rel subsystem: pretend we miss relations involving
    3875             :            * the primes generating the class group (and only those). */
    3876        1391 :           cache.missing = n;
    3877        6904 :           for ( ; n > 0; n--) mael(cache.basis, F.perm[n], F.perm[n]) = 0;
    3878             :         }
    3879       52677 :         j = done_small % (F.KC+1);
    3880       52677 :         if (j == 0) p0 = NULL;
    3881             :         else
    3882             :         {
    3883       40245 :           p0 = gel(F.LP, j);
    3884       40245 :           if (!A)
    3885             :           { /* Prevent considering both P_iP_j and P_jP_i in small_norm */
    3886             :             /* Not all elements end up in F.L_jid (eliminated by hnfspec/add or
    3887             :              * by trim_list): keep track of which ideals are being considered
    3888             :              * at each run. */
    3889       31507 :             long mj = small_multiplier[j];
    3890     1435127 :             for (i = k = 1; i < lg(F.L_jid); i++)
    3891     1403620 :               if (F.L_jid[i] > mj)
    3892             :               {
    3893     1207220 :                 small_multiplier[F.L_jid[i]] = j;
    3894     1207220 :                 F.L_jid[k++] = F.L_jid[i];
    3895             :               }
    3896       31507 :             setlg(F.L_jid, k);
    3897             :           }
    3898             :         }
    3899       52677 :         if (lg(F.L_jid) > 1)
    3900       52397 :           small_norm(&cache, &F, nf, Nrelid, M_sn, fact, p0);
    3901       52677 :         F.L_jid = F.perm; set_avma(av3);
    3902       52677 :         if (!A && cache.last != last) small_fail = 0; else small_fail++;
    3903       52677 :         if (LIE)
    3904             :         { /* restore add_rel subsystem: undo above lie */
    3905        1391 :           long n = lg(W) - 1;
    3906        6904 :           for ( ; n > 0; n--) mael(cache.basis, F.perm[n], F.perm[n]) = 1;
    3907        1391 :           cache.missing = 0;
    3908             :         }
    3909       52677 :         cache.end = cache.last;
    3910       52677 :         done_small++;
    3911       52677 :         need = F.sfb_chg = 0;
    3912             :       }
    3913       78257 :       if (need > 0)
    3914             :       { /* Random relations */
    3915       25486 :         if (++nreldep > F.MAXDEPSIZESFB) {
    3916         120 :           if (++sfb_trials > SFB_MAX && LIMC < LIMCMAX/6) goto START;
    3917          99 :           F.sfb_chg = sfb_INCREASE;
    3918          99 :           nreldep = 0;
    3919             :         }
    3920       25366 :         else if (!(nreldep % F.MAXDEPSFB))
    3921        3067 :           F.sfb_chg = sfb_CHANGE;
    3922       25465 :         if (F.sfb_chg && !subFB_change(&F)) goto START;
    3923       25465 :         rnd_rel(&cache, &F, nf, fact);
    3924       25465 :         F.L_jid = F.perm;
    3925             :       }
    3926       78236 :       if (DEBUGLEVEL) timer_start(&T);
    3927       78236 :       if (precpb)
    3928             :       {
    3929             :         REL_t *rel;
    3930          99 :         if (DEBUGLEVEL)
    3931             :         {
    3932           0 :           char str[64]; sprintf(str,"Buchall_param (%s)",precpb);
    3933           0 :           pari_warn(warnprec,str,PREC);
    3934             :         }
    3935          99 :         nf = _nfnewprec(nf, PREC, &nfisclone);
    3936          99 :         precdouble++; precpb = NULL;
    3937             : 
    3938          99 :         if (flag)
    3939             :         { /* recompute embs only, no need to redo HNF */
    3940          31 :           long j, le = lg(embs), lC = lg(C);
    3941          31 :           GEN E, M = nf_get_M(nf);
    3942          31 :           set_avma(av4);
    3943       10029 :           for (rel = cache.base+1, i = 1; i < le; i++,rel++)
    3944        9998 :             gel(embs,i) = rel_embed(rel, &F, embs, i, M, RU, R1, PREC);
    3945          31 :           E = RgM_ZM_mul(embs, rowslice(C, RU+1, nbrows(C)));
    3946       10029 :           for (j = 1; j < lC; j++)
    3947       38348 :             for (i = 1; i <= RU; i++) gcoeff(C,i,j) = gcoeff(E,i,j);
    3948          31 :           av4 = avma;
    3949             :         }
    3950             :         else
    3951             :         { /* recompute embs + HNF */
    3952       18837 :           for(i = 1; i < lg(PERM); i++) F.perm[i] = PERM[i];
    3953          68 :           cache.chk = cache.base;
    3954          68 :           W = NULL;
    3955             :         }
    3956          99 :         if (DEBUGLEVEL) timer_printf(&T, "increasing accuracy");
    3957             :       }
    3958       78236 :       set_avma(av4);
    3959       78236 :       if (cache.chk != cache.last)
    3960             :       { /* Reduce relation matrices */
    3961       31266 :         long l = cache.last - cache.chk + 1, j;
    3962       31266 :         GEN mat = cgetg(l, t_MAT);
    3963             :         REL_t *rel;
    3964             : 
    3965      264397 :         for (j=1,rel = cache.chk + 1; j < l; rel++,j++) gel(mat,j) = rel->R;
    3966       31266 :         if (!flag || W)
    3967             :         {
    3968       29621 :           embs = get_embs(&F, &cache, nf, embs, PREC);
    3969       29621 :           if (DEBUGLEVEL && timer_get(&T) > 1)
    3970           0 :             timer_printf(&T, "floating point embeddings");
    3971             :         }
    3972       31266 :         if (!W)
    3973             :         { /* never reduced before */
    3974       12150 :           C = flag? matbotid(&cache): embs;
    3975       12150 :           W = hnfspec_i(mat, F.perm, &dep, &B, &C, F.subFB ? lg(F.subFB)-1:0);
    3976       12150 :           if (DEBUGLEVEL)
    3977           0 :             timer_printf(&T, "hnfspec [%ld x %ld]", lg(F.perm)-1, l-1);
    3978       12150 :           if (flag)
    3979             :           {
    3980        1645 :             PREC += nbits2extraprec(gexpo(C));
    3981        1645 :             if (nf_get_prec(nf) < PREC) nf = _nfnewprec(nf, PREC, &nfisclone);
    3982        1645 :             embs = get_embs(&F, &cache, nf, embs, PREC);
    3983        1645 :             C = vconcat(RgM_ZM_mul(embs, C), C);
    3984             :           }
    3985       12150 :           if (DEBUGLEVEL)
    3986           0 :             timer_printf(&T, "hnfspec floating points");
    3987             :         }
    3988             :         else
    3989             :         {
    3990       19116 :           long k = lg(embs);
    3991       19116 :           GEN E = vecslice(embs, k-l+1,k-1);
    3992       19116 :           if (flag)
    3993             :           {
    3994       10286 :             E = matbotidembs(&cache, E);
    3995       10286 :             matenlarge(C, cache.last - cache.chk);
    3996             :           }
    3997       19116 :           W = hnfadd_i(W, F.perm, &dep, &B, &C, mat, E);
    3998       19116 :           if (DEBUGLEVEL)
    3999           0 :             timer_printf(&T, "hnfadd (%ld + %ld)", l-1, lg(dep)-1);
    4000             :         }
    4001       31266 :         gerepileall(av2, 5, &W,&C,&B,&dep,&embs);
    4002       31266 :         cache.chk = cache.last;
    4003             :       }
    4004       46970 :       else if (!W)
    4005             :       {
    4006           0 :         need = old_need;
    4007           0 :         F.L_jid = vecslice(F.perm, 1, need);
    4008           0 :         continue;
    4009             :       }
    4010       78236 :       need = F.KC - (lg(W)-1) - (lg(B)-1);
    4011       78236 :       if (!need && cache.missing)
    4012             :       { /* The test above will never be true except if 27449|class number.
    4013             :          * Ensure that if we have maximal rank for the ideal lattice, then
    4014             :          * cache.missing == 0. */
    4015          14 :         for (i = 1; cache.missing; i++)
    4016           7 :           if (!mael(cache.basis, i, i))
    4017             :           {
    4018             :             long j;
    4019           7 :             cache.missing--; mael(cache.basis, i, i) = 1;
    4020         427 :             for (j = i+1; j <= F.KC; j++) mael(cache.basis, j, i) = 0;
    4021             :           }
    4022             :       }
    4023       78236 :       zc = (lg(C)-1) - (lg(B)-1) - (lg(W)-1);
    4024       78236 :       if (RU-1-zc > 0) need = minss(need + RU-1-zc, F.KC); /* for units */
    4025       78236 :       if (need)
    4026             :       { /* dependent rows */
    4027       52293 :         F.L_jid = vecslice(F.perm, 1, need);
    4028       52293 :         vecsmall_sort(F.L_jid);
    4029       52293 :         if (need != old_need) { nreldep = 0; old_need = need; }
    4030             :       }
    4031             :       else
    4032             :       { /* If the relation lattice is too small, check will be > 1 and we will
    4033             :          * do a new run of small_norm/rnd_rel asking for 1 relation. This often
    4034             :          * gives a relation involving L_jid[1]. We rotate the first element of
    4035             :          * L_jid in order to increase the probability of finding relations that
    4036             :          * increases the lattice. */
    4037       25943 :         long j, n = lg(W) - 1;
    4038       30912 :         if (n > 1 && squash_index % n)
    4039             :         {
    4040        4969 :           F.L_jid = leafcopy(F.perm);
    4041       28689 :           for (j = 1; j <= n; j++)
    4042       23720 :             F.L_jid[j] = F.perm[1 + (j + squash_index - 1) % n];
    4043             :         }
    4044             :         else
    4045       20974 :           F.L_jid = F.perm;
    4046       25943 :         squash_index++;
    4047             :       }
    4048             :     }
    4049       78236 :     while (need);
    4050             : 
    4051       25943 :     if (!A)
    4052             :     {
    4053       12075 :       small_fail = old_need = 0;
    4054       12075 :       fail_limit = maxss(F.KC / FAIL_DIVISOR, MINFAIL);
    4055             :     }
    4056       25943 :     A = vecslice(C, 1, zc); /* cols corresponding to units */
    4057       25943 :     if (flag) A = rowslice(A, 1, RU);
    4058       25943 :     Ar = real_i(A);
    4059       25943 :     R = compute_multiple_of_R(Ar, RU, N, &need, &bit, &lambda);
    4060       25943 :     if (need < old_need) small_fail = 0;
    4061             : #if 0 /* A good idea if we are indeed stuck but needs tuning */
    4062             :     /* we have computed way more relations than should be necessary */
    4063             :     if (TRIES < 3 && LIMC < LIMCMAX / 24 &&
    4064             :                      cache.last - cache.base > 10 * F.KC) goto START;
    4065             : #endif
    4066       25943 :     old_need = need;
    4067       25943 :     if (!lambda)
    4068          81 :     { precpb = "bestappr"; PREC = myprecdbl(PREC, flag? C: NULL); continue; }
    4069       25862 :     if (!R)
    4070             :     { /* not full rank for units */
    4071       10261 :       if (!need)
    4072           0 :       { precpb = "regulator"; PREC = myprecdbl(PREC, flag? C: NULL); }
    4073       10261 :       continue;
    4074             :     }
    4075       15601 :     h = ZM_det_triangular(W);
    4076       15601 :     if (DEBUGLEVEL) err_printf("\n#### Tentative class number: %Ps\n", h);
    4077       15601 :     i = compute_R(lambda, mulir(h,invhr), flag? 0: RgM_bit(C, bit), &L, &R);
    4078       15601 :     if (DEBUGLEVEL)
    4079             :     {
    4080           0 :       err_printf("\n");
    4081           0 :       timer_printf(&T, "computing regulator and check");
    4082             :     }
    4083       15601 :     switch(i)
    4084             :     {
    4085        3522 :       case fupb_RELAT:
    4086        3522 :         need = 1; /* not enough relations */
    4087        3522 :         continue;
    4088          17 :       case fupb_PRECI: /* prec problem unless we cheat on Bach constant */
    4089          17 :         if ((precdouble&7) == 7 && LIMC<=LIMCMAX/6) goto START;
    4090          17 :         precpb = "compute_R"; PREC = myprecdbl(PREC, flag? C: NULL);
    4091          17 :         continue;
    4092             :     }
    4093             :     /* DONE */
    4094             : 
    4095       12062 :     if (F.KCZ2 > F.KCZ)
    4096             :     {
    4097          28 :       if (F.sfb_chg && !subFB_change(&F)) goto START;
    4098          28 :       if (!be_honest(&F, nf, auts, fact)) goto START;
    4099           7 :       if (DEBUGLEVEL) timer_printf(&T, "to be honest");
    4100             :     }
    4101       12041 :     F.KCZ2 = 0; /* be honest only once */
    4102             : 
    4103             :     /* fundamental units */
    4104             :     {
    4105       12041 :       GEN AU, CU, U, v = extract_full_lattice(L); /* L may be large */
    4106       12041 :       CU = NULL;
    4107       12041 :       if (v) { A = vecpermute(A, v); L = vecpermute(L, v); }
    4108             :       /* arch. components of fund. units */
    4109       12041 :       U = ZM_lll(L, 0.99, LLL_IM);
    4110       12041 :       U = ZM_mul(U, lll(RgM_ZM_mul(real_i(A), U)));
    4111       12041 :       AU = RgM_ZM_mul(A, U);
    4112       12041 :       A = cleanarch(AU, N, PREC);
    4113       12041 :       if (DEBUGLEVEL) timer_printf(&T, "units LLL + cleanarch");
    4114       12041 :       if (!A || (lg(A) > 1 && gprecision(A) <= 2))
    4115             :       {
    4116           0 :         long add = nbits2extraprec( gexpo(AU) + 64 ) - gprecision(AU);
    4117           0 :         precpb = "cleanarch"; PREC += maxss(add, 1); continue;
    4118             :       }
    4119       12041 :       if (flag)
    4120             :       {
    4121        1638 :         long l = lgcols(C) - RU;
    4122             :         REL_t *rel;
    4123        1638 :         SUnits = cgetg(l, t_COL);
    4124       80438 :         for (rel = cache.base+1, i = 1; i < l; i++,rel++)
    4125       78800 :           set_rel_alpha(rel, auts, SUnits, i);
    4126        1638 :         if (RU > 1)
    4127             :         {
    4128        1351 :           GEN c = v? vecpermute(C,v): vecslice(C,1,zc);
    4129        1351 :           CU = ZM_mul(rowslice(c, RU+1, nbrows(c)), U);
    4130             :         }
    4131             :       }
    4132       12041 :       if (DEBUGLEVEL) err_printf("\n#### Computing fundamental units\n");
    4133       12041 :       fu = getfu(nf, &A, CU? &U: NULL, PREC);
    4134       12041 :       CU = CU? ZM_mul(CU, U): cgetg(1, t_MAT);
    4135       12041 :       if (DEBUGLEVEL) timer_printf(&T, "getfu");
    4136       12041 :       Ce = vecslice(C, zc+1, lg(C)-1);
    4137       12041 :       if (flag) SUnits = mkvec4(SUnits, CU, rowslice(Ce, RU+1, nbrows(Ce)),
    4138             :                                 utoipos(LIMC));
    4139             :     }
    4140             :     /* class group generators */
    4141       12041 :     if (flag) Ce = rowslice(Ce, 1, RU);
    4142       12041 :     C0 = Ce; Ce = cleanarch(Ce, N, PREC);
    4143       12041 :     if (!Ce) {
    4144           1 :       long add = nbits2extraprec( gexpo(C0) + 64 ) - gprecision(C0);
    4145           1 :       precpb = "cleanarch"; PREC += maxss(add, 1);
    4146             :     }
    4147       12041 :     if (DEBUGLEVEL) timer_printf(&T, "cleanarch");
    4148       25922 :   } while (need || precpb);
    4149             : 
    4150       12040 :   Vbase = vecpermute(F.LP, F.perm);
    4151       12040 :   if (!fu) fu = cgetg(1, t_MAT);
    4152       12040 :   if (!SUnits) SUnits = gen_1;
    4153       12040 :   clg1 = class_group_gen(nf,W,Ce,Vbase,PREC, &clg2);
    4154       12040 :   res = mkvec5(clg1, R, SUnits, zu, fu);
    4155       12040 :   res = buchall_end(nf,res,clg2,W,B,A,Ce,Vbase);
    4156       12040 :   delete_FB(&F);
    4157       12040 :   res = gerepilecopy(av0, res);
    4158       12040 :   if (flag) obj_insert_shallow(res, MATAL, cgetg(1,t_VEC));
    4159       12040 :   if (precdouble) gunclone(nf);
    4160       12040 :   delete_cache(&cache);
    4161       12040 :   free_GRHcheck(&GRHcheck);
    4162       12040 :   return res;
    4163             : }

Generated by: LCOV version 1.13