Line data Source code
1 : /* Copyright (C) 2000 The PARI group.
2 :
3 : This file is part of the PARI/GP package.
4 :
5 : PARI/GP is free software; you can redistribute it and/or modify it under the
6 : terms of the GNU General Public License as published by the Free Software
7 : Foundation; either version 2 of the License, or (at your option) any later
8 : version. It is distributed in the hope that it will be useful, but WITHOUT
9 : ANY WARRANTY WHATSOEVER.
10 :
11 : Check the License for details. You should have received a copy of it, along
12 : with the package; see the file 'COPYING'. If not, write to the Free Software
13 : Foundation, Inc., 51 Franklin Street, Fifth Floor, Boston, MA 02110-1301 USA. */
14 : #include "pari.h"
15 : #include "paripriv.h"
16 :
17 : #define DEBUGLEVEL DEBUGLEVEL_bnf
18 :
19 : /*******************************************************************/
20 : /* */
21 : /* CLASS GROUP AND REGULATOR (McCURLEY, BUCHMANN) */
22 : /* GENERAL NUMBER FIELDS */
23 : /* */
24 : /*******************************************************************/
25 : /* get_random_ideal */
26 : static const long RANDOM_BITS = 4;
27 : /* Buchall */
28 : static const long RELSUP = 5;
29 : static const long FAIL_DIVISOR = 32;
30 : static const long MINFAIL = 10;
31 : /* small_norm */
32 : static const long BNF_RELPID = 4;
33 : static const long MAXTRY_FACT = 500;
34 : /* rnd_rel */
35 : static const long RND_REL_RELPID = 1;
36 : /* random relations */
37 : static const long MINSFB = 3;
38 : static const long SFB_MAX = 3;
39 : static const long DEPSIZESFBMULT = 16;
40 : static const long DEPSFBDIV = 10;
41 : /* add_rel_i */
42 : static const ulong mod_p = 27449UL;
43 : /* be_honest */
44 : static const long maxtry_HONEST = 50;
45 :
46 : typedef struct FACT {
47 : long pr, ex;
48 : } FACT;
49 :
50 : typedef struct subFB_t {
51 : GEN subFB;
52 : struct subFB_t *old;
53 : } subFB_t;
54 :
55 : /* a factor base contains only noninert primes
56 : * KC = # of P in factor base (p <= n, NP <= n2)
57 : * KC2= # of P assumed to generate class group (NP <= n2)
58 : *
59 : * KCZ = # of rational primes under ideals counted by KC
60 : * KCZ2= same for KC2 */
61 :
62 : typedef struct FB_t {
63 : GEN FB; /* FB[i] = i-th rational prime used in factor base */
64 : GEN LP; /* vector of all prime ideals in FB, by increasing norm */
65 : GEN LV; /* LV[p] = vector of P|p, NP <= n2
66 : * isclone() is set for LV[p] iff all P|p are in FB
67 : * LV[i], i not prime or i > n2, is undefined! */
68 : GEN iLP; /* iLP[p] = i such that LV[p] = [LP[i],...] */
69 : GEN L_jid; /* indexes of "useful" prime ideals for rnd_rel */
70 : long KC, KCZ, KCZ2;
71 : GEN prodZ; /* product of the primes in KCZ*/
72 : GEN subFB; /* LP o subFB = part of FB used to build random relations */
73 : int sfb_chg; /* need to change subFB ? */
74 : GEN perm; /* permutation of LP used to represent relations [updated by
75 : hnfspec/hnfadd: dense rows come first] */
76 : GEN idealperm; /* permutation of ideals under field automorphisms */
77 : GEN minidx; /* minidx[i] min ideal in orbit of LP[i] under field autom */
78 : subFB_t *allsubFB; /* all subFB's used */
79 : GEN embperm; /* permutations of the complex embeddings */
80 : long MAXDEPSIZESFB; /* # trials before increasing subFB */
81 : long MAXDEPSFB; /* MAXDEPSIZESFB / DEPSFBDIV, # trials befor rotating subFB */
82 : double ballvol;
83 : } FB_t;
84 :
85 : enum { sfb_CHANGE = 1, sfb_INCREASE = 2 };
86 :
87 : typedef struct REL_t {
88 : GEN R; /* relation vector as t_VECSMALL; clone */
89 : long nz; /* index of first nonzero elt in R (hash) */
90 : GEN m; /* pseudo-minimum yielding the relation; clone */
91 : long relorig; /* relation this one is an image of */
92 : long relaut; /* automorphim used to compute this relation from the original */
93 : GEN emb; /* archimedean embeddings */
94 : GEN junk[2]; /*make sure sizeof(struct) is a power of two.*/
95 : } REL_t;
96 :
97 : typedef struct RELCACHE_t {
98 : REL_t *chk; /* last checkpoint */
99 : REL_t *base; /* first rel found */
100 : REL_t *last; /* last rel found so far */
101 : REL_t *end; /* target for last relation. base <= last <= end */
102 : size_t len; /* number of rels pre-allocated in base */
103 : long relsup; /* how many linearly dependent relations we allow */
104 : GEN basis; /* mod p basis (generating family actually) */
105 : ulong missing; /* missing vectors in generating family above */
106 : } RELCACHE_t;
107 :
108 : typedef struct FP_t {
109 : double **q, *v, *y, *z;
110 : GEN x;
111 : } FP_t;
112 :
113 : static void
114 0 : wr_rel(GEN e)
115 : {
116 0 : long i, l = lg(e);
117 0 : for (i = 1; i < l; i++)
118 0 : if (e[i]) err_printf("%ld^%ld ",i,e[i]);
119 0 : }
120 : static void
121 0 : dbg_newrel(RELCACHE_t *cache)
122 : {
123 0 : if (DEBUGLEVEL > 1)
124 : {
125 0 : err_printf("\n++++ cglob = %ld\nrel = ", cache->last - cache->base);
126 0 : wr_rel(cache->last->R);
127 0 : err_printf("\n");
128 : }
129 : else
130 0 : err_printf("%ld ", cache->last - cache->base);
131 0 : }
132 :
133 : static void
134 64074 : delete_cache(RELCACHE_t *M)
135 : {
136 : REL_t *rel;
137 1057370 : for (rel = M->base+1; rel <= M->last; rel++)
138 : {
139 993295 : gunclone(rel->R);
140 993294 : if (rel->m) gunclone(rel->m);
141 : }
142 64075 : pari_free((void*)M->base); M->base = NULL;
143 64074 : }
144 :
145 : static void
146 66251 : delete_FB(FB_t *F)
147 : {
148 : subFB_t *s, *sold;
149 133258 : for (s = F->allsubFB; s; s = sold) { sold = s->old; pari_free(s); }
150 66251 : gunclone(F->minidx);
151 66251 : gunclone(F->idealperm);
152 66251 : }
153 :
154 : static void
155 64072 : reallocate(RELCACHE_t *M, long len)
156 : {
157 64072 : M->len = len;
158 64072 : if (!M->base)
159 64074 : M->base = (REL_t*)pari_malloc((len+1) * sizeof(REL_t));
160 : else
161 : {
162 0 : size_t l = M->last - M->base, c = M->chk - M->base, e = M->end - M->base;
163 0 : pari_realloc_ip((void**)&M->base, (len+1) * sizeof(REL_t));
164 0 : M->last = M->base + l;
165 0 : M->chk = M->base + c;
166 0 : M->end = M->base + e;
167 : }
168 64074 : }
169 :
170 : #define pr_get_smallp(pr) gel(pr,1)[2]
171 :
172 : /* don't take P|p all other Q|p are already there */
173 : static int
174 307568 : bad_subFB(FB_t *F, long t)
175 : {
176 307568 : GEN LP, P = gel(F->LP,t);
177 307568 : long p = pr_get_smallp(P);
178 307568 : LP = gel(F->LV,p);
179 307568 : return (isclone(LP) && t == F->iLP[p] + lg(LP)-1);
180 : }
181 :
182 : static void
183 67007 : assign_subFB(FB_t *F, GEN yes, long iyes)
184 : {
185 67007 : long i, lv = sizeof(subFB_t) + iyes*sizeof(long); /* for struct + GEN */
186 67007 : subFB_t *s = (subFB_t *)pari_malloc(lv);
187 67007 : s->subFB = (GEN)&s[1];
188 67007 : s->old = F->allsubFB; F->allsubFB = s;
189 288827 : for (i = 0; i < iyes; i++) s->subFB[i] = yes[i];
190 67007 : F->subFB = s->subFB;
191 67007 : F->MAXDEPSIZESFB = (iyes-1) * DEPSIZESFBMULT;
192 67007 : F->MAXDEPSFB = F->MAXDEPSIZESFB / DEPSFBDIV;
193 67007 : }
194 :
195 : /* Determine the permutation of the ideals made by each field automorphism */
196 : static GEN
197 66250 : FB_aut_perm(FB_t *F, GEN auts, GEN cyclic)
198 : {
199 66250 : long i, j, m, KC = F->KC, nauts = lg(auts)-1;
200 66250 : GEN minidx, perm = zero_Flm_copy(KC, nauts);
201 :
202 66249 : if (!nauts) { F->minidx = gclone(identity_zv(KC)); return cgetg(1,t_MAT); }
203 41959 : minidx = zero_Flv(KC);
204 91730 : for (m = 1; m < lg(cyclic); m++)
205 : {
206 49771 : GEN thiscyc = gel(cyclic, m);
207 49771 : long k0 = thiscyc[1];
208 49771 : GEN aut = gel(auts, k0), permk0 = gel(perm, k0), ppermk;
209 49771 : i = 1;
210 214408 : while (i <= KC)
211 : {
212 164637 : pari_sp av2 = avma;
213 164637 : GEN seen = zero_Flv(KC), P = gel(F->LP, i);
214 164637 : long imin = i, p, f, l;
215 164637 : p = pr_get_smallp(P);
216 164637 : f = pr_get_f(P);
217 : do
218 : {
219 482847 : if (++i > KC) break;
220 433075 : P = gel(F->LP, i);
221 : }
222 433075 : while (p == pr_get_smallp(P) && f == pr_get_f(P));
223 647471 : for (j = imin; j < i; j++)
224 : {
225 482845 : GEN img = ZM_ZC_mul(aut, pr_get_gen(gel(F->LP, j)));
226 1672180 : for (l = imin; l < i; l++)
227 1672180 : if (!seen[l] && ZC_prdvd(img, gel(F->LP, l)))
228 : {
229 482833 : seen[l] = 1; permk0[j] = l; break;
230 : }
231 : }
232 164626 : set_avma(av2);
233 : }
234 68955 : for (ppermk = permk0, i = 2; i < lg(thiscyc); i++)
235 : {
236 19184 : GEN permk = gel(perm, thiscyc[i]);
237 384169 : for (j = 1; j <= KC; j++) permk[j] = permk0[ppermk[j]];
238 19184 : ppermk = permk;
239 : }
240 : }
241 311048 : for (j = 1; j <= KC; j++)
242 : {
243 269089 : if (minidx[j]) continue;
244 129399 : minidx[j] = j;
245 361849 : for (i = 1; i <= nauts; i++) minidx[coeff(perm, j, i)] = j;
246 : }
247 41959 : F->minidx = gclone(minidx); return perm;
248 : }
249 :
250 : /* set subFB.
251 : * Fill F->perm (if != NULL): primes ideals sorted by increasing norm (except
252 : * the ones in subFB come first [dense rows for hnfspec]) */
253 : static void
254 66247 : subFBgen(FB_t *F, GEN auts, GEN cyclic, double PROD, long minsFB)
255 : {
256 : GEN y, perm, yes, no;
257 66247 : long i, j, k, iyes, ino, lv = F->KC + 1;
258 : double prod;
259 : pari_sp av;
260 :
261 66247 : F->LP = cgetg(lv, t_VEC);
262 66247 : F->L_jid = F->perm = cgetg(lv, t_VECSMALL);
263 66247 : av = avma;
264 66247 : y = cgetg(lv,t_COL); /* Norm P */
265 313405 : for (k=0, i=1; i <= F->KCZ; i++)
266 : {
267 247155 : GEN LP = gel(F->LV,F->FB[i]);
268 247155 : long l = lg(LP);
269 714684 : for (j = 1; j < l; j++)
270 : {
271 467526 : GEN P = gel(LP,j);
272 467526 : k++;
273 467526 : gel(y,k) = pr_norm(P);
274 467529 : gel(F->LP,k) = P;
275 : }
276 : }
277 : /* perm sorts LP by increasing norm */
278 66250 : perm = indexsort(y);
279 66250 : no = cgetg(lv, t_VECSMALL); ino = 1;
280 66250 : yes = cgetg(lv, t_VECSMALL); iyes = 1;
281 66250 : prod = 1.0;
282 303468 : for (i = 1; i < lv; i++)
283 : {
284 273494 : long t = perm[i];
285 273494 : if (bad_subFB(F, t)) { no[ino++] = t; continue; }
286 :
287 153012 : yes[iyes++] = t;
288 153012 : prod *= (double)itos(gel(y,t));
289 153012 : if (iyes > minsFB && prod > PROD) break;
290 : }
291 66250 : setlg(yes, iyes);
292 219259 : for (j=1; j<iyes; j++) F->perm[j] = yes[j];
293 186732 : for (i=1; i<ino; i++, j++) F->perm[j] = no[i];
294 260302 : for ( ; j<lv; j++) F->perm[j] = perm[j];
295 66249 : F->allsubFB = NULL;
296 66249 : F->idealperm = gclone(FB_aut_perm(F, auts, cyclic));
297 66251 : if (iyes) assign_subFB(F, yes, iyes);
298 66251 : set_avma(av);
299 66251 : }
300 : static int
301 26419 : subFB_change(FB_t *F)
302 : {
303 26419 : long i, iyes, minsFB = lg(F->subFB)-1, lv = F->KC + 1;
304 26419 : pari_sp av = avma;
305 26419 : GEN yes, L_jid = F->L_jid, present = zero_zv(lv-1);
306 :
307 26419 : if (F->sfb_chg == sfb_INCREASE) minsFB++;
308 :
309 26419 : yes = cgetg(minsFB+1, t_VECSMALL); iyes = 1;
310 26419 : if (L_jid)
311 : {
312 33994 : for (i = 1; i < lg(L_jid); i++)
313 : {
314 33861 : long l = L_jid[i];
315 33861 : if (bad_subFB(F, l)) continue;
316 31636 : yes[iyes++] = l;
317 31636 : present[l] = 1;
318 31636 : if (iyes > minsFB) break;
319 : }
320 : }
321 0 : else i = 1;
322 26419 : if (iyes <= minsFB)
323 : {
324 294 : for ( ; i < lv; i++)
325 : {
326 221 : long l = F->perm[i];
327 221 : if (present[l] || bad_subFB(F, l)) continue;
328 67 : yes[iyes++] = l;
329 67 : if (iyes > minsFB) break;
330 : }
331 133 : if (i == lv) return 0;
332 : }
333 26346 : if (zv_equal(F->subFB, yes))
334 : {
335 25590 : if (DEBUGLEVEL) err_printf("\n*** NOT Changing sub factor base\n");
336 : }
337 : else
338 : {
339 756 : if (DEBUGLEVEL) err_printf("\n*** Changing sub factor base\n");
340 756 : assign_subFB(F, yes, iyes);
341 : }
342 26346 : F->sfb_chg = 0; return gc_bool(av, 1);
343 : }
344 :
345 : /* make sure enough room to store n more relations */
346 : static void
347 122448 : pre_allocate(RELCACHE_t *cache, size_t n)
348 : {
349 122448 : size_t len = (cache->last - cache->base) + n;
350 122448 : if (len >= cache->len) reallocate(cache, len << 1);
351 122448 : }
352 :
353 : void
354 134434 : init_GRHcheck(GRHcheck_t *S, long N, long R1, double LOGD)
355 : {
356 134434 : const double c1 = M_PI*M_PI/2;
357 134434 : const double c2 = 3.663862376709;
358 134434 : const double c3 = 3.801387092431; /* Euler + log(8*Pi)*/
359 134434 : S->clone = 0;
360 134434 : S->cN = R1*c2 + N*c1;
361 134434 : S->cD = LOGD - N*c3 - R1*M_PI/2;
362 134434 : S->maxprimes = 16000; /* sufficient for LIMC=176081*/
363 134434 : S->primes = (GRHprime_t*)pari_malloc(S->maxprimes*sizeof(*S->primes));
364 134434 : S->nprimes = 0;
365 134434 : S->limp = 0;
366 134434 : u_forprime_init(&S->P, 2, ULONG_MAX);
367 134433 : }
368 :
369 : void
370 134434 : free_GRHcheck(GRHcheck_t *S)
371 : {
372 134434 : if (S->clone)
373 : {
374 64001 : long i = S->nprimes;
375 : GRHprime_t *pr;
376 7578954 : for (pr = S->primes, i = S->nprimes; i > 0; pr++, i--) gunclone(pr->dec);
377 : }
378 134436 : pari_free(S->primes);
379 134434 : }
380 :
381 : int
382 1530860 : GRHok(GRHcheck_t *S, double L, double SA, double SB)
383 : {
384 1530860 : return (S->cD + (S->cN + 2*SB) / L - 2*SA < -1e-8);
385 : }
386 :
387 : /* Return factorization pattern of p: [f,n], where n[i] primes of
388 : * residue degree f[i] */
389 : static GEN
390 7511614 : get_fs(GEN nf, GEN P, GEN index, ulong p)
391 : {
392 : long j, k, f, n, l;
393 : GEN fs, ns;
394 :
395 7511614 : if (umodiu(index, p))
396 : { /* easy case: p does not divide index */
397 7472701 : GEN F = Flx_degfact(ZX_to_Flx(P,p), p);
398 7473552 : fs = gel(F,1); l = lg(fs);
399 : }
400 : else
401 : {
402 38507 : GEN F = idealprimedec(nf, utoipos(p));
403 38714 : l = lg(F);
404 38714 : fs = cgetg(l, t_VECSMALL);
405 121177 : for (j = 1; j < l; j++) fs[j] = pr_get_f(gel(F,j));
406 : }
407 7512266 : ns = cgetg(l, t_VECSMALL);
408 7509838 : f = fs[1]; n = 1;
409 13911797 : for (j = 2, k = 1; j < l; j++)
410 6401959 : if (fs[j] == f)
411 4677324 : n++;
412 : else
413 : {
414 1724635 : ns[k] = n; fs[k] = f; k++;
415 1724635 : f = fs[j]; n = 1;
416 : }
417 7509838 : ns[k] = n; fs[k] = f; k++;
418 7509838 : setlg(fs, k);
419 7508889 : setlg(ns, k); return mkvec2(fs,ns);
420 : }
421 :
422 : /* cache data for all rational primes up to the LIM */
423 : static void
424 921200 : cache_prime_dec(GRHcheck_t *S, ulong LIM, GEN nf)
425 : {
426 921200 : pari_sp av = avma;
427 : GRHprime_t *pr;
428 : GEN index, P;
429 : double nb;
430 :
431 921200 : if (S->limp >= LIM) return;
432 329614 : S->clone = 1;
433 329614 : nb = primepi_upper_bound((double)LIM); /* #{p <= LIM} <= nb */
434 329620 : GRH_ensure(S, nb+1); /* room for one extra prime */
435 329621 : P = nf_get_pol(nf);
436 329620 : index = nf_get_index(nf);
437 329619 : for (pr = S->primes + S->nprimes;;)
438 7182384 : {
439 7512003 : ulong p = u_forprime_next(&(S->P));
440 7511671 : pr->p = p;
441 7511671 : pr->logp = log((double)p);
442 7511671 : pr->dec = gclone(get_fs(nf, P, index, p));
443 7512093 : S->nprimes++;
444 7512093 : pr++;
445 7512093 : set_avma(av);
446 : /* store up to nextprime(LIM) included */
447 7512007 : if (p >= LIM) { S->limp = p; break; }
448 : }
449 : }
450 :
451 : static double
452 2258439 : tailresback(long R1, long R2, double rK, long C, double C2, double C3, double r1K, double r2K, double logC, double logC2, double logC3)
453 : {
454 2258439 : const double rQ = 1.83787706641;
455 2258439 : const double r1Q = 1.98505372441;
456 2258439 : const double r2Q = 1.07991541347;
457 4516878 : return fabs((R1+R2-1)*(12*logC3+4*logC2-9*logC-6)/(2*C*logC3)
458 2258439 : + (rK-rQ)*(6*logC2 + 5*logC + 2)/(C*logC3)
459 2258439 : - R2*(6*logC2+11*logC+6)/(C2*logC2)
460 2258439 : - 2*(r1K-r1Q)*(3*logC2 + 4*logC + 2)/(C2*logC3)
461 2258439 : + (R1+R2-1)*(12*logC3+40*logC2+45*logC+18)/(6*C3*logC3)
462 2258439 : + (r2K-r2Q)*(2*logC2 + 3*logC + 2)/(C3*logC3));
463 : }
464 :
465 : static double
466 1129233 : tailres(long R1, long R2, double al2K, double rKm, double rKM, double r1Km,
467 : double r1KM, double r2Km, double r2KM, double C, long i)
468 : {
469 : /* C >= 3*2^i, lower bound for eint1(log(C)/2) */
470 : /* for(i=0,30,print(eint1(log(3*2^i)/2))) */
471 : static double tab[] = {
472 : 0.50409264803,
473 : 0.26205336997,
474 : 0.14815491171,
475 : 0.08770540561,
476 : 0.05347651832,
477 : 0.03328934284,
478 : 0.02104510690,
479 : 0.01346475900,
480 : 0.00869778586,
481 : 0.00566279855,
482 : 0.00371111950,
483 : 0.00244567837,
484 : 0.00161948049,
485 : 0.00107686891,
486 : 0.00071868750,
487 : 0.00048119961,
488 : 0.00032312188,
489 : 0.00021753772,
490 : 0.00014679818,
491 : 9.9272855581E-5,
492 : 6.7263969995E-5,
493 : 4.5656812967E-5,
494 : 3.1041124593E-5,
495 : 2.1136011590E-5,
496 : 1.4411645381E-5,
497 : 9.8393304088E-6,
498 : 6.7257395409E-6,
499 : 4.6025878272E-6,
500 : 3.1529719271E-6,
501 : 2.1620490021E-6,
502 : 1.4839266071E-6
503 : };
504 1129233 : const double logC = log(C), logC2 = logC*logC, logC3 = logC*logC2;
505 1129233 : const double C2 = C*C, C3 = C*C2;
506 1129233 : double E1 = i >30? 0: tab[i];
507 1129233 : return al2K*((33*logC2+22*logC+8)/(8*logC3*sqrt(C))+15*E1/16)
508 1129233 : + maxdd(tailresback(rKm,r1KM,r2Km, C,C2,C3,R1,R2,logC,logC2,logC3),
509 1129239 : tailresback(rKM,r1Km,r2KM, C,C2,C3,R1,R2,logC,logC2,logC3))/2
510 1129239 : + ((R1+R2-1)*4*C+R2)*(C2+6*logC)/(4*C2*C2*logC2);
511 : }
512 :
513 : static long
514 64000 : primeneeded(long N, long R1, long R2, double LOGD)
515 : {
516 64000 : const double lim = 0.25; /* should be log(2)/2 == 0.34657... */
517 64000 : const double al2K = 0.3526*LOGD - 0.8212*N + 4.5007;
518 64000 : const double rKm = -1.0155*LOGD + 2.1042*N - 8.3419;
519 64000 : const double rKM = -0.5 *LOGD + 1.2076*N + 1;
520 64000 : const double r1Km = - LOGD + 1.4150*N;
521 64000 : const double r1KM = - LOGD + 1.9851*N;
522 64000 : const double r2Km = - LOGD + 0.9151*N;
523 64000 : const double r2KM = - LOGD + 1.0800*N;
524 64000 : long Cmin = 3, Cmax = 3, i = 0;
525 574127 : while (tailres(R1, R2, al2K, rKm, rKM, r1Km, r1KM, r2Km, r2KM, Cmax, i) > lim)
526 : {
527 510127 : Cmin = Cmax;
528 510127 : Cmax *= 2;
529 510127 : i++;
530 : }
531 63994 : i--;
532 619195 : while (Cmax - Cmin > 1)
533 : {
534 555200 : long t = (Cmin + Cmax)/2;
535 555200 : if (tailres(R1, R2, al2K, rKm, rKM, r1Km, r1KM, r2Km, r2KM, t, i) > lim)
536 344300 : Cmin = t;
537 : else
538 210901 : Cmax = t;
539 : }
540 63995 : return Cmax;
541 : }
542 :
543 : /* ~ 1 / Res(s = 1, zeta_K) */
544 : static GEN
545 64001 : compute_invres(GRHcheck_t *S, long LIMC)
546 : {
547 64001 : pari_sp av = avma;
548 64001 : double loginvres = 0.;
549 : GRHprime_t *pr;
550 : long i;
551 64001 : double logLIMC = log((double)LIMC);
552 64001 : double logLIMC2 = logLIMC*logLIMC, denc;
553 : double c0, c1, c2;
554 64001 : denc = 1/(pow((double)LIMC, 3.) * logLIMC * logLIMC2);
555 64001 : c2 = ( logLIMC2 + 3 * logLIMC / 2 + 1) * denc;
556 64001 : denc *= LIMC;
557 64001 : c1 = (3 * logLIMC2 + 4 * logLIMC + 2) * denc;
558 64001 : denc *= LIMC;
559 64001 : c0 = (3 * logLIMC2 + 5 * logLIMC / 2 + 1) * denc;
560 7522275 : for (pr = S->primes, i = S->nprimes; i > 0; pr++, i--)
561 : {
562 : GEN dec, fs, ns;
563 : long addpsi;
564 : double addpsi1, addpsi2;
565 7514364 : double logp = pr->logp, NPk;
566 7514364 : long j, k, limp = logLIMC/logp;
567 7514364 : ulong p = pr->p, p2 = p*p;
568 7514364 : if (limp < 1) break;
569 7458274 : dec = pr->dec;
570 7458274 : fs = gel(dec, 1); ns = gel(dec, 2);
571 7458274 : loginvres += 1./p;
572 : /* NB: limp = 1 nearly always and limp > 2 for very few primes */
573 8824624 : for (k=2, NPk = p; k <= limp; k++) { NPk *= p; loginvres += 1/(k * NPk); }
574 7458274 : addpsi = limp;
575 7458274 : addpsi1 = p *(pow((double)p , (double)limp)-1)/(p -1);
576 7458274 : addpsi2 = p2*(pow((double)p2, (double)limp)-1)/(p2-1);
577 7458274 : j = lg(fs);
578 16630148 : while (--j > 0)
579 : {
580 : long f, nb, kmax;
581 : double NP, NP2, addinvres;
582 9171874 : f = fs[j]; if (f > limp) continue;
583 3983460 : nb = ns[j];
584 3983460 : NP = pow((double)p, (double)f);
585 3983460 : addinvres = 1/NP;
586 3983460 : kmax = limp / f;
587 4860310 : for (k=2, NPk = NP; k <= kmax; k++) { NPk *= NP; addinvres += 1/(k*NPk); }
588 3983460 : NP2 = NP*NP;
589 3983460 : loginvres -= nb * addinvres;
590 3983460 : addpsi -= nb * f * kmax;
591 3983460 : addpsi1 -= nb*(f*NP *(pow(NP ,(double)kmax)-1)/(NP -1));
592 3983460 : addpsi2 -= nb*(f*NP2*(pow(NP2,(double)kmax)-1)/(NP2-1));
593 : }
594 7458274 : loginvres -= (addpsi*c0 - addpsi1*c1 + addpsi2*c2)*logp;
595 : }
596 64001 : return gc_leaf(av, mpexp(dbltor(loginvres)));
597 : }
598 :
599 : static long
600 63999 : nthideal(GRHcheck_t *S, GEN nf, long n)
601 : {
602 63999 : pari_sp av = avma;
603 63999 : GEN P = nf_get_pol(nf);
604 63999 : ulong p = 0, *vecN = (ulong*)const_vecsmall(n, LONG_MAX);
605 63999 : long i, N = poldegree(P, -1);
606 63999 : for (i = 0; ; i++)
607 230578 : {
608 : GRHprime_t *pr;
609 : GEN fs;
610 294577 : cache_prime_dec(S, p+1, nf);
611 294578 : pr = S->primes + i;
612 294578 : fs = gel(pr->dec, 1);
613 294578 : p = pr->p;
614 294578 : if (fs[1] != N)
615 : {
616 198011 : GEN ns = gel(pr->dec, 2);
617 198011 : long k, l, j = lg(fs);
618 443720 : while (--j > 0)
619 : {
620 245709 : ulong NP = upowuu(p, fs[j]);
621 : long nf;
622 245709 : if (!NP) continue;
623 754375 : for (k = 1; k <= n; k++) if (vecN[k] > NP) break;
624 245317 : if (k > n) continue;
625 : /* vecN[k] <= NP */
626 158866 : nf = ns[j]; /*#{primes of norme NP} = nf, insert them here*/
627 355595 : for (l = k+nf; l <= n; l++) vecN[l] = vecN[l-nf];
628 401542 : for (l = 0; l < nf && k+l <= n; l++) vecN[k+l] = NP;
629 365199 : while (l <= k) vecN[l++] = NP;
630 : }
631 : }
632 294578 : if (p > vecN[n]) break;
633 : }
634 64000 : return gc_long(av, vecN[n]);
635 : }
636 :
637 : /* volume of unit ball in R^n: \pi^{n/2} / \Gamma(n/2 + 1) */
638 : static double
639 66247 : ballvol(long n)
640 : {
641 66247 : double v = odd(n)? 2: 1;
642 151516 : for (; n > 1; n -= 2) v *= (2 * M_PI) / n;
643 66247 : return v;
644 : }
645 :
646 : /* Compute FB, LV, iLP + KC*. Reset perm
647 : * C2: bound for norm of tested prime ideals (includes be_honest())
648 : * C1: bound for p, such that P|p (NP <= C2) used to build relations */
649 : static void
650 66250 : FBgen(FB_t *F, GEN nf, long N, ulong C1, ulong C2, GRHcheck_t *S)
651 : {
652 : GRHprime_t *pr;
653 : long i, ip;
654 : GEN prim;
655 66250 : const double L = log((double)C2 + 0.5);
656 :
657 66250 : cache_prime_dec(S, C2, nf);
658 66250 : pr = S->primes;
659 66250 : F->sfb_chg = 0;
660 66250 : F->FB = cgetg(C2+1, t_VECSMALL);
661 66250 : F->iLP = cgetg(C2+1, t_VECSMALL);
662 66250 : F->LV = zerovec(C2);
663 :
664 66250 : prim = icopy(gen_1);
665 66252 : i = ip = 0;
666 66252 : F->KC = F->KCZ = 0;
667 438150 : for (;; pr++) /* p <= C2 */
668 438150 : {
669 504402 : ulong p = pr->p;
670 : long k, l, m;
671 : GEN LP, nb, f;
672 :
673 504402 : if (!F->KC && p > C1) { F->KCZ = i; F->KC = ip; }
674 504402 : if (p > C2) break;
675 :
676 467194 : if (DEBUGLEVEL>1) err_printf(" %ld",p);
677 :
678 467195 : f = gel(pr->dec, 1); nb = gel(pr->dec, 2);
679 467195 : if (f[1] == N)
680 : {
681 146468 : if (p == C2) break;
682 137928 : continue; /* p inert */
683 : }
684 320727 : l = (long)(L/pr->logp); /* p^f <= C2 <=> f <= l */
685 584608 : for (k=0, m=1; m < lg(f) && f[m]<=l; m++) k += nb[m];
686 320727 : if (!k)
687 : { /* too inert to appear in FB */
688 73555 : if (p == C2) break;
689 72701 : continue;
690 : }
691 247172 : prim[2] = p; LP = idealprimedec_limit_f(nf,prim, l);
692 : /* keep noninert ideals with Norm <= C2 */
693 247170 : if (m == lg(f)) setisclone(LP); /* flag it: all prime divisors in FB */
694 247170 : F->FB[++i]= p;
695 247170 : gel(F->LV,p) = LP;
696 247170 : F->iLP[p] = ip; ip += k;
697 247170 : if (p == C2) break;
698 : }
699 66251 : if (!F->KC) { F->KCZ = i; F->KC = ip; }
700 : /* Note F->KC > 0 otherwise GRHchk is false */
701 66251 : setlg(F->FB, F->KCZ+1); F->KCZ2 = i;
702 66250 : F->prodZ = zv_prod_Z(F->FB);
703 66248 : if (DEBUGLEVEL>1)
704 : {
705 0 : err_printf("\n");
706 0 : if (DEBUGLEVEL>6)
707 : {
708 0 : err_printf("########## FACTORBASE ##########\n\n");
709 0 : err_printf("KC2=%ld, KC=%ld, KCZ=%ld, KCZ2=%ld\n",
710 : ip, F->KC, F->KCZ, F->KCZ2);
711 0 : for (i=1; i<=F->KCZ; i++) err_printf("++ LV[%ld] = %Ps",i,gel(F->LV,F->FB[i]));
712 : }
713 : }
714 66248 : F->perm = NULL; F->L_jid = NULL;
715 66248 : F->ballvol = ballvol(nf_get_degree(nf));
716 66247 : }
717 :
718 : static int
719 496385 : GRHchk(GEN nf, GRHcheck_t *S, ulong LIMC)
720 : {
721 496385 : double logC = log((double)LIMC), SA = 0, SB = 0;
722 496385 : GRHprime_t *pr = S->primes;
723 :
724 496385 : cache_prime_dec(S, LIMC, nf);
725 496384 : for (pr = S->primes;; pr++)
726 3084591 : {
727 3580975 : ulong p = pr->p;
728 : GEN dec, fs, ns;
729 : double logCslogp;
730 : long j;
731 :
732 3580975 : if (p > LIMC) break;
733 3190972 : dec = pr->dec; fs = gel(dec, 1); ns = gel(dec,2);
734 3190972 : logCslogp = logC/pr->logp;
735 5022085 : for (j = 1; j < lg(fs); j++)
736 : {
737 3929194 : long f = fs[j], M, nb;
738 : double logNP, q, A, B;
739 3929194 : if (f > logCslogp) break;
740 1831114 : logNP = f * pr->logp;
741 1831114 : q = 1/sqrt((double)upowuu(p, f));
742 1831113 : A = logNP * q; B = logNP * A; M = (long)(logCslogp/f);
743 1831113 : if (M > 1)
744 : {
745 376640 : double inv1_q = 1 / (1-q);
746 376640 : A *= (1 - pow(q, (double)M)) * inv1_q;
747 376640 : B *= (1 - pow(q, (double)M)*(M+1 - M*q)) * inv1_q * inv1_q;
748 : }
749 1831113 : nb = ns[j];
750 1831113 : SA += nb * A;
751 1831113 : SB += nb * B;
752 : }
753 3190971 : if (p == LIMC) break;
754 : }
755 496383 : return GRHok(S, logC, SA, SB);
756 : }
757 :
758 : /* SMOOTH IDEALS */
759 : static void
760 9461780 : store(long i, long e, FACT *fact)
761 : {
762 9461780 : ++fact[0].pr;
763 9461780 : fact[fact[0].pr].pr = i; /* index */
764 9461780 : fact[fact[0].pr].ex = e; /* exponent */
765 9461780 : }
766 :
767 : /* divide out m by all P|p, k = v_p(Nm) */
768 : static int
769 2343 : divide_p_elt(GEN LP, long ip, long k, GEN m, FACT *fact)
770 : {
771 2343 : long j, l = lg(LP);
772 3188 : for (j=1; j<l; j++)
773 : {
774 3188 : GEN P = gel(LP,j);
775 3188 : long v = ZC_nfval(m, P);
776 3188 : if (!v) continue;
777 2771 : store(ip + j, v, fact); /* v = v_P(m) > 0 */
778 2771 : k -= v * pr_get_f(P);
779 2771 : if (!k) return 1;
780 : }
781 0 : return 0;
782 : }
783 : /* divide out I by all P|p, k = v_p(NI) */
784 : static int
785 185478 : divide_p_id(GEN LP, long ip, long k, GEN nf, GEN I, FACT *fact)
786 : {
787 185478 : long j, l = lg(LP);
788 278162 : for (j=1; j<l; j++)
789 : {
790 270301 : GEN P = gel(LP,j);
791 270301 : long v = idealval(nf,I, P);
792 270300 : if (!v) continue;
793 181115 : store(ip + j, v, fact); /* v = v_P(I) > 0 */
794 181117 : k -= v * pr_get_f(P);
795 181117 : if (!k) return 1;
796 : }
797 7861 : return 0;
798 : }
799 : /* divide out m/I by all P|p, k = v_p(Nm/NI) */
800 : static int
801 5436028 : divide_p_quo(GEN LP, long ip, long k, GEN nf, GEN I, GEN m, FACT *fact)
802 : {
803 5436028 : long j, l = lg(LP);
804 17749507 : for (j=1; j<l; j++)
805 : {
806 17660797 : GEN P = gel(LP,j);
807 17660797 : long v = ZC_nfval(m, P);
808 17659669 : if (!v) continue;
809 8106003 : v -= idealval(nf,I, P);
810 8107344 : if (!v) continue;
811 6887672 : store(ip + j, v, fact); /* v = v_P(m / I) > 0 */
812 6887707 : k -= v * pr_get_f(P);
813 6887630 : if (!k) return 1;
814 : }
815 88710 : return 0;
816 : }
817 :
818 : static int
819 5623840 : divide_p(FB_t *F, long p, long k, GEN nf, GEN I, GEN m, FACT *fact)
820 : {
821 5623840 : GEN LP = gel(F->LV,p);
822 5623840 : long ip = F->iLP[p];
823 5623840 : if (!m) return divide_p_id (LP,ip,k,nf,I,fact);
824 5438362 : if (!I) return divide_p_elt(LP,ip,k,m,fact);
825 5436019 : return divide_p_quo(LP,ip,k,nf,I,m,fact);
826 : }
827 :
828 : /* Let x = m if I == NULL,
829 : * I if m == NULL,
830 : * m/I otherwise.
831 : * Can we factor the integral primitive ideal x ? |N| = Norm x > 0 */
832 : static long
833 18355117 : can_factor(FB_t *F, GEN nf, GEN I, GEN m, GEN N, FACT *fact)
834 : {
835 : GEN f, p, e;
836 : long i, l;
837 18355117 : fact[0].pr = 0;
838 18355117 : if (is_pm1(N)) return 1;
839 17376301 : if (!is_pm1(Z_ppo(N, F->prodZ))) return 0;
840 2841684 : f = absZ_factor(N); p = gel(f,1); e = gel(f,2); l = lg(p);
841 8369311 : for (i = 1; i < l; i++)
842 5623828 : if (!divide_p(F, itou(gel(p,i)), itou(gel(e,i)), nf, I, m, fact))
843 : {
844 96304 : if (DEBUGLEVEL > 1) err_printf(".");
845 96304 : return 0;
846 : }
847 2745483 : return 1;
848 : }
849 :
850 : /* can we factor m/I ? [m in I from idealpseudomin_nonscalar], NI = norm I */
851 : static long
852 16942547 : factorgen(FB_t *F, GEN nf, GEN I, GEN NI, GEN m, FACT *fact)
853 : {
854 : long e;
855 16942547 : GEN Nm = embed_norm(RgM_RgC_mul(nf_get_M(nf),m), nf_get_r1(nf));
856 16942600 : GEN N = grndtoi(NI? divri(Nm, NI): Nm, &e); /* ~ N(m/I) */
857 16942531 : if (e > -32)
858 : {
859 0 : if (DEBUGLEVEL > 1) err_printf("+");
860 0 : return 0;
861 : }
862 16942531 : return can_factor(F, nf, I, m, N, fact);
863 : }
864 :
865 : /* FUNDAMENTAL UNITS */
866 :
867 : /* a, y real. Return (Re(x) + a) + I * (Im(x) % y) */
868 : static GEN
869 6906846 : addRe_modIm(GEN x, GEN a, GEN y, GEN iy)
870 : {
871 : GEN z;
872 6906846 : if (typ(x) == t_COMPLEX)
873 : {
874 4906342 : GEN re, im = modRr_i(gel(x,2), y, iy);
875 4906303 : if (!im) return NULL;
876 4906303 : re = gadd(gel(x,1), a);
877 4906280 : z = gequal0(im)? re: mkcomplex(re, im);
878 : }
879 : else
880 2000504 : z = gadd(x, a);
881 6906779 : return z;
882 : }
883 : static GEN
884 204620 : modIm(GEN x, GEN y, GEN iy)
885 : {
886 204620 : if (typ(x) == t_COMPLEX)
887 : {
888 191445 : GEN im = modRr_i(gel(x,2), y, iy);
889 191442 : if (!im) return NULL;
890 191442 : x = gequal0(im)? gel(x,1): mkcomplex(gel(x,1), im);
891 : }
892 204617 : return x;
893 : }
894 :
895 : /* clean archimedean components. ipi = 2^n / pi (n arbitrary); its
896 : * exponent may be modified */
897 : static GEN
898 3071479 : cleanarch(GEN x, long N, GEN ipi, long prec)
899 : {
900 : long i, l, R1, RU;
901 3071479 : GEN s, y = cgetg_copy(x, &l);
902 :
903 3071476 : if (!ipi) ipi = invr(mppi(prec));
904 3071474 : if (typ(x) == t_MAT)
905 : {
906 529385 : for (i = 1; i < l; i++)
907 465235 : if (!(gel(y,i) = cleanarch(gel(x,i), N, ipi, prec))) return NULL;
908 64150 : return y;
909 : }
910 3007320 : RU = l-1; R1 = (RU<<1) - N;
911 3007320 : s = gdivgs(RgV_sum(real_i(x)), -N); /* -log |norm(x)| / N */
912 3007295 : i = 1;
913 3007295 : if (R1)
914 : {
915 2523880 : GEN pi2 = Pi2n(1, prec);
916 2523887 : setexpo(ipi, -3); /* 1/(2pi) */
917 7790474 : for (; i <= R1; i++)
918 5266616 : if (!(gel(y,i) = addRe_modIm(gel(x,i), s, pi2, ipi))) return NULL;
919 : }
920 3007273 : if (i <= RU)
921 : {
922 1080298 : GEN pi4 = Pi2n(2, prec), s2 = gmul2n(s, 1);
923 1080312 : setexpo(ipi, -4); /* 1/(4pi) */
924 2720580 : for (; i <= RU; i++)
925 1640238 : if (!(gel(y,i) = addRe_modIm(gel(x,i), s2, pi4, ipi))) return NULL;
926 : }
927 3007317 : return y;
928 : }
929 : GEN
930 197126 : nf_cxlog_normalize(GEN nf, GEN x, long prec)
931 : {
932 197126 : long N = nf_get_degree(nf);
933 197126 : return cleanarch(x, N, NULL, prec);
934 : }
935 :
936 : /* clean unit archimedean components. ipi = 2^n / pi (n arbitrary); its
937 : * exponent may be modified */
938 : static GEN
939 133868 : cleanarchunit(GEN x, long N, GEN ipi, long prec)
940 : {
941 : long i, l, R1, RU;
942 133868 : GEN y = cgetg_copy(x, &l);
943 :
944 133868 : if (!ipi) ipi = invr(mppi(prec));
945 133868 : if (typ(x) == t_MAT)
946 : {
947 133861 : for (i = 1; i < l; i++)
948 69860 : if (!(gel(y,i) = cleanarchunit(gel(x,i), N, ipi, prec))) return NULL;
949 64001 : return y;
950 : }
951 69860 : if (gexpo(RgV_sum(real_i(x))) > -10) return NULL;
952 69852 : RU = l-1; R1 = (RU<<1) - N;
953 69852 : i = 1;
954 69852 : if (R1)
955 : {
956 55376 : GEN pi2 = Pi2n(1, prec);
957 55376 : setexpo(ipi, -3); /* 1/(2pi) */
958 189321 : for (; i <= R1; i++)
959 133949 : if (!(gel(y,i) = modIm(gel(x,i), pi2, ipi))) return NULL;
960 : }
961 69848 : if (i <= RU)
962 : {
963 34447 : GEN pi4 = Pi2n(2, prec);
964 34447 : setexpo(ipi, -4); /* 1/(4pi) */
965 105124 : for (; i <= RU; i++)
966 70672 : if (!(gel(y,i) = modIm(gel(x,i), pi4, ipi))) return NULL;
967 : }
968 69853 : return y;
969 : }
970 :
971 : static GEN
972 396 : not_given(long reason)
973 : {
974 396 : if (DEBUGLEVEL)
975 0 : switch(reason)
976 : {
977 0 : case fupb_LARGE:
978 0 : pari_warn(warner,"fundamental units too large, not given");
979 0 : break;
980 0 : case fupb_PRECI:
981 0 : pari_warn(warner,"insufficient precision for fundamental units, not given");
982 0 : break;
983 : }
984 396 : return NULL;
985 : }
986 :
987 : /* check whether exp(x) will 1) get too big (real(x) large), 2) require
988 : * large accuracy for argument reduction (imag(x) large) */
989 : static long
990 2840326 : expbitprec(GEN x, long *e)
991 : {
992 : GEN re, im;
993 2840326 : if (typ(x) != t_COMPLEX) re = x;
994 : else
995 : {
996 1778704 : im = gel(x,2); *e = maxss(*e, expo(im) + 5 - bit_prec(im));
997 1778704 : re = gel(x,1);
998 : }
999 2840326 : return (expo(re) <= 20);
1000 :
1001 : }
1002 : static long
1003 1237088 : RgC_expbitprec(GEN x)
1004 : {
1005 1237088 : long l = lg(x), i, e = - (long)HIGHEXPOBIT;
1006 3872328 : for (i = 1; i < l; i++)
1007 2635814 : if (!expbitprec(gel(x,i), &e)) return LONG_MAX;
1008 1236514 : return e;
1009 : }
1010 : static long
1011 48790 : RgM_expbitprec(GEN x)
1012 : {
1013 48790 : long i, j, I, J, e = - (long)HIGHEXPOBIT;
1014 48790 : RgM_dimensions(x, &I,&J);
1015 118573 : for (j = 1; j <= J; j++)
1016 274295 : for (i = 1; i <= I; i++)
1017 204512 : if (!expbitprec(gcoeff(x,i,j), &e)) return LONG_MAX;
1018 48727 : return e;
1019 : }
1020 :
1021 : static GEN
1022 1379 : FlxqX_chinese_unit(GEN X, GEN U, GEN invzk, GEN D, GEN T, ulong p)
1023 : {
1024 1379 : long i, lU = lg(U), lX = lg(X), d = lg(invzk)-1;
1025 1379 : GEN M = cgetg(lU, t_MAT);
1026 1379 : if (D)
1027 : {
1028 1272 : D = Flv_inv(D, p);
1029 69740 : for (i = 1; i < lX; i++)
1030 68468 : if (uel(D, i) != 1)
1031 56623 : gel(X,i) = Flx_Fl_mul(gel(X,i), uel(D,i), p);
1032 : }
1033 3878 : for (i = 1; i < lU; i++)
1034 : {
1035 2499 : GEN H = FlxqV_factorback(X, gel(U, i), T, p);
1036 2499 : gel(M, i) = Flm_Flc_mul(invzk, Flx_to_Flv(H, d), p);
1037 : }
1038 1379 : return M;
1039 : }
1040 :
1041 : static GEN
1042 274 : chinese_unit_slice(GEN A, GEN U, GEN B, GEN D, GEN C, GEN P, GEN *mod)
1043 : {
1044 274 : pari_sp av = avma;
1045 274 : long i, n = lg(P)-1, v = varn(C);
1046 : GEN H, T;
1047 274 : if (n == 1)
1048 : {
1049 0 : ulong p = uel(P,1);
1050 0 : GEN a = ZXV_to_FlxV(A, p), b = ZM_to_Flm(B, p), c = ZX_to_Flx(C, p);
1051 0 : GEN d = D ? ZV_to_Flv(D, p): NULL;
1052 0 : GEN Hp = FlxqX_chinese_unit(a, U, b, d, c, p);
1053 0 : H = gc_upto(av, Flm_to_ZM(Hp));
1054 0 : *mod = utoi(p);
1055 0 : return H;
1056 : }
1057 274 : T = ZV_producttree(P);
1058 274 : A = ZXC_nv_mod_tree(A, P, T, v);
1059 274 : B = ZM_nv_mod_tree(B, P, T);
1060 274 : D = D ? ZV_nv_mod_tree(D, P, T): NULL;
1061 274 : C = ZX_nv_mod_tree(C, P, T);
1062 :
1063 274 : H = cgetg(n+1, t_VEC);
1064 1653 : for(i=1; i <= n; i++)
1065 : {
1066 1379 : ulong p = P[i];
1067 1379 : GEN a = gel(A,i), b = gel(B,i), c = gel(C,i), d = D ? gel(D,i): NULL;
1068 1379 : gel(H,i) = FlxqX_chinese_unit(a, U, b, d, c, p);
1069 : }
1070 274 : H = nmV_chinese_center_tree_seq(H, P, T, ZV_chinesetree(P, T));
1071 274 : *mod = gmael(T, lg(T)-1, 1); return gc_all(av, 2, &H, mod);
1072 : }
1073 :
1074 : GEN
1075 274 : chinese_unit_worker(GEN P, GEN A, GEN U, GEN B, GEN D, GEN C)
1076 : {
1077 274 : GEN V = cgetg(3, t_VEC);
1078 274 : gel(V,1) = chinese_unit_slice(A, U, B, isintzero(D) ? NULL: D, C, P, &gel(V,2));
1079 274 : return V;
1080 : }
1081 :
1082 : /* Let x = \prod X[i]^E[i] = u, return u.
1083 : * If dX != NULL, X[i] = nX[i] / dX[i] where nX[i] is a ZX, dX[i] in Z */
1084 : static GEN
1085 94 : chinese_unit(GEN nf, GEN nX, GEN dX, GEN U, ulong bnd)
1086 : {
1087 94 : pari_sp av = avma;
1088 94 : GEN f = nf_get_index(nf), T = nf_get_pol(nf), invzk = nf_get_invzk(nf);
1089 : GEN H, mod;
1090 : forprime_t S;
1091 94 : GEN worker = snm_closure(is_entry("_chinese_unit_worker"),
1092 : mkcol5(nX, U, invzk, dX? dX: gen_0, T));
1093 94 : init_modular_big(&S);
1094 94 : H = gen_crt("chinese_units", worker, &S, f, bnd, 0, &mod, nmV_chinese_center, FpM_center);
1095 94 : settyp(H, t_VEC); return gc_GEN(av, H);
1096 : }
1097 :
1098 : /* *pE a ZM */
1099 : static void
1100 164 : ZM_remove_unused(GEN *pE, GEN *pX)
1101 : {
1102 164 : long j, k, l = lg(*pX);
1103 164 : GEN E = *pE, v = cgetg(l, t_VECSMALL);
1104 16395 : for (j = k = 1; j < l; j++)
1105 16231 : if (!ZMrow_equal0(E, j)) v[k++] = j;
1106 164 : if (k < l)
1107 : {
1108 164 : setlg(v, k);
1109 164 : *pX = vecpermute(*pX,v);
1110 164 : *pE = rowpermute(E,v);
1111 : }
1112 164 : }
1113 :
1114 : /* s = -log|norm(x)|/N */
1115 : static GEN
1116 1306928 : fixarch(GEN x, GEN s, long R1)
1117 : {
1118 : long i, l;
1119 1306928 : GEN y = cgetg_copy(x, &l);
1120 3644835 : for (i = 1; i <= R1; i++) gel(y,i) = gadd(s, gel(x,i));
1121 1810044 : for ( ; i < l; i++) gel(y,i) = gadd(s, gmul2n(gel(x,i),-1));
1122 1306940 : return y;
1123 : }
1124 :
1125 : static GEN
1126 64001 : getfu(GEN nf, GEN *ptA, GEN *ptU, long prec)
1127 : {
1128 64001 : GEN U, y, matep, A, T = nf_get_pol(nf), M = nf_get_M(nf);
1129 64001 : long e, j, R1, RU, N = degpol(T);
1130 :
1131 64001 : R1 = nf_get_r1(nf); RU = (N+R1) >> 1;
1132 64001 : if (RU == 1) return cgetg(1,t_VEC);
1133 :
1134 48790 : A = *ptA;
1135 48790 : matep = cgetg(RU,t_MAT);
1136 118642 : for (j = 1; j < RU; j++)
1137 : {
1138 69853 : GEN Aj = gel(A,j), s = gdivgs(RgV_sum(real_i(Aj)), -N);
1139 69850 : gel(matep,j) = fixarch(Aj, s, R1);
1140 : }
1141 48789 : U = lll(real_i(matep));
1142 48790 : if (lg(U) < RU) return not_given(fupb_PRECI);
1143 48790 : if (ptU) { *ptU = U; *ptA = A = RgM_ZM_mul(A,U); }
1144 48790 : y = RgM_ZM_mul(matep,U);
1145 48790 : e = RgM_expbitprec(y);
1146 48790 : if (e >= 0) return not_given(e == LONG_MAX? fupb_LARGE: fupb_PRECI);
1147 48727 : if (prec <= 0) prec = gprecision(A);
1148 48727 : y = RgM_solve_realimag(M, gexp(y,prec));
1149 48727 : if (!y) return not_given(fupb_PRECI);
1150 48727 : y = grndtoi(y, &e); if (e >= 0) return not_given(fupb_PRECI);
1151 48400 : settyp(y, t_VEC);
1152 :
1153 48400 : if (!ptU) *ptA = A = RgM_ZM_mul(A, U);
1154 117505 : for (j = 1; j < RU; j++)
1155 : { /* y[i] are hopefully unit generators. Normalize: smallest T2 norm */
1156 69111 : GEN u = gel(y,j), v = zk_inv(nf, u);
1157 69113 : if (!v || !is_pm1(Q_denom(v)) || ZV_isscalar(u))
1158 8 : return not_given(fupb_PRECI);
1159 69105 : if (gcmp(RgC_fpnorml2(v,DEFAULTPREC), RgC_fpnorml2(u,DEFAULTPREC)) < 0)
1160 : {
1161 29885 : gel(A,j) = RgC_neg(gel(A,j));
1162 29885 : if (ptU) gel(U,j) = ZC_neg(gel(U,j));
1163 29885 : u = v;
1164 : }
1165 69105 : gel(y,j) = nf_to_scalar_or_alg(nf, u);
1166 : }
1167 48394 : return y;
1168 : }
1169 :
1170 : static void
1171 0 : err_units() { pari_err_PREC("makeunits [cannot get units, use bnfinit(,1)]"); }
1172 :
1173 : /* bound for log2 |sigma(u)|, sigma complex embedding, u fundamental unit
1174 : * attached to bnf_get_logfu */
1175 : static double
1176 94 : log2fubound(GEN bnf)
1177 : {
1178 94 : GEN LU = bnf_get_logfu(bnf);
1179 94 : long i, j, l = lg(LU), r1 = nf_get_r1(bnf_get_nf(bnf));
1180 94 : double e = 0.0;
1181 330 : for (j = 1; j < l; j++)
1182 : {
1183 236 : GEN u = gel(LU,j);
1184 624 : for (i = 1; i <= r1; i++)
1185 : {
1186 388 : GEN E = real_i(gel(u,i));
1187 388 : e = maxdd(e, gtodouble(E));
1188 : }
1189 842 : for ( ; i <= l; i++)
1190 : {
1191 606 : GEN E = real_i(gel(u,i));
1192 606 : e = maxdd(e, gtodouble(E) / 2);
1193 : }
1194 : }
1195 94 : return e / M_LN2;
1196 : }
1197 : /* bound for log2(|RgM_solve_realimag(M, y)|_oo / |y|_oo)*/
1198 : static double
1199 94 : log2Mbound(GEN nf)
1200 : {
1201 94 : GEN G = nf_get_G(nf), D = nf_get_disc(nf);
1202 94 : long r2 = nf_get_r2(nf), l = lg(G), i;
1203 94 : double e, d = dbllog2(D)/2 - r2 * M_LN2; /* log2 |det(split_realimag(M))| */
1204 94 : e = log2(nf_get_degree(nf));
1205 535 : for (i = 2; i < l; i++) e += dbllog2(gnorml2(gel(G,i))); /* Hadamard bound */
1206 94 : return e / 2 - d;
1207 : }
1208 :
1209 : static GEN
1210 94 : vec_chinese_units(GEN bnf)
1211 : {
1212 94 : GEN nf = bnf_get_nf(bnf), SUnits = bnf_get_sunits(bnf);
1213 94 : double bnd = ceil(log2Mbound(nf) + log2fubound(bnf));
1214 94 : GEN X, dX, Y, U, f = nf_get_index(nf);
1215 94 : long j, l, v = nf_get_varn(nf);
1216 94 : if (!SUnits) err_units(); /* no compact units */
1217 94 : Y = gel(SUnits,1);
1218 94 : U = gel(SUnits,2);
1219 94 : ZM_remove_unused(&U, &Y); l = lg(Y); X = cgetg(l, t_VEC);
1220 94 : if (is_pm1(f)) f = dX = NULL; else dX = cgetg(l, t_VEC);
1221 5142 : for (j = 1; j < l; j++)
1222 : {
1223 5048 : GEN t = nf_to_scalar_or_alg(nf, gel(Y,j));
1224 5048 : if (f)
1225 : {
1226 : GEN den;
1227 4202 : t = Q_remove_denom(t, &den);
1228 4202 : gel(dX,j) = den ? den: gen_1;
1229 : }
1230 5048 : gel(X,j) = typ(t) == t_INT? scalarpol_shallow(t,v): t;
1231 : }
1232 94 : if (dblexpo(bnd) >= BITS_IN_LONG)
1233 0 : pari_err_OVERFLOW("vec_chinese_units [units too large]");
1234 94 : return chinese_unit(nf, X, dX, U, (ulong)bnd);
1235 : }
1236 :
1237 : static GEN
1238 24914 : makeunits(GEN bnf)
1239 : {
1240 24914 : GEN nf = bnf_get_nf(bnf), fu = bnf_get_fu_nocheck(bnf);
1241 24914 : GEN tu = nf_to_scalar_or_basis(nf, bnf_get_tuU(bnf));
1242 24914 : fu = (typ(fu) == t_MAT)? vec_chinese_units(bnf): matalgtobasis(nf, fu);
1243 24915 : return vec_prepend(fu, tu);
1244 : }
1245 :
1246 : /*******************************************************************/
1247 : /* */
1248 : /* PRINCIPAL IDEAL ALGORITHM (DISCRETE LOG) */
1249 : /* */
1250 : /*******************************************************************/
1251 :
1252 : /* G: prime ideals, E: vector of nonnegative exponents.
1253 : * C = possible extra prime (^1) or NULL
1254 : * Return Norm (product) */
1255 : static GEN
1256 69 : get_norm_fact_primes(GEN G, GEN E, GEN C)
1257 : {
1258 69 : pari_sp av=avma;
1259 69 : GEN N = gen_1, P, p;
1260 69 : long i, c = lg(E);
1261 69 : for (i=1; i<c; i++)
1262 : {
1263 0 : GEN ex = gel(E,i);
1264 0 : long s = signe(ex);
1265 0 : if (!s) continue;
1266 :
1267 0 : P = gel(G,i); p = pr_get_p(P);
1268 0 : N = mulii(N, powii(p, mului(pr_get_f(P), ex)));
1269 : }
1270 69 : if (C) N = mulii(N, pr_norm(C));
1271 69 : return gc_INT(av, N);
1272 : }
1273 :
1274 : /* gen: HNF ideals */
1275 : static GEN
1276 1231475 : get_norm_fact(GEN gen, GEN ex, GEN *pd)
1277 : {
1278 1231475 : long i, c = lg(ex);
1279 : GEN d,N,I,e,n,ne,de;
1280 1231475 : d = N = gen_1;
1281 1528233 : for (i=1; i<c; i++)
1282 296758 : if (signe(gel(ex,i)))
1283 : {
1284 175574 : I = gel(gen,i); e = gel(ex,i); n = ZM_det_triangular(I);
1285 175574 : ne = powii(n,e);
1286 175574 : de = equalii(n, gcoeff(I,1,1))? ne: powii(gcoeff(I,1,1), e);
1287 175574 : N = mulii(N, ne);
1288 175574 : d = mulii(d, de);
1289 : }
1290 1231475 : *pd = d; return N;
1291 : }
1292 :
1293 : static GEN
1294 1392614 : get_pr_lists(GEN FB, long N, int list_pr)
1295 : {
1296 : GEN pr, L;
1297 1392614 : long i, l = lg(FB), p, pmax;
1298 :
1299 1392614 : pmax = 0;
1300 9558023 : for (i=1; i<l; i++)
1301 : {
1302 8165409 : pr = gel(FB,i); p = pr_get_smallp(pr);
1303 8165409 : if (p > pmax) pmax = p;
1304 : }
1305 1392614 : L = const_vec(pmax, NULL);
1306 1392614 : if (list_pr)
1307 : {
1308 0 : for (i=1; i<l; i++)
1309 : {
1310 0 : pr = gel(FB,i); p = pr_get_smallp(pr);
1311 0 : if (!L[p]) gel(L,p) = vectrunc_init(N+1);
1312 0 : vectrunc_append(gel(L,p), pr);
1313 : }
1314 0 : for (p=1; p<=pmax; p++)
1315 0 : if (L[p]) gen_sort_inplace(gel(L,p), (void*)&cmp_prime_over_p,
1316 : &cmp_nodata, NULL);
1317 : }
1318 : else
1319 : {
1320 9558020 : for (i=1; i<l; i++)
1321 : {
1322 8165407 : pr = gel(FB,i); p = pr_get_smallp(pr);
1323 8165407 : if (!L[p]) gel(L,p) = vecsmalltrunc_init(N+1);
1324 8165406 : vecsmalltrunc_append(gel(L,p), i);
1325 : }
1326 : }
1327 1392613 : return L;
1328 : }
1329 :
1330 : /* recover FB, LV, iLP, KCZ from Vbase */
1331 : static GEN
1332 1392614 : recover_partFB(FB_t *F, GEN Vbase, long N)
1333 : {
1334 1392614 : GEN FB, LV, iLP, L = get_pr_lists(Vbase, N, 0);
1335 1392613 : long l = lg(L), p, ip, i;
1336 :
1337 1392613 : i = ip = 0;
1338 1392613 : FB = cgetg(l, t_VECSMALL);
1339 1392614 : iLP= cgetg(l, t_VECSMALL);
1340 1392612 : LV = cgetg(l, t_VEC);
1341 20446863 : for (p = 2; p < l; p++)
1342 : {
1343 19054249 : if (!L[p]) continue;
1344 4455878 : FB[++i] = p;
1345 4455878 : gel(LV,p) = vecpermute(Vbase, gel(L,p));
1346 4455879 : iLP[p]= ip; ip += lg(gel(L,p))-1;
1347 : }
1348 1392614 : F->KCZ = i;
1349 1392614 : F->KC = ip;
1350 1392614 : F->FB = FB; setlg(FB, i+1);
1351 1392612 : F->prodZ = zv_prod_Z(F->FB);
1352 1392612 : F->LV = LV;
1353 1392612 : F->iLP= iLP; return L;
1354 : }
1355 :
1356 : /* add v^e to factorization */
1357 : static void
1358 2919764 : add_to_fact(long v, long e, FACT *fact)
1359 : {
1360 2919764 : long i, n = fact[0].pr;
1361 9910470 : for (i=1; i<=n; i++)
1362 7520247 : if (fact[i].pr == v) { fact[i].ex += e; return; }
1363 2390223 : store(v, e, fact);
1364 : }
1365 : static void
1366 0 : inv_fact(FACT *fact)
1367 : {
1368 0 : long i, n = fact[0].pr;
1369 0 : for (i=1; i<=n; i++) fact[i].ex = -fact[i].ex;
1370 0 : }
1371 :
1372 : /* L (small) list of primes above the same p including pr. Return pr index */
1373 : static int
1374 3461 : pr_index(GEN L, GEN pr)
1375 : {
1376 3461 : long j, l = lg(L);
1377 3461 : GEN al = pr_get_gen(pr);
1378 3461 : for (j=1; j<l; j++)
1379 3461 : if (ZV_equal(al, pr_get_gen(gel(L,j)))) return j;
1380 0 : pari_err_BUG("codeprime");
1381 : return 0; /* LCOV_EXCL_LINE */
1382 : }
1383 :
1384 : static long
1385 3461 : Vbase_to_FB(FB_t *F, GEN pr)
1386 : {
1387 3461 : long p = pr_get_smallp(pr);
1388 3461 : return F->iLP[p] + pr_index(gel(F->LV,p), pr);
1389 : }
1390 :
1391 : /* x, y 2 extended ideals whose first component is an integral HNF and second
1392 : * a famat */
1393 : static GEN
1394 3537 : idealHNF_mulred(GEN nf, GEN x, GEN y)
1395 : {
1396 3537 : GEN A = idealHNF_mul(nf, gel(x,1), gel(y,1));
1397 3537 : GEN F = famat_mul_shallow(gel(x,2), gel(y,2));
1398 3537 : return idealred(nf, mkvec2(A, F));
1399 : }
1400 : /* idealred(x * pr^n), n > 0 is small, x extended ideal. Reduction in order to
1401 : * avoid prec pb: don't let id become too large as lgsub increases */
1402 : static GEN
1403 4663 : idealmulpowprime2(GEN nf, GEN x, GEN pr, ulong n)
1404 : {
1405 4663 : GEN A = idealmulpowprime(nf, gel(x,1), pr, utoipos(n));
1406 4663 : return mkvec2(A, gel(x,2));
1407 : }
1408 : static GEN
1409 65805 : init_famat(GEN x) { return mkvec2(x, trivial_fact()); }
1410 : /* optimized idealfactorback + reduction; z = init_famat() */
1411 : static GEN
1412 28784 : genback(GEN z, GEN nf, GEN P, GEN E)
1413 : {
1414 28784 : long i, l = lg(E);
1415 28784 : GEN I = NULL;
1416 76689 : for (i = 1; i < l; i++)
1417 47905 : if (signe(gel(E,i)))
1418 : {
1419 : GEN J;
1420 32321 : gel(z,1) = gel(P,i);
1421 32321 : J = idealpowred(nf, z, gel(E,i));
1422 32321 : I = I? idealHNF_mulred(nf, I, J): J;
1423 : }
1424 28784 : return I; /* != NULL since a generator */
1425 : }
1426 :
1427 : static GEN
1428 1252848 : SPLIT_i(FB_t *F, GEN nf, GEN G, GEN x, GEN xred, GEN Nx, FACT *fact)
1429 : {
1430 1252848 : pari_sp av = avma;
1431 1252848 : GEN L = idealpseudominvec(xred, G);
1432 1252849 : long k, l = lg(L);
1433 1336537 : for(k = 1; k < l; k++)
1434 1320449 : if (factorgen(F, nf, x, Nx, gel(L,k), fact)) return gel(L,k);
1435 16088 : return gc_NULL(av);
1436 : }
1437 : /* return famat y (principal ideal) such that y / x is smooth [wrt Vbase] */
1438 : static GEN
1439 1408956 : SPLIT(FB_t *F, GEN nf, GEN x, GEN Vbase, FACT *fact)
1440 : {
1441 1408956 : GEN vecG, ex, y, x0, Nx = ZM_det_triangular(x);
1442 : long nbtest_lim, nbtest, i, j, ru, lgsub;
1443 : pari_sp av;
1444 :
1445 : /* try without reduction if x is small */
1446 2817704 : if (expi(gcoeff(x,1,1)) < 100 &&
1447 1580945 : can_factor(F, nf, x, NULL, Nx, fact)) return NULL;
1448 1236759 : if ((y = SPLIT_i(F, nf, nf_get_roundG(nf), x, x, Nx, fact))) return y;
1449 :
1450 : /* reduce in various directions */
1451 8866 : ru = lg(nf_get_roots(nf));
1452 8866 : vecG = cgetg(ru, t_VEC);
1453 14471 : for (j=1; j<ru; j++)
1454 : {
1455 12667 : gel(vecG,j) = nf_get_Gtwist1(nf, j);
1456 12667 : if ((y = SPLIT_i(F, nf, gel(vecG,j), x, x, Nx, fact))) return y;
1457 : }
1458 :
1459 : /* tough case, multiply by random products */
1460 1804 : lgsub = 3; nbtest = 1; nbtest_lim = 4;
1461 1804 : ex = cgetg(lgsub, t_VECSMALL);
1462 1804 : x0 = init_famat(x);
1463 : for(;;)
1464 629 : {
1465 2433 : GEN Ired, I, NI, id = x0;
1466 2433 : av = avma;
1467 2433 : if (DEBUGLEVEL>2) err_printf("# ideals tried = %ld\n",nbtest);
1468 7411 : for (i=1; i<lgsub; i++)
1469 : {
1470 4978 : ex[i] = random_bits(RANDOM_BITS);
1471 4978 : if (ex[i]) id = idealmulpowprime2(nf, id, gel(Vbase,i), ex[i]);
1472 : }
1473 2433 : if (id == x0) continue;
1474 : /* I^(-1) * \prod Vbase[i]^ex[i] = (id[2]) / x */
1475 :
1476 2412 : I = gel(id,1); NI = ZM_det_triangular(I);
1477 2412 : if (can_factor(F, nf, I, NULL, NI, fact))
1478 : {
1479 0 : inv_fact(fact); /* I^(-1) */
1480 0 : for (i=1; i<lgsub; i++)
1481 0 : if (ex[i]) add_to_fact(Vbase_to_FB(F,gel(Vbase,i)), ex[i], fact);
1482 0 : return gel(id,2);
1483 : }
1484 2412 : Ired = ru == 2? I: ZM_lll(I, 0.99, LLL_INPLACE);
1485 4029 : for (j=1; j<ru; j++)
1486 3421 : if ((y = SPLIT_i(F, nf, gel(vecG,j), I, Ired, NI, fact)))
1487 : {
1488 5433 : for (i=1; i<lgsub; i++)
1489 3629 : if (ex[i]) add_to_fact(Vbase_to_FB(F,gel(Vbase,i)), ex[i], fact);
1490 1804 : return famat_mul_shallow(gel(id,2), y);
1491 : }
1492 608 : set_avma(av);
1493 608 : if (++nbtest > nbtest_lim)
1494 : {
1495 21 : nbtest = 0;
1496 21 : if (++lgsub < minss(8, lg(Vbase)-1))
1497 : {
1498 21 : nbtest_lim <<= 1;
1499 21 : ex = cgetg(lgsub, t_VECSMALL);
1500 : }
1501 0 : else nbtest_lim = LONG_MAX; /* don't increase further */
1502 21 : if (DEBUGLEVEL>2) err_printf("SPLIT: increasing factor base [%ld]\n",lgsub);
1503 : }
1504 : }
1505 : }
1506 :
1507 : INLINE GEN
1508 1397525 : bnf_get_W(GEN bnf) { return gel(bnf,1); }
1509 : INLINE GEN
1510 2785115 : bnf_get_B(GEN bnf) { return gel(bnf,2); }
1511 : INLINE GEN
1512 2819258 : bnf_get_C(GEN bnf) { return gel(bnf,4); }
1513 : INLINE GEN
1514 1392632 : bnf_get_vbase(GEN bnf) { return gel(bnf,5); }
1515 : INLINE GEN
1516 1392549 : bnf_get_Ur(GEN bnf) { return gmael(bnf,9,1); }
1517 : INLINE GEN
1518 271827 : bnf_get_ga(GEN bnf) { return gmael(bnf,9,2); }
1519 : INLINE GEN
1520 276783 : bnf_get_GD(GEN bnf) { return gmael(bnf,9,3); }
1521 :
1522 : /* Return y (as an elt of K or a t_MAT representing an elt in Z[K])
1523 : * such that x / (y) is smooth and store the exponents of its factorization
1524 : * on g_W and g_B in Wex / Bex; return NULL for y = 1 */
1525 : static GEN
1526 1392549 : split_ideal(GEN bnf, GEN x, GEN *pWex, GEN *pBex)
1527 : {
1528 1392549 : GEN L, y, Vbase = bnf_get_vbase(bnf);
1529 1392549 : GEN Wex, W = bnf_get_W(bnf);
1530 1392549 : GEN Bex, B = bnf_get_B(bnf);
1531 : long p, j, i, l, nW, nB;
1532 : FACT *fact;
1533 : FB_t F;
1534 :
1535 1392549 : L = recover_partFB(&F, Vbase, lg(x)-1);
1536 1392549 : fact = (FACT*)stack_malloc((F.KC+1)*sizeof(FACT));
1537 1392547 : y = SPLIT(&F, bnf_get_nf(bnf), x, Vbase, fact);
1538 1392549 : nW = lg(W)-1; *pWex = Wex = zero_zv(nW);
1539 1392549 : nB = lg(B)-1; *pBex = Bex = zero_zv(nB); l = lg(F.FB);
1540 1392549 : p = j = 0; /* -Wall */
1541 2092934 : for (i = 1; i <= fact[0].pr; i++)
1542 : { /* decode index C = ip+j --> (p,j) */
1543 700385 : long a, b, t, C = fact[i].pr;
1544 1953386 : for (t = 1; t < l; t++)
1545 : {
1546 1865910 : long q = F.FB[t], k = C - F.iLP[q];
1547 1865910 : if (k <= 0) break;
1548 1253001 : p = q;
1549 1253001 : j = k;
1550 : }
1551 700385 : a = gel(L, p)[j];
1552 700385 : b = a - nW;
1553 700385 : if (b <= 0) Wex[a] = y? -fact[i].ex: fact[i].ex;
1554 546764 : else Bex[b] = y? -fact[i].ex: fact[i].ex;
1555 : }
1556 1392549 : return y;
1557 : }
1558 :
1559 : GEN
1560 1107813 : init_red_mod_units(GEN bnf, long prec)
1561 : {
1562 1107813 : GEN s = gen_0, p1,s1,mat, logfu = bnf_get_logfu(bnf);
1563 1107813 : long i,j, RU = lg(logfu);
1564 :
1565 1107813 : if (RU == 1) return NULL;
1566 1107813 : mat = cgetg(RU,t_MAT);
1567 2507125 : for (j=1; j<RU; j++)
1568 : {
1569 1399314 : p1 = cgetg(RU+1,t_COL); gel(mat,j) = p1;
1570 1399314 : s1 = gen_0;
1571 3465874 : for (i=1; i<RU; i++)
1572 : {
1573 2066562 : gel(p1,i) = real_i(gcoeff(logfu,i,j));
1574 2066562 : s1 = mpadd(s1, mpsqr(gel(p1,i)));
1575 : }
1576 1399312 : gel(p1,RU) = gen_0; if (mpcmp(s1,s) > 0) s = s1;
1577 : }
1578 1107811 : s = gsqrt(gmul2n(s,RU),prec);
1579 1107813 : if (expo(s) < 27) s = utoipos(1UL << 27);
1580 1107811 : return mkvec2(mat, s);
1581 : }
1582 :
1583 : /* z computed above. Return unit exponents that would reduce col (arch) */
1584 : GEN
1585 1107810 : red_mod_units(GEN col, GEN z)
1586 : {
1587 : long i,RU;
1588 : GEN x,mat,N2;
1589 :
1590 1107810 : if (!z) return NULL;
1591 1107810 : mat= gel(z,1);
1592 1107810 : N2 = gel(z,2);
1593 1107810 : RU = lg(mat); x = cgetg(RU+1,t_COL);
1594 2507123 : for (i=1; i<RU; i++) gel(x,i) = real_i(gel(col,i));
1595 1107811 : gel(x,RU) = N2;
1596 1107811 : x = lll(shallowconcat(mat,x));
1597 1107813 : if (typ(x) != t_MAT || lg(x) <= RU) return NULL;
1598 1107813 : x = gel(x,RU);
1599 1107813 : if (signe(gel(x,RU)) < 0) x = gneg_i(x);
1600 1107813 : if (!gequal1(gel(x,RU))) pari_err_BUG("red_mod_units");
1601 1107813 : setlg(x,RU); return x;
1602 : }
1603 :
1604 : static GEN
1605 2250484 : add(GEN a, GEN t) { return a = a? RgC_add(a,t): t; }
1606 :
1607 : /* [x] archimedian components, A column vector. return [x] A */
1608 : static GEN
1609 2057279 : act_arch(GEN A, GEN x)
1610 : {
1611 : GEN a;
1612 2057279 : long i,l = lg(A), tA = typ(A);
1613 2057279 : if (tA == t_MAT)
1614 : { /* assume lg(x) >= l */
1615 192290 : a = cgetg(l, t_MAT);
1616 282152 : for (i=1; i<l; i++) gel(a,i) = act_arch(gel(A,i), x);
1617 192291 : return a;
1618 : }
1619 1864989 : if (l==1) return cgetg(1, t_COL);
1620 1864989 : a = NULL;
1621 1864989 : if (tA == t_VECSMALL)
1622 : {
1623 7202978 : for (i=1; i<l; i++)
1624 : {
1625 5971506 : long c = A[i];
1626 5971506 : if (c) a = add(a, gmulsg(c, gel(x,i)));
1627 : }
1628 : }
1629 : else
1630 : { /* A a t_COL of t_INT. Assume lg(A)==lg(x) */
1631 1384260 : for (i=1; i<l; i++)
1632 : {
1633 750746 : GEN c = gel(A,i);
1634 750746 : if (signe(c)) a = add(a, gmul(c, gel(x,i)));
1635 : }
1636 : }
1637 1864986 : return a? a: zerocol(lgcols(x)-1);
1638 : }
1639 : /* act_arch(matdiagonal(v), x) */
1640 : static GEN
1641 64096 : diagact_arch(GEN v, GEN x)
1642 : {
1643 64096 : long i, l = lg(v);
1644 64096 : GEN a = cgetg(l, t_MAT);
1645 92951 : for (i = 1; i < l; i++) gel(a,i) = gmul(gel(x,i), gel(v,i));
1646 64097 : return a;
1647 : }
1648 :
1649 : static long
1650 1410518 : prec_arch(GEN bnf)
1651 : {
1652 1410518 : GEN a = bnf_get_C(bnf);
1653 1410518 : long i, l = lg(a), prec;
1654 :
1655 1410518 : for (i=1; i<l; i++)
1656 1410434 : if ( (prec = gprecision(gel(a,i))) ) return prec;
1657 84 : return DEFAULTPREC;
1658 : }
1659 :
1660 : static long
1661 3857 : needed_bitprec(GEN x)
1662 : {
1663 3857 : long i, e = 0, l = lg(x);
1664 22552 : for (i = 1; i < l; i++)
1665 : {
1666 18695 : GEN c = gel(x,i);
1667 18695 : long f = gexpo(c) - gprecision(c);
1668 18695 : if (f > e) e = f;
1669 : }
1670 3857 : return e;
1671 : }
1672 :
1673 : /* col = archimedian components of x, Nx its norm, dx a multiple of its
1674 : * denominator. Return x or NULL (fail) */
1675 : GEN
1676 1237084 : isprincipalarch(GEN bnf, GEN col, GEN kNx, GEN e, GEN dx, long *pe)
1677 : {
1678 : GEN nf, x, y, logfu, s, M;
1679 1237084 : long N, prec = gprecision(col);
1680 1237086 : bnf = checkbnf(bnf); nf = bnf_get_nf(bnf); M = nf_get_M(nf);
1681 1237086 : if (!prec) prec = prec_arch(bnf);
1682 1237086 : *pe = 128;
1683 1237086 : logfu = bnf_get_logfu(bnf);
1684 1237086 : N = nf_get_degree(nf);
1685 1237086 : if (!(col = cleanarch(col,N,NULL,prec))) return NULL;
1686 1237087 : if (lg(col) > 2)
1687 : { /* reduce mod units */
1688 1107813 : GEN u, z = init_red_mod_units(bnf,prec);
1689 1107810 : if (!(u = red_mod_units(col,z))) return NULL;
1690 1107813 : col = RgC_add(col, RgM_RgC_mul(logfu, u));
1691 1107813 : if (!(col = cleanarch(col,N,NULL,prec))) return NULL;
1692 : }
1693 1237085 : s = divru(mulir(e, glog(kNx,prec)), N);
1694 1237079 : col = fixarch(col, s, nf_get_r1(nf));
1695 1237088 : if (RgC_expbitprec(col) >= 0) return NULL;
1696 1236514 : col = gexp(col, prec);
1697 : /* d.alpha such that x = alpha \prod gj^ej */
1698 1236514 : x = RgM_solve_realimag(M,col); if (!x) return NULL;
1699 1236511 : x = RgC_Rg_mul(x, dx);
1700 1236510 : y = grndtoi(x, pe);
1701 1236512 : if (*pe > -5) { *pe = needed_bitprec(x); return NULL; }
1702 1232655 : return RgC_Rg_div(y, dx);
1703 : }
1704 :
1705 : /* y = C \prod g[i]^e[i] ? */
1706 : static int
1707 1228546 : fact_ok(GEN nf, GEN y, GEN C, GEN g, GEN e)
1708 : {
1709 1228546 : pari_sp av = avma;
1710 1228546 : long i, c = lg(e);
1711 1228546 : GEN z = C? C: gen_1;
1712 1505577 : for (i=1; i<c; i++)
1713 277031 : if (signe(gel(e,i))) z = idealmul(nf, z, idealpow(nf, gel(g,i), gel(e,i)));
1714 1228546 : if (typ(z) != t_MAT) z = idealhnf_shallow(nf,z);
1715 1228546 : if (typ(y) != t_MAT) y = idealhnf_shallow(nf,y);
1716 1228546 : return gc_bool(av, ZM_equal(y,z));
1717 : }
1718 : static GEN
1719 1392550 : ZV_divrem(GEN A, GEN B, GEN *pR)
1720 : {
1721 1392550 : long i, l = lg(A);
1722 1392550 : GEN Q = cgetg(l, t_COL), R = cgetg(l, t_COL);
1723 1899433 : for (i = 1; i < l; i++) gel(Q,i) = truedvmdii(gel(A,i), gel(B,i), &gel(R,i));
1724 1392549 : *pR = R; return Q;
1725 : }
1726 :
1727 : static GEN
1728 1392549 : Ur_ZC_mul(GEN bnf, GEN v)
1729 : {
1730 1392549 : GEN w, U = bnf_get_Ur(bnf);
1731 1392549 : long i, l = lg(bnf_get_cyc(bnf)); /* may be < lgcols(U) */
1732 :
1733 1392549 : w = cgetg(l, t_COL);
1734 1899435 : for (i = 1; i < l; i++) gel(w,i) = ZMrow_ZC_mul(U, v, i);
1735 1392551 : return w;
1736 : }
1737 :
1738 : static GEN
1739 7325 : ZV_mul(GEN x, GEN y)
1740 : {
1741 7325 : long i, l = lg(x);
1742 7325 : GEN z = cgetg(l, t_COL);
1743 31948 : for (i = 1; i < l; i++) gel(z,i) = mulii(gel(x,i), gel(y,i));
1744 7325 : return z;
1745 : }
1746 : static int
1747 1228033 : dump_gen(GEN SUnits, GEN x, long flag)
1748 : {
1749 : GEN d;
1750 : long e;
1751 1228033 : if (!(flag & nf_GENMAT) || !SUnits) return 0;
1752 272953 : e = gexpo(gel(SUnits,2)); if (e > 64) return 0; /* U large */
1753 272857 : x = Q_remove_denom(x, &d);
1754 272855 : return (d && expi(d) > 32) || gexpo(x) > 32;
1755 : }
1756 :
1757 : /* assume x in HNF; cf class_group_gen for notations. Return NULL iff
1758 : * flag & nf_FORCE and computation of principal ideal generator fails */
1759 : static GEN
1760 1408721 : isprincipalall(GEN bnf, GEN x, long *pprec, long flag)
1761 : {
1762 : GEN xar, Wex, Bex, gen, xc, col, A, Q, R, UA, SUnits;
1763 1408721 : GEN C = bnf_get_C(bnf), nf = bnf_get_nf(bnf), cyc = bnf_get_cyc(bnf);
1764 : long nB, nW, e;
1765 :
1766 1408721 : if (lg(cyc) == 1 && !(flag & (nf_GEN|nf_GENMAT|nf_GEN_IF_PRINCIPAL)))
1767 4732 : return cgetg(1,t_COL);
1768 1403989 : if (lg(x) == 2)
1769 : { /* nf = Q */
1770 84 : col = gel(x,1);
1771 84 : if (flag & nf_GENMAT) col = Q_to_famat(gel(col,1));
1772 84 : return (flag & nf_GEN_IF_PRINCIPAL)? col: mkvec2(cgetg(1,t_COL), col);
1773 : }
1774 :
1775 1403905 : x = Q_primitive_part(x, &xc);
1776 1403903 : if (equali1(gcoeff(x,1,1))) /* trivial ideal */
1777 : {
1778 11354 : R = zerocol(lg(cyc)-1);
1779 11354 : if (!(flag & (nf_GEN|nf_GENMAT|nf_GEN_IF_PRINCIPAL))) return R;
1780 11305 : if (flag & nf_GEN_IF_PRINCIPAL)
1781 6468 : return scalarcol_shallow(xc? xc: gen_1, nf_get_degree(nf));
1782 4837 : if (flag & nf_GENMAT)
1783 2191 : col = xc? Q_to_famat(xc): trivial_fact();
1784 : else
1785 2646 : col = scalarcol_shallow(xc? xc: gen_1, nf_get_degree(nf));
1786 4837 : return mkvec2(R, col);
1787 : }
1788 1392549 : xar = split_ideal(bnf, x, &Wex, &Bex);
1789 : /* x = g_W Wex + g_B Bex + [xar] = g_W (Wex - B*Bex) + [xar] + [C_B]Bex */
1790 1392549 : A = zc_to_ZC(Wex); nB = lg(Bex)-1;
1791 1392549 : if (nB) A = ZC_sub(A, ZM_zc_mul(bnf_get_B(bnf), Bex));
1792 1392549 : UA = Ur_ZC_mul(bnf, A);
1793 1392550 : Q = ZV_divrem(UA, cyc, &R);
1794 : /* g_W (Wex - B*Bex) = G Ur A - [ga]A = G R + [GD]Q - [ga]A
1795 : * Finally: x = G R + [xar] + [C_B]Bex + [GD]Q - [ga]A */
1796 1392549 : if (!(flag & (nf_GEN|nf_GENMAT|nf_GEN_IF_PRINCIPAL))) return R;
1797 1232078 : if ((flag & nf_GEN_IF_PRINCIPAL) && !ZV_equal0(R)) return gen_0;
1798 :
1799 1232071 : nW = lg(Wex)-1;
1800 1232071 : gen = bnf_get_gen(bnf);
1801 1232071 : col = NULL;
1802 1232071 : SUnits = bnf_get_sunits(bnf);
1803 1232071 : if (lg(R) == 1
1804 272424 : || abscmpiu(gel(R,vecindexmax(R)), 4 * (*pprec)) < 0)
1805 : { /* q = N (x / prod gj^ej) = N(alpha), denom(alpha) | d */
1806 1231474 : GEN d, q = gdiv(ZM_det_triangular(x), get_norm_fact(gen, R, &d));
1807 1231476 : col = xar? nf_cxlog(nf, xar, *pprec): NULL;
1808 1231476 : if (nB) col = add(col, act_arch(Bex, nW? vecslice(C,nW+1,lg(C)-1): C));
1809 1231472 : if (nW) col = add(col, RgC_sub(act_arch(Q, bnf_get_GD(bnf)),
1810 : act_arch(A, bnf_get_ga(bnf))));
1811 1231472 : col = isprincipalarch(bnf, col, q, gen_1, d, &e);
1812 1231476 : if (col && (dump_gen(SUnits, col, flag)
1813 1228031 : || !fact_ok(nf,x, col,gen,R))) col = NULL;
1814 : }
1815 1232069 : if (!col && (flag & nf_GENMAT))
1816 : {
1817 8066 : if (SUnits)
1818 : {
1819 7577 : GEN X = gel(SUnits,1), U = gel(SUnits,2), C = gel(SUnits,3);
1820 7577 : GEN v = gel(bnf,9), Ge = gel(v,4), M1 = gel(v,5), M2 = gel(v,6);
1821 7577 : GEN z = NULL, F = NULL;
1822 7577 : if (nB)
1823 : {
1824 7577 : GEN C2 = nW? vecslice(C, nW+1, lg(C)-1): C;
1825 7577 : z = ZM_zc_mul(C2, Bex);
1826 : }
1827 7577 : if (nW)
1828 : { /* [GD]Q - [ga]A = ([X]M1 - [Ge]D) Q - ([X]M2 - [Ge]Ur) A */
1829 7325 : GEN C1 = vecslice(C, 1, nW);
1830 7325 : GEN v = ZC_sub(ZM_ZC_mul(M1,Q), ZM_ZC_mul(M2,A));
1831 7325 : z = add(z, ZM_ZC_mul(C1, v));
1832 7325 : F = famat_reduce(famatV_factorback(Ge, ZC_sub(UA, ZV_mul(cyc,Q))));
1833 7325 : if (lgcols(F) == 1) F = NULL;
1834 : }
1835 : /* reduce modulo units and Q^* */
1836 7577 : if (lg(U) != 1) z = ZC_sub(z, ZM_ZC_mul(U, RgM_Babai(U,z)));
1837 7577 : col = mkmat2(X, z);
1838 7577 : if (F) col = famat_mul_shallow(col, F);
1839 7577 : col = famat_remove_trivial(col);
1840 7577 : if (xar) col = famat_mul_shallow(col, xar);
1841 : }
1842 489 : else if (!ZV_equal0(R))
1843 : { /* in case isprincipalfact calls bnfinit() due to prec trouble...*/
1844 483 : GEN y = isprincipalfact(bnf, x, gen, ZC_neg(R), flag);
1845 483 : if (typ(y) != t_VEC) return y;
1846 483 : col = gel(y,2);
1847 : }
1848 : }
1849 1232069 : if (col)
1850 : { /* add back missing content */
1851 1231979 : if (typ(col) == t_MAT)
1852 8060 : { if (xc) col = famat_mul_shallow(col, xc); }
1853 1223919 : else if (flag & nf_GENMAT)
1854 : {
1855 : GEN c;
1856 1210220 : if (RgV_isscalar(col))
1857 3651 : col = Q_to_famat(mul_content(xc, gel(col,1)));
1858 : else
1859 : {
1860 1206568 : col = Q_primitive_part(col, &c);
1861 1206569 : col = to_famat_shallow(col, gen_1);
1862 1206570 : xc = mul_content(xc, c);
1863 1206570 : if (xc) col = famat_mul(col, Q_to_famat(xc));
1864 : }
1865 : }
1866 : else
1867 13699 : { if (xc) col = RgC_Rg_mul(col,xc); }
1868 : }
1869 : else
1870 : {
1871 90 : if (e < 0) e = 0;
1872 90 : *pprec += nbits2extraprec(e + 128);
1873 90 : if (flag & nf_FORCE)
1874 : {
1875 76 : if (DEBUGLEVEL)
1876 0 : pari_warn(warner,"precision too low for generators, e = %ld",e);
1877 76 : return NULL;
1878 : }
1879 14 : pari_warn(warner,"precision too low for generators, not given");
1880 14 : col = cgetg(1, t_COL);
1881 : }
1882 1231995 : return (flag & nf_GEN_IF_PRINCIPAL)? col: mkvec2(R, col);
1883 : }
1884 :
1885 : static GEN
1886 464660 : triv_gen(GEN bnf, GEN x, long flag)
1887 : {
1888 464660 : pari_sp av = avma;
1889 464660 : GEN nf = bnf_get_nf(bnf);
1890 : long c;
1891 464660 : if (flag & nf_GEN_IF_PRINCIPAL)
1892 : {
1893 7 : if (!(flag & nf_GENMAT)) return algtobasis(nf,x);
1894 7 : x = nf_to_scalar_or_basis(nf,x);
1895 7 : if (typ(x) == t_INT && is_pm1(x)) return trivial_fact();
1896 0 : return gc_GEN(av, to_famat_shallow(x, gen_1));
1897 : }
1898 464653 : c = lg(bnf_get_cyc(bnf)) - 1;
1899 464653 : if (flag & nf_GENMAT)
1900 455049 : retmkvec2(zerocol(c), to_famat_shallow(algtobasis(nf,x), gen_1));
1901 9604 : if (flag & nf_GEN)
1902 28 : retmkvec2(zerocol(c), algtobasis(nf,x));
1903 9576 : return zerocol(c);
1904 : }
1905 :
1906 : GEN
1907 1841332 : bnfisprincipal0(GEN bnf,GEN x,long flag)
1908 : {
1909 1841332 : pari_sp av = avma;
1910 : GEN c, nf;
1911 : long pr;
1912 :
1913 1841332 : bnf = checkbnf(bnf);
1914 1841333 : nf = bnf_get_nf(bnf);
1915 1841333 : switch( idealtyp(&x, NULL) )
1916 : {
1917 59080 : case id_PRINCIPAL:
1918 59080 : if (gequal0(x)) pari_err_DOMAIN("bnfisprincipal","ideal","=",gen_0,x);
1919 59080 : return triv_gen(bnf, x, flag);
1920 1758592 : case id_PRIME:
1921 1758592 : if (pr_is_inert(x)) return triv_gen(bnf, pr_get_p(x), flag);
1922 1353019 : x = pr_hnf(nf, x);
1923 1353022 : break;
1924 23660 : case id_MAT:
1925 23660 : if (lg(x)==1) pari_err_DOMAIN("bnfisprincipal","ideal","=",gen_0,x);
1926 23660 : if (nf_get_degree(nf) != lg(x)-1)
1927 0 : pari_err_TYPE("idealtyp [dimension != degree]", x);
1928 : }
1929 1376683 : pr = prec_arch(bnf); /* precision of unit matrix */
1930 1376681 : c = getrand();
1931 : for (;;)
1932 6 : {
1933 1376688 : pari_sp av1 = avma;
1934 1376688 : GEN y = isprincipalall(bnf,x,&pr,flag);
1935 1376686 : if (y) return gc_GEN(av, y);
1936 :
1937 6 : if (DEBUGLEVEL) pari_warn(warnprec,"isprincipal",pr);
1938 6 : set_avma(av1); bnf = bnfnewprec_shallow(bnf,pr); setrand(c);
1939 : }
1940 : }
1941 : GEN
1942 174779 : isprincipal(GEN bnf,GEN x) { return bnfisprincipal0(bnf,x,0); }
1943 :
1944 : /* FIXME: OBSOLETE */
1945 : GEN
1946 0 : isprincipalgen(GEN bnf,GEN x)
1947 0 : { return bnfisprincipal0(bnf,x,nf_GEN); }
1948 : GEN
1949 0 : isprincipalforce(GEN bnf,GEN x)
1950 0 : { return bnfisprincipal0(bnf,x,nf_FORCE); }
1951 : GEN
1952 0 : isprincipalgenforce(GEN bnf,GEN x)
1953 0 : { return bnfisprincipal0(bnf,x,nf_GEN | nf_FORCE); }
1954 :
1955 : /* lg(u) > 1 */
1956 : static int
1957 105 : RgV_is1(GEN u) { return isint1(gel(u,1)) && RgV_isscalar(u); }
1958 : static GEN
1959 31963 : add_principal_part(GEN nf, GEN u, GEN v, long flag)
1960 : {
1961 31963 : if (flag & nf_GENMAT)
1962 14540 : return (typ(u) == t_COL && RgV_is1(u))? v: famat_mul_shallow(v,u);
1963 : else
1964 17423 : return nfmul(nf, v, u);
1965 : }
1966 :
1967 : #if 0
1968 : /* compute C prod P[i]^e[i], e[i] >=0 for all i. C may be NULL (omitted)
1969 : * e destroyed ! */
1970 : static GEN
1971 : expand(GEN nf, GEN C, GEN P, GEN e)
1972 : {
1973 : long i, l = lg(e), done = 1;
1974 : GEN id = C;
1975 : for (i=1; i<l; i++)
1976 : {
1977 : GEN ei = gel(e,i);
1978 : if (signe(ei))
1979 : {
1980 : if (mod2(ei)) id = id? idealmul(nf, id, gel(P,i)): gel(P,i);
1981 : ei = shifti(ei,-1);
1982 : if (signe(ei)) done = 0;
1983 : gel(e,i) = ei;
1984 : }
1985 : }
1986 : if (id != C) id = idealred(nf, id);
1987 : if (done) return id;
1988 : return idealmulred(nf, id, idealsqr(nf, expand(nf,id,P,e)));
1989 : }
1990 : /* C is an extended ideal, possibly with C[1] = NULL */
1991 : static GEN
1992 : expandext(GEN nf, GEN C, GEN P, GEN e)
1993 : {
1994 : long i, l = lg(e), done = 1;
1995 : GEN A = gel(C,1);
1996 : for (i=1; i<l; i++)
1997 : {
1998 : GEN ei = gel(e,i);
1999 : if (signe(ei))
2000 : {
2001 : if (mod2(ei)) A = A? idealmul(nf, A, gel(P,i)): gel(P,i);
2002 : ei = shifti(ei,-1);
2003 : if (signe(ei)) done = 0;
2004 : gel(e,i) = ei;
2005 : }
2006 : }
2007 : if (A == gel(C,1))
2008 : A = C;
2009 : else
2010 : A = idealred(nf, mkvec2(A, gel(C,2)));
2011 : if (done) return A;
2012 : return idealmulred(nf, A, idealsqr(nf, expand(nf,A,P,e)));
2013 : }
2014 : #endif
2015 :
2016 : static GEN
2017 0 : expand(GEN nf, GEN C, GEN P, GEN e)
2018 : {
2019 0 : long i, l = lg(e);
2020 0 : GEN B, A = C;
2021 0 : for (i=1; i<l; i++) /* compute prod P[i]^e[i] */
2022 0 : if (signe(gel(e,i)))
2023 : {
2024 0 : B = idealpowred(nf, gel(P,i), gel(e,i));
2025 0 : A = A? idealmulred(nf,A,B): B;
2026 : }
2027 0 : return A;
2028 : }
2029 : static GEN
2030 31984 : expandext(GEN nf, GEN C, GEN P, GEN e)
2031 : {
2032 31984 : long i, l = lg(e);
2033 31984 : GEN B, A = gel(C,1), C1 = A;
2034 94800 : for (i=1; i<l; i++) /* compute prod P[i]^e[i] */
2035 62816 : if (signe(gel(e,i)))
2036 : {
2037 35006 : gel(C,1) = gel(P,i);
2038 35006 : B = idealpowred(nf, C, gel(e,i));
2039 35006 : A = A? idealmulred(nf,A,B): B;
2040 : }
2041 31984 : return A == C1? C: A;
2042 : }
2043 :
2044 : /* isprincipal for C * \prod P[i]^e[i] (C omitted if NULL) */
2045 : GEN
2046 31984 : isprincipalfact(GEN bnf, GEN C, GEN P, GEN e, long flag)
2047 : {
2048 31984 : const long gen = flag & (nf_GEN|nf_GENMAT|nf_GEN_IF_PRINCIPAL);
2049 : long prec;
2050 31984 : pari_sp av = avma;
2051 31984 : GEN C0, Cext, c, id, nf = bnf_get_nf(bnf);
2052 :
2053 31984 : if (gen)
2054 : {
2055 14547 : Cext = (flag & nf_GENMAT)? trivial_fact()
2056 31984 : : mkpolmod(gen_1,nf_get_pol(nf));
2057 31984 : C0 = mkvec2(C, Cext);
2058 31984 : id = expandext(nf, C0, P, e);
2059 : } else {
2060 0 : Cext = NULL;
2061 0 : C0 = C;
2062 0 : id = expand(nf, C, P, e);
2063 : }
2064 31984 : if (id == C0) /* e = 0 */
2065 : {
2066 12477 : if (!C) return bnfisprincipal0(bnf, gen_1, flag);
2067 12463 : switch(typ(C))
2068 : {
2069 7 : case t_INT: case t_FRAC: case t_POL: case t_POLMOD: case t_COL:
2070 7 : return triv_gen(bnf, C, flag);
2071 : }
2072 12456 : C = idealhnf_shallow(nf,C);
2073 : }
2074 : else
2075 : {
2076 19507 : if (gen) { C = gel(id,1); Cext = gel(id,2); } else C = id;
2077 : }
2078 31963 : prec = prec_arch(bnf);
2079 31963 : c = getrand();
2080 : for (;;)
2081 70 : {
2082 32033 : pari_sp av1 = avma;
2083 32033 : GEN y = isprincipalall(bnf, C, &prec, flag);
2084 32033 : if (y)
2085 : {
2086 31963 : if (flag & nf_GEN_IF_PRINCIPAL)
2087 : {
2088 20818 : if (typ(y) == t_INT) return gc_NULL(av);
2089 20818 : y = add_principal_part(nf, y, Cext, flag);
2090 : }
2091 : else
2092 : {
2093 11145 : GEN u = gel(y,2);
2094 11145 : if (!gen || typ(y) != t_VEC) return gc_upto(av,y);
2095 11145 : if (lg(u) != 1) gel(y,2) = add_principal_part(nf, u, Cext, flag);
2096 : }
2097 31962 : return gc_GEN(av, y);
2098 : }
2099 70 : if (DEBUGLEVEL) pari_warn(warnprec,"isprincipal",prec);
2100 70 : set_avma(av1); bnf = bnfnewprec_shallow(bnf,prec); setrand(c);
2101 : }
2102 : }
2103 : GEN
2104 0 : isprincipalfact_or_fail(GEN bnf, GEN C, GEN P, GEN e)
2105 : {
2106 0 : const long flag = nf_GENMAT|nf_FORCE;
2107 : long prec;
2108 0 : pari_sp av = avma;
2109 0 : GEN u, y, id, C0, Cext, nf = bnf_get_nf(bnf);
2110 :
2111 0 : Cext = trivial_fact();
2112 0 : C0 = mkvec2(C, Cext);
2113 0 : id = expandext(nf, C0, P, e);
2114 0 : if (id == C0) /* e = 0 */
2115 0 : C = idealhnf_shallow(nf,C);
2116 : else {
2117 0 : C = gel(id,1); Cext = gel(id,2);
2118 : }
2119 0 : prec = prec_arch(bnf);
2120 0 : y = isprincipalall(bnf, C, &prec, flag);
2121 0 : if (!y) return gc_utoipos(av, prec);
2122 0 : u = gel(y,2);
2123 0 : if (lg(u) != 1) gel(y,2) = add_principal_part(nf, u, Cext, flag);
2124 0 : return gc_GEN(av, y);
2125 : }
2126 :
2127 : GEN
2128 162325 : nfsign_from_logarch(GEN LA, GEN invpi, GEN archp)
2129 : {
2130 162325 : long l = lg(archp), i;
2131 162325 : GEN y = cgetg(l, t_VECSMALL);
2132 162327 : pari_sp av = avma;
2133 :
2134 310066 : for (i=1; i<l; i++)
2135 : {
2136 147738 : GEN c = ground( gmul(imag_i(gel(LA,archp[i])), invpi) );
2137 147741 : y[i] = mpodd(c)? 1: 0;
2138 : }
2139 162328 : set_avma(av); return y;
2140 : }
2141 :
2142 : GEN
2143 239979 : nfsign_tu(GEN bnf, GEN archp)
2144 : {
2145 : long n;
2146 239979 : if (bnf_get_tuN(bnf) != 2) return cgetg(1, t_VECSMALL);
2147 172876 : n = archp? lg(archp) - 1: nf_get_r1(bnf_get_nf(bnf));
2148 172876 : return const_vecsmall(n, 1);
2149 : }
2150 : GEN
2151 241155 : nfsign_fu(GEN bnf, GEN archp)
2152 : {
2153 241155 : GEN invpi, y, A = bnf_get_logfu(bnf), nf = bnf_get_nf(bnf);
2154 241155 : long j = 1, RU = lg(A);
2155 :
2156 241155 : if (!archp) archp = identity_perm( nf_get_r1(nf) );
2157 241155 : invpi = invr( mppi(nf_get_prec(nf)) );
2158 241150 : y = cgetg(RU,t_MAT);
2159 403388 : for (j = 1; j < RU; j++)
2160 162227 : gel(y,j) = nfsign_from_logarch(gel(A,j), invpi, archp);
2161 241161 : return y;
2162 : }
2163 : GEN
2164 35 : nfsign_units(GEN bnf, GEN archp, int add_zu)
2165 : {
2166 35 : GEN sfu = nfsign_fu(bnf, archp);
2167 35 : return add_zu? vec_prepend(sfu, nfsign_tu(bnf, archp)): sfu;
2168 : }
2169 :
2170 : /* obsolete */
2171 : GEN
2172 7 : signunits(GEN bnf)
2173 : {
2174 : pari_sp av;
2175 : GEN S, y, nf;
2176 : long i, j, r1, r2;
2177 :
2178 7 : bnf = checkbnf(bnf); nf = bnf_get_nf(bnf);
2179 7 : nf_get_sign(nf, &r1,&r2);
2180 7 : S = zeromatcopy(r1, r1+r2-1); av = avma;
2181 7 : y = nfsign_fu(bnf, NULL);
2182 14 : for (j = 1; j < lg(y); j++)
2183 : {
2184 7 : GEN Sj = gel(S,j), yj = gel(y,j);
2185 21 : for (i = 1; i <= r1; i++) gel(Sj,i) = yj[i]? gen_m1: gen_1;
2186 : }
2187 7 : set_avma(av); return S;
2188 : }
2189 :
2190 : static GEN
2191 726091 : get_log_embed(REL_t *rel, GEN M, long RU, long R1, long prec)
2192 : {
2193 726091 : GEN arch, C, z = rel->m;
2194 : long i;
2195 726091 : arch = typ(z) == t_COL? RgM_RgC_mul(M, z): const_col(nbrows(M), z);
2196 726098 : C = cgetg(RU+1, t_COL); arch = glog(arch, prec);
2197 1654551 : for (i=1; i<=R1; i++) gel(C,i) = gel(arch,i);
2198 1565321 : for ( ; i<=RU; i++) gel(C,i) = gmul2n(gel(arch,i), 1);
2199 726079 : return C;
2200 : }
2201 : static GEN
2202 1017064 : rel_embed(REL_t *rel, FB_t *F, GEN embs, long ind, GEN M, long RU, long R1,
2203 : long prec)
2204 : {
2205 : GEN C, D, perm;
2206 : long i, n;
2207 1017064 : if (!rel->relaut) return get_log_embed(rel, M, RU, R1, prec);
2208 : /* image of another relation by automorphism */
2209 290973 : C = gel(embs, ind - rel->relorig);
2210 290973 : perm = gel(F->embperm, rel->relaut);
2211 290973 : D = cgetg_copy(C, &n);
2212 1217812 : for (i = 1; i < n; i++)
2213 : {
2214 926831 : long v = perm[i];
2215 926831 : gel(D,i) = (v > 0)? gel(C,v): conj_i(gel(C,-v));
2216 : }
2217 290981 : return D;
2218 : }
2219 : static GEN
2220 122410 : get_embs(FB_t *F, RELCACHE_t *cache, GEN nf, GEN embs, long PREC)
2221 : {
2222 122410 : long ru, j, k, l = cache->last - cache->chk + 1, r1 = nf_get_r1(nf);
2223 122410 : GEN M = nf_get_M(nf), nembs = cgetg(cache->last - cache->base+1, t_MAT);
2224 : REL_t *rel;
2225 :
2226 1517102 : for (k = 1; k <= cache->chk - cache->base; k++) gel(nembs,k) = gel(embs,k);
2227 122409 : embs = nembs; ru = nbrows(M);
2228 1126890 : for (j=1,rel = cache->chk + 1; j < l; rel++,j++,k++)
2229 1004488 : gel(embs,k) = rel_embed(rel, F, embs, k, M, ru, r1, PREC);
2230 122402 : return embs;
2231 : }
2232 : static void
2233 948646 : set_rel_alpha(REL_t *rel, GEN auts, GEN vA, long ind)
2234 : {
2235 : GEN u;
2236 948646 : if (!rel->relaut)
2237 674823 : u = rel->m;
2238 : else
2239 273823 : u = ZM_ZC_mul(gel(auts, rel->relaut), gel(vA, ind - rel->relorig));
2240 948643 : gel(vA, ind) = u;
2241 948643 : }
2242 : static GEN
2243 2315386 : set_fact(FB_t *F, FACT *fact, GEN e, long *pnz)
2244 : {
2245 2315386 : long i, n = fact[0].pr, nz = F->KC + 1;
2246 2315386 : GEN c = zero_Flv(F->KC);
2247 10861413 : for (i = 1; i <= n; i++)
2248 : {
2249 8546022 : long p = fact[i].pr;
2250 8546022 : if (p < nz) nz = p;
2251 8546022 : c[p] = fact[i].ex;
2252 : }
2253 2315391 : if (e)
2254 : {
2255 114517 : long l = lg(e);
2256 334252 : for (i = 1; i < l; i++)
2257 219735 : if (e[i]) { long v = F->subFB[i]; c[v] += e[i]; if (v < nz) nz = v; }
2258 : }
2259 2315391 : *pnz = nz; return c;
2260 : }
2261 :
2262 : /* Is cols already in the cache ? bs = index of first non zero coeff in cols
2263 : * General check for colinearity useless since exceedingly rare */
2264 : static int
2265 2968710 : already_known(RELCACHE_t *cache, long bs, GEN cols)
2266 : {
2267 : REL_t *r;
2268 2968710 : long l = lg(cols);
2269 220206453 : for (r = cache->last; r > cache->base; r--)
2270 217727466 : if (bs == r->nz)
2271 : {
2272 34413572 : GEN coll = r->R;
2273 34413572 : long b = bs;
2274 124138000 : while (b < l && cols[b] == coll[b]) b++;
2275 34413572 : if (b == l) return 1;
2276 : }
2277 2478987 : return 0;
2278 : }
2279 :
2280 : /* Add relation R to cache, nz = index of first non zero coeff in R.
2281 : * If relation is a linear combination of the previous ones, return 0.
2282 : * Otherwise, update basis and return > 0. Compute mod p (much faster)
2283 : * so some kernel vector might not be genuine. */
2284 : static int
2285 2972852 : add_rel_i(RELCACHE_t *cache, GEN R, long nz, GEN m, long orig, long aut, REL_t **relp, long in_rnd_rel)
2286 : {
2287 2972852 : long i, k, n = lg(R)-1;
2288 :
2289 2972852 : if (nz == n+1) { k = 0; goto ADD_REL; }
2290 2968705 : if (already_known(cache, nz, R)) return -1;
2291 2479042 : if (cache->last >= cache->base + cache->len) return 0;
2292 2479042 : if (DEBUGLEVEL>6)
2293 : {
2294 0 : err_printf("adding vector = %Ps\n",R);
2295 0 : err_printf("generators =\n%Ps\n", cache->basis);
2296 : }
2297 2479059 : if (cache->missing)
2298 : {
2299 2082767 : GEN a = leafcopy(R), basis = cache->basis;
2300 2082763 : k = lg(a);
2301 127975974 : do --k; while (!a[k]);
2302 7310661 : while (k)
2303 : {
2304 5697540 : GEN c = gel(basis, k);
2305 5697540 : if (c[k])
2306 : {
2307 5227898 : long ak = a[k];
2308 269126373 : for (i=1; i < k; i++) if (c[i]) a[i] = (a[i] + ak*(mod_p-c[i])) % mod_p;
2309 5227898 : a[k] = 0;
2310 132735985 : do --k; while (!a[k]); /* k cannot go below 0: codeword is a sentinel */
2311 : }
2312 : else
2313 : {
2314 469642 : ulong invak = Fl_inv(uel(a,k), mod_p);
2315 : /* Cleanup a */
2316 14002440 : for (i = k; i-- > 1; )
2317 : {
2318 13532798 : long j, ai = a[i];
2319 13532798 : c = gel(basis, i);
2320 13532798 : if (!ai || !c[i]) continue;
2321 237847 : ai = mod_p-ai;
2322 4508750 : for (j = 1; j < i; j++) if (c[j]) a[j] = (a[j] + ai*c[j]) % mod_p;
2323 237847 : a[i] = 0;
2324 : }
2325 : /* Insert a/a[k] as k-th column */
2326 469642 : c = gel(basis, k);
2327 14002441 : for (i = 1; i<k; i++) if (a[i]) c[i] = (a[i] * invak) % mod_p;
2328 469642 : c[k] = 1; a = c;
2329 : /* Cleanup above k */
2330 13796494 : for (i = k+1; i<n; i++)
2331 : {
2332 : long j, ck;
2333 13326852 : c = gel(basis, i);
2334 13326852 : ck = c[k];
2335 13326852 : if (!ck) continue;
2336 2764952 : ck = mod_p-ck;
2337 100012552 : for (j = 1; j < k; j++) if (a[j]) c[j] = (c[j] + ck*a[j]) % mod_p;
2338 2764952 : c[k] = 0;
2339 : }
2340 469642 : cache->missing--;
2341 469642 : break;
2342 : }
2343 : }
2344 : }
2345 : else
2346 396292 : k = (cache->last - cache->base) + 1;
2347 2479055 : if (k || cache->relsup > 0 || (m && in_rnd_rel))
2348 : {
2349 : REL_t *rel;
2350 :
2351 989128 : ADD_REL:
2352 993275 : rel = ++cache->last;
2353 993275 : if (!k && cache->relsup && nz < n+1)
2354 : {
2355 122849 : cache->relsup--;
2356 122849 : k = (rel - cache->base) + cache->missing;
2357 : }
2358 993275 : rel->R = gclone(R);
2359 993285 : rel->m = m ? gclone(m) : NULL;
2360 993277 : rel->nz = nz;
2361 993277 : if (aut)
2362 : {
2363 288377 : rel->relorig = (rel - cache->base) - orig;
2364 288377 : rel->relaut = aut;
2365 : }
2366 : else
2367 704900 : rel->relaut = 0;
2368 993277 : if (relp) *relp = rel;
2369 993277 : if (DEBUGLEVEL) dbg_newrel(cache);
2370 : }
2371 2483200 : return k;
2372 : }
2373 :
2374 : /* m a t_INT or primitive t_COL */
2375 : static int
2376 2488410 : add_rel(RELCACHE_t *cache, FB_t *F, GEN R, long nz, GEN m, long in_rnd_rel)
2377 : {
2378 : REL_t *rel;
2379 : long k, l, reln;
2380 2488410 : const long lauts = lg(F->idealperm), KC = F->KC;
2381 :
2382 2488410 : k = add_rel_i(cache, R, nz, m, 0, 0, &rel, in_rnd_rel);
2383 2488463 : if (k > 0 && typ(m) != t_INT)
2384 : {
2385 531393 : GEN Rl = cgetg(KC+1, t_VECSMALL);
2386 531395 : reln = rel - cache->base;
2387 1015855 : for (l = 1; l < lauts; l++)
2388 : {
2389 484452 : GEN perml = gel(F->idealperm, l);
2390 484452 : long i, nzl = perml[nz];
2391 :
2392 20547061 : for (i = 1; i <= KC; i++) Rl[i] = 0;
2393 18319287 : for (i = nz; i <= KC; i++)
2394 17834835 : if (R[i])
2395 : {
2396 1288183 : long v = perml[i];
2397 :
2398 1288183 : if (v < nzl) nzl = v;
2399 1288183 : Rl[v] = R[i];
2400 : }
2401 484452 : (void)add_rel_i(cache, Rl, nzl, NULL, reln, l, NULL, in_rnd_rel);
2402 : }
2403 : }
2404 2488473 : return k;
2405 : }
2406 :
2407 : INLINE void
2408 28415394 : step(GEN x, double *y, GEN inc, long k)
2409 : {
2410 28415394 : if (!y[k])
2411 2228806 : x[k]++; /* leading coeff > 0 */
2412 : else
2413 : {
2414 26186588 : long i = inc[k];
2415 26186588 : x[k] += i;
2416 26186588 : inc[k] = (i > 0)? -1-i: 1-i;
2417 : }
2418 28415394 : }
2419 :
2420 : /* degree n >= 2 */
2421 : static double
2422 266571 : Fincke_Pohst_bound(double T, GEN r)
2423 : {
2424 266571 : pari_sp av = avma;
2425 266571 : GEN zT = dbltor(T * T), p = gmael(r,1,1), B = NULL;
2426 266564 : long i, n = lg(r)-1;
2427 266564 : double g = 0.;
2428 611379 : for (i = 2; i <= n; i++)
2429 : {
2430 611365 : p = gmul(p, gmael(r,i,i));
2431 611388 : B = sqrtnr(gmul(zT,p), i);
2432 611372 : if (i == n || cmprr(B, gmael(r,i+1,i+1)) < 0) break;
2433 : }
2434 266574 : if (!gisdouble(B,&g)) g = 0.;
2435 266565 : return gc_double(av, g);
2436 : }
2437 :
2438 : static void
2439 2052113 : fact_update(GEN R, FB_t *F, long ipr, GEN c)
2440 : {
2441 2052113 : GEN pr = gel(F->LP,ipr), p = pr_get_p(pr);
2442 2052111 : long v = Z_lval(c, itou(p));
2443 2052127 : if (v) R[ipr] -= pr_get_e(pr) * v;
2444 2052126 : }
2445 :
2446 : static long
2447 266572 : Fincke_Pohst_ideal_both(RELCACHE_t *cache, GEN V, FB_t *F, GEN nf, GEN I, GEN NI,
2448 : FACT *fact, long max_FACT, long Nrelid, FP_t *fp, GEN rex, long jid, long jid0, long e0,
2449 : long *Nsmall, long *Nfact)
2450 : {
2451 : pari_sp av;
2452 266572 : GEN G = nf_get_G(nf), G0 = nf_get_roundG(nf), r, u, gx, cgx, inc, ideal;
2453 266569 : long prec = nf_get_prec(nf), N = nf_get_degree(nf);
2454 266566 : long j, k, skipfirst, relid = 0, try_factor = 0, rel = 1;
2455 266566 : long try_elt = 0, maxtry_ELEMENT = 4*max_FACT*max_FACT;
2456 : double BOUND, B1, B2;
2457 :
2458 266566 : inc = const_vecsmall(N, 1);
2459 266565 : u = ZM_lll(ZM_mul(G0, I), 0.99, LLL_IM);
2460 266573 : ideal = ZM_mul(I,u); /* approximate T2-LLL reduction */
2461 266551 : r = gaussred_from_QR(RgM_mul(G, ideal), prec); /* Cholesky for T2 | ideal */
2462 266575 : if (!r) pari_err_BUG("small_norm (precision too low)");
2463 :
2464 1184216 : for (k=1; k<=N; k++)
2465 : {
2466 917645 : if (!gisdouble(gcoeff(r,k,k),&(fp->v[k]))) return 0;
2467 2724822 : for (j=1; j<k; j++) if (!gisdouble(gcoeff(r,j,k),&(fp->q[j][k]))) return 0;
2468 917641 : if (DEBUGLEVEL>3) err_printf("v[%ld]=%.4g ",k,fp->v[k]);
2469 : }
2470 266571 : B1 = fp->v[1]; /* T2(ideal[1]) */
2471 266571 : B2 = fp->v[2] + B1 * fp->q[1][2] * fp->q[1][2]; /* T2(ideal[2]) */
2472 266571 : skipfirst = ZV_isscalar(gel(ideal,1));
2473 266571 : BOUND = maxdd(2*B2, Fincke_Pohst_bound(4 * max_FACT / F->ballvol, r));
2474 266563 : if (DEBUGLEVEL>1)
2475 : {
2476 0 : if (DEBUGLEVEL>3) err_printf("\n");
2477 0 : err_printf("BOUND = %.4g\n",BOUND);
2478 : }
2479 :
2480 266564 : k = N; fp->y[N] = fp->z[N] = 0; fp->x[N] = 0;
2481 18816107 : for (av = avma;; set_avma(av), step(fp->x,fp->y,inc,k))
2482 18549192 : {
2483 : GEN R;
2484 : long nz;
2485 : do
2486 : { /* look for primitive element of small norm, cf minim00 */
2487 23670214 : int fl = 0;
2488 : double p;
2489 23670214 : if (k > 1)
2490 : {
2491 5121086 : long l = k-1;
2492 5121086 : fp->z[l] = 0;
2493 45542105 : for (j=k; j<=N; j++) fp->z[l] += fp->q[l][j]*fp->x[j];
2494 5121086 : p = (double)fp->x[k] + fp->z[k];
2495 5121086 : fp->y[l] = fp->y[k] + p*p*fp->v[k];
2496 5121086 : if (l <= skipfirst && !fp->y[1]) fl = 1;
2497 5121086 : fp->x[l] = (long)floor(-fp->z[l] + 0.5);
2498 5121086 : k = l;
2499 : }
2500 4476289 : for(;; step(fp->x,fp->y,inc,k))
2501 : {
2502 28146460 : if (!fl)
2503 : {
2504 28106878 : if (++try_elt > maxtry_ELEMENT) goto END_Fincke_Pohst_ideal;
2505 28104085 : p = (double)fp->x[k] + fp->z[k];
2506 28104085 : if (fp->y[k] + p*p*fp->v[k] <= BOUND) break;
2507 :
2508 5389828 : step(fp->x,fp->y,inc,k);
2509 :
2510 5390035 : p = (double)fp->x[k] + fp->z[k];
2511 5390035 : if (fp->y[k] + p*p*fp->v[k] <= BOUND) break;
2512 : }
2513 4479082 : fl = 0; inc[k] = 1;
2514 4479082 : if (++k > N) goto END_Fincke_Pohst_ideal;
2515 : }
2516 23667585 : } while (k > 1);
2517 :
2518 : /* element complete */
2519 33923847 : if (zv_content(fp->x) !=1) continue; /* not primitive */
2520 15662669 : gx = ZM_zc_mul(ideal,fp->x);
2521 15662598 : if (ZV_isscalar(gx)) continue;
2522 15786161 : if (++try_factor > max_FACT) break;
2523 :
2524 15622085 : if (DEBUGLEVEL && Nsmall) (*Nsmall)++;
2525 15622085 : if (!factorgen(F,nf,I,NI,gx,fact)) continue;
2526 2413444 : if (!Nrelid) return 1;
2527 2313735 : if (jid == jid0)
2528 28724 : add_to_fact(jid, 1 + e0, fact);
2529 : else
2530 : {
2531 2285011 : add_to_fact(jid, 1, fact);
2532 2285238 : if (jid0) add_to_fact(jid0, e0, fact);
2533 : }
2534 :
2535 : /* smooth element */
2536 2313962 : R = set_fact(F, fact, rex, &nz);
2537 2313973 : cgx = Z_content(gx);
2538 2313934 : if (cgx)
2539 : { /* relatively rare, compute relation attached to gx/cgx */
2540 502269 : long i, n = fact[0].pr;
2541 502269 : gx = Q_div_to_int(gx, cgx);
2542 2504689 : for (i = 1; i <= n; i++) fact_update(R, F, fact[i].pr, cgx);
2543 502276 : if (rex)
2544 : {
2545 32891 : long l = lg(rex);
2546 108205 : for (i = 1; i < l; i++)
2547 75314 : if (rex[i])
2548 : {
2549 73300 : long t, ipr = F->subFB[i];
2550 235280 : for (t = 1; t <= n; t++)
2551 185596 : if (fact[t].pr == ipr) break;
2552 73300 : if (t > n) fact_update(R, F, ipr, cgx);
2553 : }
2554 : }
2555 : }
2556 2313941 : if (DEBUGLEVEL && Nfact) (*Nfact)++;
2557 2313941 : if (cache)
2558 : {
2559 : /* make sure we get maximal rank first, then allow all relations */
2560 2313941 : if (add_rel(cache, F, R, nz, gx, rex? 1: 0) <= 0)
2561 : { /* probably Q-dependent from previous ones: forget it */
2562 1783030 : if (DEBUGLEVEL>1) err_printf("*");
2563 1783030 : continue;
2564 : }
2565 530966 : if (cache->last >= cache->end) return 1; /* we have enough */
2566 : } else
2567 : {
2568 0 : gel(V,rel++) = gerepilecopy(av, mkvec3(R, stoi(nz), gx));
2569 8 : av = avma;
2570 : }
2571 431279 : if (++relid == Nrelid) break;
2572 : }
2573 166869 : END_Fincke_Pohst_ideal:
2574 166869 : return 0;
2575 : }
2576 :
2577 : static long
2578 266572 : Fincke_Pohst_ideal(RELCACHE_t *cache, FB_t *F, GEN nf, GEN I, GEN NI,
2579 : FACT *fact, long Nrelid, long max_FACT, FP_t *fp, GEN rex, long jid, long jid0, long e0,
2580 : long *Nsmall, long *Nfact)
2581 : {
2582 266572 : return Fincke_Pohst_ideal_both(cache, NULL,
2583 : F, nf, I, NI, fact, max_FACT, Nrelid, fp, rex, jid, jid0, e0, Nsmall, Nfact);
2584 : }
2585 :
2586 : static long
2587 0 : Fincke_Pohst_ideal_par(GEN V, FB_t *F, GEN nf, GEN I, GEN NI,
2588 : FACT *fact, long max_FACT, FP_t *fp, GEN rex, long jid, long jid0, long e0,
2589 : long *Nsmall, long *Nfact)
2590 : {
2591 0 : return Fincke_Pohst_ideal_both(NULL, V,
2592 : F, nf, I, NI, fact, max_FACT, max_FACT, fp, rex, jid, jid0, e0, Nsmall, Nfact);
2593 : }
2594 :
2595 : static GEN
2596 0 : pack_FB(FB_t *F, long s)
2597 : {
2598 0 : return mkvecn(s ? 8: 7, F->FB, F->LP, F->LV, F->iLP, F->idealperm, F->prodZ,
2599 : mkvecsmall3(F->KC,F->KCZ,F->KCZ2), F->subFB);
2600 : }
2601 :
2602 : static void
2603 0 : unpack_FB(FB_t *F, GEN P)
2604 : {
2605 0 : F->FB = gel(P,1);
2606 0 : F->LP = gel(P,2);
2607 0 : F->LV = gel(P,3);
2608 0 : F->iLP = gel(P,4);
2609 0 : F->idealperm = gel(P,5);
2610 0 : F->prodZ = gel(P,6);
2611 0 : F->KC = mael(P,7,1);
2612 0 : F->KCZ = mael(P,7,2);
2613 0 : F->KCZ2 = mael(P,7,3);
2614 0 : if (lg(P) > 8)
2615 0 : F->subFB = gel(P,8);
2616 0 : }
2617 :
2618 : GEN
2619 0 : bnfinit_FP_worker(GEN INI, GEN PF, GEN nf, long max_FACT, GEN rex, long jid0, long e0)
2620 : {
2621 0 : pari_sp av = avma;
2622 : FB_t F;
2623 : FP_t fp;
2624 : FACT * fact;
2625 0 : long Nsmall = 0, Nfact = 0, res, N = nf_get_degree(nf), jid = itos(gel(INI,3));
2626 0 : GEN vec = zerovec(max_FACT);
2627 0 : GEN I = gel(INI,1), NI = gel(INI,2);
2628 0 : if (isintzero(rex)) rex = NULL;
2629 0 : unpack_FB(&F, PF);
2630 0 : F.ballvol = ballvol(N);
2631 0 : minim_alloc(N+1, &fp.q, &fp.x, &fp.y, &fp.z, &fp.v);
2632 0 : fact = (FACT*)stack_malloc((F.KC+1)*sizeof(FACT));
2633 0 : res = Fincke_Pohst_ideal_par(vec, &F, nf, I, NI, fact, max_FACT, &fp, rex, jid, jid0, e0, &Nsmall, &Nfact);
2634 0 : return gerepilecopy(av, mkvec2(vec, mkvecsmall3(res, Nsmall, Nfact)));
2635 : }
2636 :
2637 : static void
2638 0 : small_norm_par(RELCACHE_t *cache, FB_t *F, GEN nf, long max_FACT, long idex, long nbthr, long j0)
2639 : {
2640 0 : GEN L_jid = F->L_jid, Np0 = NULL, p0 = j0? gel(F->LP,j0): NULL;
2641 0 : long Nsmall, Nfact, n = lg(L_jid)-1, e0 = 0;
2642 : pari_timer T;
2643 0 : long nt = nbthr? nbthr: mt_nbthreads();
2644 0 : GEN worker, vec = cgetg(nt+1, t_VEC);
2645 :
2646 0 : if (DEBUGLEVEL)
2647 : {
2648 0 : timer_start(&T);
2649 0 : err_printf("#### Look for %ld relations in %ld ideals (small_norm)\n",
2650 0 : cache->end - cache->last, lg(L_jid)-1);
2651 0 : if (p0) err_printf("Look in p0 = %Ps\n", vecslice(p0,1,4));
2652 : }
2653 0 : Nsmall = Nfact = 0;
2654 0 : if (p0)
2655 : {
2656 0 : GEN N0 = pr_norm(p0);
2657 0 : e0 = idex ? idex: logint0(sqri(pr_norm(veclast(F->LP))), N0, NULL);
2658 0 : p0 = idealpows(nf, p0, e0);
2659 0 : Np0 = powiu(N0, e0);
2660 : }
2661 0 : worker = snm_closure(is_entry("_bnfinit_FP_worker"),
2662 : mkcol6(pack_FB(F,0), nf, stoi(max_FACT), gen_0, stoi(j0), stoi(e0)));
2663 0 : while(n)
2664 : {
2665 : GEN VB;
2666 0 : long k, m = minss(n, nt);
2667 0 : for (k = 1; k <= m; k++, n--)
2668 : {
2669 0 : long j = L_jid[n];
2670 0 : GEN id = gel(F->LP,j), Nid;
2671 0 : if (p0)
2672 0 : { Nid = mulii(Np0, pr_norm(id)); id = idealmul(nf, p0, id); }
2673 : else
2674 0 : { Nid = pr_norm(id); id = pr_hnf(nf, id); }
2675 0 : gel(vec,k) = mkvec3(id, Nid, stoi(j));
2676 : }
2677 0 : setlg(vec,k);
2678 0 : VB = gen_parapply(worker,vec);
2679 0 : for (k = 1; k <= m; k++)
2680 : {
2681 0 : GEN B = gel(VB,k), B1 = gel(B,1), B2 = gel(B,2);
2682 0 : long i, lB = lg(B1);
2683 0 : Nsmall += B2[2]; Nfact += B2[3];
2684 0 : for (i = 1; i<lB && !isintzero(gel(B1,i)); i++)
2685 : {
2686 0 : GEN Bi = gel(B1,i), R = gel(Bi,1), gx = gel(Bi,3);
2687 0 : long nz = itos(gel(Bi,2));
2688 0 : if (cache->last < cache->end)
2689 0 : add_rel(cache, F, R, nz, gx, 0);
2690 : }
2691 0 : if (cache->last >= cache->end) { n = 0; break; }
2692 : }
2693 : }
2694 0 : if (DEBUGLEVEL && Nsmall)
2695 : {
2696 0 : if (DEBUGLEVEL == 1)
2697 0 : { if (Nfact) err_printf("\n"); }
2698 : else
2699 0 : err_printf(" \nnb. fact./nb. small norm = %ld/%ld = %.3f\n",
2700 0 : Nfact,Nsmall,((double)Nfact)/Nsmall);
2701 0 : if (timer_get(&T)>1) timer_printf(&T,"small_norm");
2702 : }
2703 0 : }
2704 :
2705 : static void
2706 66762 : small_norm_seq(RELCACHE_t *cache, FB_t *F, GEN nf, long Nrelid, long max_fact, long idex, FACT *fact, long j0)
2707 : {
2708 66762 : const long N = nf_get_degree(nf);
2709 : FP_t fp;
2710 : pari_sp av;
2711 66762 : GEN L_jid = F->L_jid, Np0 = NULL, p0 = j0? gel(F->LP,j0): NULL;
2712 66762 : long Nsmall, Nfact, n = lg(L_jid), e0 = 0;
2713 : pari_timer T;
2714 :
2715 66762 : if (DEBUGLEVEL)
2716 : {
2717 0 : timer_start(&T);
2718 0 : err_printf("#### Look for %ld relations in %ld ideals (small_norm)\n",
2719 0 : cache->end - cache->last, lg(L_jid)-1);
2720 0 : if (p0) err_printf("Look in p0 = %Ps\n", vecslice(p0,1,4));
2721 : }
2722 66762 : Nsmall = Nfact = 0;
2723 66762 : minim_alloc(N+1, &fp.q, &fp.x, &fp.y, &fp.z, &fp.v);
2724 66762 : if (p0)
2725 : {
2726 26701 : GEN n = pr_norm(p0);
2727 26701 : e0 = idex ? idex: logint0(sqri(pr_norm(veclast(F->LP))), n, NULL);
2728 26701 : p0 = idealpows(nf, p0, e0);
2729 26701 : Np0 = powiu(n,e0);
2730 : }
2731 165088 : for (av = avma; --n; set_avma(av))
2732 : {
2733 164670 : long j = L_jid[n];
2734 164670 : GEN id = gel(F->LP, j), Nid;
2735 164670 : if (DEBUGLEVEL>1)
2736 0 : err_printf("\n*** Ideal no %ld: %Ps\n", j, vecslice(id,1,4));
2737 164670 : if (p0)
2738 : {
2739 32492 : if (j == j0)
2740 : { /* avoid trivial relation */
2741 3894 : long e = pr_get_e(id);
2742 3894 : if ((e0 + 1) % e == 0 && e * pr_get_f(id) == N) continue;
2743 : }
2744 31849 : Nid = mulii(Np0, pr_norm(id)); id = idealmul(nf, p0, id);
2745 : }
2746 : else
2747 132178 : { Nid = pr_norm(id); id = pr_hnf(nf, id);}
2748 164028 : if (Fincke_Pohst_ideal(cache, F, nf, id, Nid, fact, Nrelid, max_fact, &fp,
2749 66346 : NULL, j, j0, e0, &Nsmall, &Nfact)) break;
2750 : }
2751 66764 : if (DEBUGLEVEL && Nsmall)
2752 : {
2753 0 : if (DEBUGLEVEL == 1)
2754 0 : { if (Nfact) err_printf("\n"); }
2755 : else
2756 0 : err_printf(" \nnb. fact./nb. small norm = %ld/%ld = %.3f\n",
2757 0 : Nfact,Nsmall,((double)Nfact)/Nsmall);
2758 0 : if (timer_get(&T)>1) timer_printf(&T,"small_norm");
2759 : }
2760 66764 : }
2761 :
2762 : static void
2763 66762 : small_norm(RELCACHE_t *cache, FB_t *F, GEN nf, long Nrelid, long max_fact, long idex, long nbthr, FACT *fact, long j0)
2764 : {
2765 66762 : if (nbthr==1)
2766 66762 : return small_norm_seq(cache, F, nf, Nrelid, max_fact, idex, fact, j0);
2767 : else
2768 0 : return small_norm_par(cache, F, nf, max_fact, idex, nbthr, j0);
2769 : }
2770 :
2771 : static GEN
2772 55612 : get_random_ideal(FB_t *F, GEN nf, GEN ex)
2773 : {
2774 55612 : long i, l = lg(ex);
2775 : for (;;)
2776 1033 : {
2777 56645 : GEN I = NULL;
2778 158690 : for (i = 1; i < l; i++)
2779 102046 : if ((ex[i] = random_bits(RANDOM_BITS)))
2780 : {
2781 97766 : GEN pr = gel(F->LP, F->subFB[i]), e = utoipos(ex[i]);
2782 97764 : I = I? idealmulpowprime(nf, I, pr, e): idealpow(nf, pr, e);
2783 : }
2784 56644 : if (I && !ZM_isscalar(I,NULL)) return I; /* != (n)Z_K */
2785 : }
2786 : }
2787 :
2788 : static void
2789 0 : rnd_rel_par(RELCACHE_t *cache, FB_t *F, GEN nf, long max_FACT)
2790 : {
2791 : pari_timer T;
2792 0 : GEN L_jid = F->L_jid, R, ex;
2793 0 : long k, l = lg(L_jid), Nfact = 0, Nsmall = 0;
2794 : GEN worker, vec, VB, NR;
2795 :
2796 0 : if (DEBUGLEVEL) {
2797 0 : timer_start(&T);
2798 0 : err_printf("#### Look for %ld relations in %ld ideals (rnd_rel)\n",
2799 0 : cache->end - cache->last, l-1);
2800 : }
2801 0 : ex = cgetg(lg(F->subFB), t_VECSMALL);
2802 0 : R = get_random_ideal(F, nf, ex); /* random product from subFB */
2803 0 : NR = ZM_det_triangular(R);
2804 0 : worker = snm_closure(is_entry("_bnfinit_FP_worker"),
2805 : mkcol6(pack_FB(F,1), nf, stoi(max_FACT), ex, gen_0, gen_0));
2806 0 : vec = cgetg(l, t_VEC);
2807 0 : for (k = 1; k < l; k++)
2808 : {
2809 0 : long j = L_jid[k];
2810 0 : GEN id = gel(F->LP,j), Nid;
2811 0 : Nid = mulii(NR, pr_norm(id)); id = idealmul(nf, R, id);
2812 0 : gel(vec,k) = mkvec3(id, Nid, stoi(j));
2813 : }
2814 0 : VB = gen_parapply(worker,vec);
2815 0 : for (k = 1; k < l; k++)
2816 : {
2817 0 : GEN B = gel(VB,k), B1 = gel(B,1), B2 = gel(B,2);
2818 0 : long i, lB = lg(B1);
2819 0 : Nsmall += B2[2]; Nfact += B2[3];
2820 0 : for (i = 1; i<lB && !isintzero(gel(B1,i)); i++)
2821 : {
2822 0 : GEN Bi = gel(B1,i), R = gel(Bi,1), gx = gel(Bi,3);
2823 0 : long nz = itos(gel(Bi,2));
2824 0 : if (cache->last < cache->end)
2825 0 : add_rel(cache, F, R, nz, gx, 1);
2826 : }
2827 0 : if (cache->last >= cache->end) break;
2828 : }
2829 :
2830 0 : if (DEBUGLEVEL)
2831 : {
2832 0 : if (DEBUGLEVEL == 1)
2833 0 : { if (Nfact) err_printf("\n"); }
2834 : else
2835 0 : err_printf(" \nnb. fact./nb. small norm = %ld/%ld = %.3f\n",
2836 0 : Nfact,Nsmall,((double)Nfact)/Nsmall);
2837 0 : if (timer_get(&T)>=0) timer_printf(&T,"rnd_rel");
2838 : }
2839 0 : }
2840 :
2841 : static void
2842 55612 : rnd_rel_seq(RELCACHE_t *cache, FB_t *F, GEN nf, long max_fact, FACT *fact)
2843 : {
2844 : pari_timer T;
2845 55612 : GEN L_jid = F->L_jid, R, NR, ex;
2846 55612 : long i, l = lg(L_jid), Nfact = 0;
2847 : FP_t fp;
2848 : pari_sp av;
2849 :
2850 55612 : if (DEBUGLEVEL) {
2851 0 : timer_start(&T);
2852 0 : err_printf("#### Look for %ld relations in %ld ideals (rnd_rel)\n",
2853 0 : cache->end - cache->last, l-1);
2854 : }
2855 55612 : ex = cgetg(lg(F->subFB), t_VECSMALL);
2856 55612 : R = get_random_ideal(F, nf, ex); /* random product from subFB */
2857 55611 : NR = ZM_det_triangular(R);
2858 55608 : minim_alloc(nf_get_degree(nf)+1, &fp.q, &fp.x, &fp.y, &fp.z, &fp.v);
2859 124797 : for (av = avma, i = 1; i < l; i++, set_avma(av))
2860 : { /* try P[j] * base */
2861 102533 : long j = L_jid[i];
2862 102533 : GEN P = gel(F->LP, j), Nid = mulii(NR, pr_norm(P));
2863 102520 : if (DEBUGLEVEL>1) err_printf("\n*** Ideal %ld: %Ps\n", j, vecslice(P,1,4));
2864 102520 : if (Fincke_Pohst_ideal(cache, F, nf, idealHNF_mul(nf, R, P), Nid, fact,
2865 33349 : RND_REL_RELPID, max_fact, &fp, ex, j, 0, 0, NULL, &Nfact)) break;
2866 : }
2867 55613 : if (DEBUGLEVEL)
2868 : {
2869 0 : if (Nfact) err_printf("\n");
2870 0 : if (timer_get(&T)>=0) timer_printf(&T,"rnd_rel");
2871 : }
2872 55613 : }
2873 : static void
2874 55612 : rnd_rel(RELCACHE_t *cache, FB_t *F, GEN nf, long max_fact, long nbthr, FACT *fact)
2875 : {
2876 55612 : if (nbthr==1)
2877 55612 : return rnd_rel_seq(cache, F, nf, max_fact, fact);
2878 : else
2879 0 : return rnd_rel_par(cache, F, nf, max_fact);
2880 : }
2881 :
2882 : static GEN
2883 64000 : automorphism_perms(GEN M, GEN auts, GEN cyclic, long r1, long r2, long N)
2884 : {
2885 64000 : long L = lgcols(M), lauts = lg(auts), lcyc = lg(cyclic), i, j, l, m;
2886 64000 : GEN Mt, perms = cgetg(lauts, t_VEC);
2887 : pari_sp av;
2888 :
2889 128833 : for (l = 1; l < lauts; l++) gel(perms, l) = cgetg(L, t_VECSMALL);
2890 64000 : av = avma;
2891 64000 : Mt = shallowtrans(gprec_w(M, LOWDEFAULTPREC));
2892 64000 : Mt = shallowconcat(Mt, conj_i(vecslice(Mt, r1+1, r1+r2)));
2893 111450 : for (l = 1; l < lcyc; l++)
2894 : {
2895 47450 : GEN thiscyc = gel(cyclic, l), thisperm, perm, prev, Nt;
2896 47450 : long k = thiscyc[1];
2897 :
2898 47450 : Nt = RgM_mul(shallowtrans(gel(auts, k)), Mt);
2899 47451 : perm = gel(perms, k);
2900 157262 : for (i = 1; i < L; i++)
2901 : {
2902 109810 : GEN v = gel(Nt, i), minD;
2903 109810 : minD = gnorml2(gsub(v, gel(Mt, 1)));
2904 109808 : perm[i] = 1;
2905 577162 : for (j = 2; j <= N; j++)
2906 : {
2907 467351 : GEN D = gnorml2(gsub(v, gel(Mt, j)));
2908 467352 : if (gcmp(D, minD) < 0) { minD = D; perm[i] = j >= L ? r2-j : j; }
2909 : }
2910 : }
2911 66065 : for (prev = perm, m = 2; m < lg(thiscyc); m++, prev = thisperm)
2912 : {
2913 18613 : thisperm = gel(perms, thiscyc[m]);
2914 94360 : for (i = 1; i < L; i++)
2915 : {
2916 75747 : long pp = labs(prev[i]);
2917 75747 : thisperm[i] = prev[i] < 0 ? -perm[pp] : perm[pp];
2918 : }
2919 : }
2920 : }
2921 64000 : set_avma(av); return perms;
2922 : }
2923 :
2924 : /* Determine the field automorphisms as matrices on the integral basis */
2925 : static GEN
2926 64061 : automorphism_matrices(GEN nf, GEN *cycp)
2927 : {
2928 64061 : pari_sp av = avma;
2929 64061 : GEN auts = galoisconj(nf, NULL), mats, cyclic, cyclicidx;
2930 64062 : long nauts = lg(auts)-1, i, j, k, l;
2931 :
2932 64062 : cyclic = cgetg(nauts+1, t_VEC);
2933 64062 : cyclicidx = zero_Flv(nauts);
2934 98478 : for (l = 1; l <= nauts; l++)
2935 : {
2936 98478 : GEN aut = gel(auts, l);
2937 98478 : if (gequalX(aut)) { swap(gel(auts, l), gel(auts, nauts)); break; }
2938 : }
2939 : /* trivial automorphism is last */
2940 192980 : for (l = 1; l <= nauts; l++) gel(auts, l) = algtobasis(nf, gel(auts, l));
2941 : /* Compute maximal cyclic subgroups */
2942 128920 : for (l = nauts; --l > 0; ) if (!cyclicidx[l])
2943 : {
2944 48956 : GEN elt = gel(auts, l), aut = elt, cyc = cgetg(nauts+1, t_VECSMALL);
2945 48956 : cyc[1] = cyclicidx[l] = l; j = 1;
2946 : do
2947 : {
2948 68112 : elt = galoisapply(nf, elt, aut);
2949 221626 : for (k = 1; k <= nauts; k++) if (gequal(elt, gel(auts, k))) break;
2950 68111 : cyclicidx[k] = l; cyc[++j] = k;
2951 : }
2952 68111 : while (k != nauts);
2953 48955 : setlg(cyc, j);
2954 48955 : gel(cyclic, l) = cyc;
2955 : }
2956 128919 : for (i = j = 1; i < nauts; i++)
2957 64859 : if (cyclicidx[i] == i) cyclic[j++] = cyclic[i];
2958 64060 : setlg(cyclic, j);
2959 64061 : mats = cgetg(nauts, t_VEC);
2960 111541 : while (--j > 0)
2961 : {
2962 47478 : GEN cyc = gel(cyclic, j);
2963 47478 : long id = cyc[1];
2964 47478 : GEN M, Mi, aut = gel(auts, id);
2965 :
2966 47478 : gel(mats, id) = Mi = M = nfgaloismatrix(nf, aut);
2967 66093 : for (i = 2; i < lg(cyc); i++) gel(mats, cyc[i]) = Mi = ZM_mul(Mi, M);
2968 : }
2969 64063 : (void)gc_all(av, 2, &mats, &cyclic);
2970 64063 : if (cycp) *cycp = cyclic;
2971 64063 : return mats;
2972 : }
2973 :
2974 : /* vP a list of maximal ideals above the same p from idealprimedec: f(P/p) is
2975 : * increasing; 1 <= j <= #vP; orbit a zc of length <= #vP; auts a vector of
2976 : * automorphisms in ZM form.
2977 : * Set orbit[i] = 1 for all vP[i], i >= j, in the orbit of pr = vP[j] wrt auts.
2978 : * N.B.1 orbit need not be initialized to 0: useful to incrementally run
2979 : * through successive orbits
2980 : * N.B.2 i >= j, so primes with index < j will be missed; run incrementally
2981 : * starting from j = 1 ! */
2982 : static void
2983 11865 : pr_orbit_fill(GEN orbit, GEN auts, GEN vP, long j)
2984 : {
2985 11865 : GEN pr = gel(vP,j), gen = pr_get_gen(pr);
2986 11865 : long i, l = lg(auts), J = lg(orbit), f = pr_get_f(pr);
2987 11865 : orbit[j] = 1;
2988 23730 : for (i = 1; i < l; i++)
2989 : {
2990 11865 : GEN g = ZM_ZC_mul(gel(auts,i), gen);
2991 : long k;
2992 11886 : for (k = j+1; k < J; k++)
2993 : {
2994 35 : GEN prk = gel(vP,k);
2995 35 : if (pr_get_f(prk) > f) break; /* f(P[k]) increases with k */
2996 : /* don't check that e matches: (almost) always 1 ! */
2997 35 : if (!orbit[k] && ZC_prdvd(g, prk)) { orbit[k] = 1; break; }
2998 : }
2999 : }
3000 11865 : }
3001 : /* remark: F->KCZ changes if be_honest() fails */
3002 : static int
3003 7 : be_honest(FB_t *F, GEN nf, GEN auts, FACT *fact)
3004 : {
3005 : long i, iz, nbtest;
3006 7 : long lgsub = lg(F->subFB), KCZ0 = F->KCZ, N = nf_get_degree(nf);
3007 : FP_t fp;
3008 : pari_sp av;
3009 :
3010 7 : if (DEBUGLEVEL) {
3011 0 : err_printf("Be honest for %ld primes from %ld to %ld\n", F->KCZ2 - F->KCZ,
3012 0 : F->FB[ F->KCZ+1 ], F->FB[ F->KCZ2 ]);
3013 : }
3014 7 : minim_alloc(N+1, &fp.q, &fp.x, &fp.y, &fp.z, &fp.v);
3015 7 : if (lg(auts) == 1) auts = NULL;
3016 7 : av = avma;
3017 14 : for (iz=F->KCZ+1; iz<=F->KCZ2; iz++, set_avma(av))
3018 : {
3019 7 : long p = F->FB[iz];
3020 7 : GEN pr_orbit, P = gel(F->LV,p);
3021 7 : long j, J = lg(P); /* > 1 */
3022 : /* the P|p, NP > C2 are assumed in subgroup generated by FB + last P
3023 : * with NP <= C2 is unramified --> check all but last */
3024 7 : if (pr_get_e(gel(P,J-1)) == 1) J--;
3025 7 : if (J == 1) continue;
3026 7 : if (DEBUGLEVEL>1) err_printf("%ld ", p);
3027 7 : pr_orbit = auts? zero_zv(J-1): NULL;
3028 28 : for (j = 1; j < J; j++)
3029 : {
3030 : GEN Nid, id, id0;
3031 21 : if (pr_orbit)
3032 : {
3033 21 : if (pr_orbit[j]) continue;
3034 : /* discard all primes in automorphism orbit simultaneously */
3035 14 : pr_orbit_fill(pr_orbit, auts, P, j);
3036 : }
3037 14 : id = id0 = pr_hnf(nf,gel(P,j));
3038 14 : Nid = pr_norm(gel(P,j));
3039 14 : for (nbtest=0;;)
3040 : {
3041 14 : if (Fincke_Pohst_ideal(NULL, F, nf, id, Nid, fact, 0, MAXTRY_FACT, &fp,
3042 14 : NULL, 0, 0, 0, NULL, NULL)) break;
3043 0 : if (++nbtest > maxtry_HONEST)
3044 : {
3045 0 : if (DEBUGLEVEL)
3046 0 : pari_warn(warner,"be_honest() failure on prime %Ps\n", gel(P,j));
3047 0 : return 0;
3048 : }
3049 : /* occurs at most once in the whole function */
3050 0 : for (i = 1, id = id0; i < lgsub; i++)
3051 : {
3052 0 : long ex = random_bits(RANDOM_BITS);
3053 0 : if (ex)
3054 : {
3055 0 : GEN pr = gel(F->LP, F->subFB[i]);
3056 0 : id = idealmulpowprime(nf, id, pr, utoipos(ex));
3057 : }
3058 : }
3059 0 : if (!equali1(gcoeff(id,N,N))) id = Q_primpart(id);
3060 0 : if (expi(gcoeff(id,1,1)) > 100) id = idealred(nf, id);
3061 0 : Nid = ZM_det_triangular(id);
3062 : }
3063 : }
3064 7 : F->KCZ++; /* SUCCESS, "enlarge" factorbase */
3065 : }
3066 7 : F->KCZ = KCZ0; return gc_bool(av,1);
3067 : }
3068 :
3069 : /* all primes with N(P) <= BOUND factor on factorbase ? */
3070 : void
3071 63 : bnftestprimes(GEN bnf, GEN BOUND)
3072 : {
3073 63 : pari_sp av0 = avma, av;
3074 63 : ulong count = 0;
3075 63 : GEN auts, p, nf = bnf_get_nf(bnf), Vbase = bnf_get_vbase(bnf);
3076 63 : GEN fb = gen_sort_shallow(Vbase, (void*)&cmp_prime_ideal, cmp_nodata);
3077 63 : ulong pmax = pr_get_smallp(veclast(fb)); /*largest p in factorbase*/
3078 : forprime_t S;
3079 : FACT *fact;
3080 : FB_t F;
3081 :
3082 63 : (void)recover_partFB(&F, Vbase, nf_get_degree(nf));
3083 63 : fact = (FACT*)stack_malloc((F.KC+1)*sizeof(FACT));
3084 63 : forprime_init(&S, gen_2, BOUND);
3085 63 : auts = automorphism_matrices(nf, NULL);
3086 63 : if (lg(auts) == 1) auts = NULL;
3087 63 : av = avma;
3088 37604 : while (( p = forprime_next(&S) ))
3089 : {
3090 : GEN pr_orbit, vP;
3091 : long j, J;
3092 :
3093 37541 : if (DEBUGLEVEL == 1 && ++count > 1000)
3094 : {
3095 0 : err_printf("passing p = %Ps / %Ps\n", p, BOUND);
3096 0 : count = 0;
3097 : }
3098 37541 : set_avma(av);
3099 37541 : vP = idealprimedec_limit_norm(nf, p, BOUND);
3100 37541 : J = lg(vP);
3101 : /* if last is unramified, all P|p in subgroup generated by FB: skip last */
3102 37541 : if (J > 1 && pr_get_e(gel(vP,J-1)) == 1) J--;
3103 37541 : if (J == 1) continue;
3104 14525 : if (DEBUGLEVEL>1) err_printf("*** p = %Ps\n",p);
3105 14525 : pr_orbit = auts? zero_zv(J-1): NULL;
3106 31549 : for (j = 1; j < J; j++)
3107 : {
3108 17024 : GEN P = gel(vP,j);
3109 17024 : long k = 0;
3110 17024 : if (pr_orbit)
3111 : {
3112 11858 : if (pr_orbit[j]) continue;
3113 : /* discard all primes in automorphism orbit simultaneously */
3114 11851 : pr_orbit_fill(pr_orbit, auts, vP, j);
3115 : }
3116 17017 : if (abscmpiu(p, pmax) > 0 || !(k = tablesearch(fb, P, &cmp_prime_ideal)))
3117 16408 : (void)SPLIT(&F, nf, pr_hnf(nf,P), Vbase, fact);
3118 17017 : if (DEBUGLEVEL>1)
3119 : {
3120 0 : err_printf(" Testing P = %Ps\n",P);
3121 0 : if (k) err_printf(" #%ld in factor base\n",k);
3122 0 : else err_printf(" is %Ps\n", isprincipal(bnf,P));
3123 : }
3124 : }
3125 : }
3126 63 : set_avma(av0);
3127 63 : }
3128 :
3129 : /* A t_MAT of complex floats, in fact reals. Extract a submatrix B
3130 : * whose columns are definitely nonzero, i.e. gexpo(A[j]) >= -2
3131 : *
3132 : * If possible precision problem (t_REAL 0 with large exponent), set
3133 : * *precpb to 1 */
3134 : static GEN
3135 93032 : clean_cols(GEN A, int *precpb)
3136 : {
3137 93032 : long l = lg(A), h, i, j, k;
3138 : GEN B;
3139 93032 : *precpb = 0;
3140 93032 : if (l == 1) return A;
3141 93032 : h = lgcols(A);;
3142 93032 : B = cgetg(l, t_MAT);
3143 1092750 : for (i = k = 1; i < l; i++)
3144 : {
3145 999718 : GEN Ai = gel(A,i);
3146 999718 : int non0 = 0;
3147 4356649 : for (j = 1; j < h; j++)
3148 : {
3149 3356931 : GEN c = gel(Ai,j);
3150 3356931 : if (gexpo(c) >= -2)
3151 : {
3152 2017454 : if (gequal0(c)) *precpb = 1; else non0 = 1;
3153 : }
3154 : }
3155 999718 : if (non0) gel(B, k++) = Ai;
3156 : }
3157 93032 : setlg(B, k); return B;
3158 : }
3159 :
3160 : static long
3161 575962 : compute_multiple_of_R_pivot(GEN X, GEN x0/*unused*/, long ix, GEN c)
3162 : {
3163 575962 : GEN x = gel(X,ix);
3164 575962 : long i, k = 0, ex = - (long)HIGHEXPOBIT, lx = lg(x);
3165 : (void)x0;
3166 2836227 : for (i=1; i<lx; i++)
3167 2260280 : if (!c[i] && !gequal0(gel(x,i)))
3168 : {
3169 729598 : long e = gexpo(gel(x,i));
3170 729584 : if (e > ex) { ex = e; k = i; }
3171 : }
3172 575947 : return (k && ex > -32)? k: lx;
3173 : }
3174 :
3175 : /* Ar = (log |sigma_i(u_j)|) for units (u_j) found so far;
3176 : * RU = R1+R2 = target rank for unit matrix, after adding [1 x r1, 2 x r2];
3177 : * N = field degree, need = unit rank defect;
3178 : * L = NULL (prec problem) or B^(-1) * A with approximate rational entries
3179 : * (as t_REAL), B a submatrix of A, with (probably) maximal rank RU */
3180 : static GEN
3181 121719 : compute_multiple_of_R(GEN Ar, long RU, long N, long *pneed, long *bit, GEN *ptL)
3182 : {
3183 : GEN T, d, mdet, Im_mdet, kR, L;
3184 121719 : long i, j, r, R1 = 2*RU - N;
3185 : int precpb;
3186 121719 : pari_sp av = avma;
3187 :
3188 121719 : if (RU == 1) { *ptL = zeromat(0, lg(Ar)-1); return gen_1; }
3189 :
3190 93032 : if (DEBUGLEVEL) err_printf("\n#### Computing regulator multiple\n");
3191 93032 : mdet = clean_cols(Ar, &precpb);
3192 : /* will cause precision to increase on later failure, but we may succeed! */
3193 93032 : *ptL = precpb? NULL: gen_1;
3194 93032 : T = cgetg(RU+1,t_COL);
3195 235149 : for (i=1; i<=R1; i++) gel(T,i) = gen_1;
3196 193194 : for ( ; i<=RU; i++) gel(T,i) = gen_2;
3197 93032 : mdet = shallowconcat(T, mdet); /* det(Span(mdet)) = N * R */
3198 :
3199 : /* could be using indexrank(), but need custom "get_pivot" function */
3200 93032 : d = RgM_pivots(mdet, &r, &compute_multiple_of_R_pivot, NULL);
3201 : /* # of independent columns = target rank ? */
3202 93032 : if (lg(mdet)-1 - r != RU)
3203 : {
3204 25361 : if (DEBUGLEVEL)
3205 0 : err_printf("Units matrix target rank = %ld < %ld\n",lg(mdet)-1 - r, RU);
3206 25361 : *pneed = RU - (lg(mdet)-1-r); return gc_NULL(av);
3207 : }
3208 :
3209 67671 : Im_mdet = cgetg(RU+1, t_MAT); /* extract independent columns */
3210 : /* N.B: d[1] = 1, corresponding to T above */
3211 67671 : gel(Im_mdet, 1) = T;
3212 257409 : for (i = j = 2; i <= RU; j++)
3213 189738 : if (d[j]) gel(Im_mdet, i++) = gel(mdet,j);
3214 :
3215 : /* integral multiple of R: the cols we picked form a Q-basis, they have an
3216 : * index in the full lattice. First column is T */
3217 67671 : kR = divru(det2(Im_mdet), N);
3218 : /* R > 0.2 uniformly */
3219 67671 : if (!signe(kR) || expo(kR) < -3)
3220 : {
3221 0 : if (DEBUGLEVEL) err_printf("Regulator is zero.\n");
3222 0 : *pneed = 0; return gc_NULL(av);
3223 : }
3224 67671 : d = det2(rowslice(vecslice(Im_mdet, 2, RU), 2, RU));
3225 67671 : setabssign(d); setabssign(kR);
3226 67671 : if (gexpo(gsub(d,kR)) - gexpo(d) > -20) { *ptL = NULL; return gc_NULL(av); }
3227 67670 : L = RgM_inv(Im_mdet);
3228 : /* estimate # of correct bits in result */
3229 67669 : if (!L || (*bit = -gexpo(RgM_Rg_sub_shallow(RgM_mul(L,Im_mdet), gen_1))) < 16)
3230 10 : { *ptL = NULL; return gc_NULL(av); }
3231 :
3232 67660 : *ptL = RgM_mul(rowslice(L,2,RU), Ar); /* approximate rational entries */
3233 67660 : return gc_all(av,2, &kR, ptL);
3234 : }
3235 :
3236 : /* leave small integer n as is, convert huge n to t_REAL (for readability) */
3237 : static GEN
3238 0 : i2print(GEN n)
3239 0 : { return lgefint(n) <= DEFAULTPREC? n: itor(n,LOWDEFAULTPREC); }
3240 :
3241 : static long
3242 96198 : bad_check(GEN c)
3243 : {
3244 96198 : long ec = gexpo(c);
3245 96198 : if (DEBUGLEVEL) err_printf("\n ***** check = %.28Pg\n",c);
3246 : /* safe check for c < 0.75 : avoid underflow in gtodouble() */
3247 96198 : if (ec < -1 || (ec == -1 && gtodouble(c) < 0.75)) return fupb_PRECI;
3248 : /* safe check for c > 1.3 : avoid overflow */
3249 96198 : if (ec > 0 || (ec == 0 && gtodouble(c) > 1.3)) return fupb_RELAT;
3250 64006 : return fupb_NONE;
3251 : }
3252 : /* Input:
3253 : * lambda = approximate rational entries: coords of units found so far on a
3254 : * sublattice of maximal rank (sublambda)
3255 : * *ptkR = regulator of sublambda = multiple of regulator of lambda
3256 : * Compute R = true regulator of lambda.
3257 : *
3258 : * If c := Rz ~ 1, by Dirichlet's formula, then lambda is the full group of
3259 : * units AND the full set of relations for the class group has been computed.
3260 : * In fact z is a very rough approximation and we only expect 0.75 < Rz < 1.3
3261 : *
3262 : * Output: *ptkR = R, *ptL = numerator(units) (in terms of lambda) */
3263 : static long
3264 96260 : compute_R(GEN lambda, GEN z, GEN *ptL, GEN *ptkR)
3265 : {
3266 96260 : pari_sp av = avma;
3267 96260 : long bit, r, reason, RU = lg(lambda) == 1? 1: lgcols(lambda);
3268 : GEN L, H, D, den, R, c;
3269 :
3270 96261 : *ptL = NULL;
3271 96261 : if (RU == 1) { *ptkR = gen_1; *ptL = lambda; return bad_check(z); }
3272 67575 : D = gmul2n(mpmul(*ptkR,z), 1); /* bound for denom(lambda) */
3273 67576 : if (expo(D) < 0 && rtodbl(D) < 0.95) return fupb_PRECI;
3274 67576 : L = bestappr(lambda,D);
3275 67576 : if (lg(L) == 1)
3276 : {
3277 0 : if (DEBUGLEVEL) err_printf("truncation error in bestappr\n");
3278 0 : return fupb_PRECI;
3279 : }
3280 67576 : den = Q_denom(L);
3281 67576 : if (mpcmp(den,D) > 0)
3282 : {
3283 15 : if (DEBUGLEVEL) err_printf("D = %Ps\nden = %Ps\n",D, i2print(den));
3284 15 : return fupb_PRECI;
3285 : }
3286 67561 : bit = -gexpo(gsub(L, lambda)); /* input accuracy */
3287 67561 : L = Q_muli_to_int(L, den);
3288 67558 : if (gexpo(L) + expi(den) > bit - 32)
3289 : {
3290 47 : if (DEBUGLEVEL) err_printf("dubious bestappr; den = %Ps\n", i2print(den));
3291 47 : return fupb_PRECI;
3292 : }
3293 67512 : H = ZM_hnf(L); r = lg(H)-1;
3294 67513 : if (!r || r != nbrows(H))
3295 0 : R = gen_0; /* wrong rank */
3296 : else
3297 67513 : R = gmul(*ptkR, gdiv(ZM_det_triangular(H), powiu(den, r)));
3298 : /* R = tentative regulator; regulator > 0.2 uniformly */
3299 67513 : if (gexpo(R) < -3) {
3300 0 : if (DEBUGLEVEL) err_printf("\n#### Tentative regulator: %.28Pg\n", R);
3301 0 : return gc_long(av, fupb_PRECI);
3302 : }
3303 67513 : c = gmul(R,z); /* should be n (= 1 if we are done) */
3304 67513 : if (DEBUGLEVEL) err_printf("\n#### Tentative regulator: %.28Pg\n", R);
3305 67513 : if ((reason = bad_check(c))) return gc_long(av, reason);
3306 48796 : *ptkR = R; *ptL = L; return fupb_NONE;
3307 : }
3308 : static GEN
3309 64096 : get_clg2(GEN cyc, GEN Ga, GEN C, GEN Ur, GEN Ge, GEN M1, GEN M2)
3310 : {
3311 64096 : GEN GD = gsub(act_arch(M1, C), diagact_arch(cyc, Ga));
3312 64097 : GEN ga = gsub(act_arch(M2, C), act_arch(Ur, Ga));
3313 64098 : return mkvecn(6, Ur, ga, GD, Ge, M1, M2);
3314 : }
3315 : /* compute class group (clg1) + data for isprincipal (clg2) */
3316 : static GEN
3317 64001 : class_group_gen(GEN nf,GEN W,GEN C,GEN Vbase,long prec, GEN *pclg2)
3318 : {
3319 : GEN M1, M2, z, G, Ga, Ge, cyc, X, Y, D, U, V, Ur, Ui, Uir;
3320 : long j, l;
3321 :
3322 64001 : D = ZM_snfall(W,&U,&V); /* UWV=D, D diagonal, G = g Ui (G=new gens, g=old) */
3323 64001 : Ui = ZM_inv(U, NULL);
3324 64001 : l = lg(D); cyc = cgetg(l, t_VEC); /* elementary divisors */
3325 92785 : for (j = 1; j < l; j++)
3326 : {
3327 30366 : gel(cyc,j) = gcoeff(D,j,j); /* strip useless components */
3328 30366 : if (is_pm1(gel(cyc,j))) break;
3329 : }
3330 64001 : l = j;
3331 64001 : Ur = ZM_hnfdivrem(U, D, &Y);
3332 64001 : Uir = ZM_hnfdivrem(Ui,W, &X);
3333 : /* {x} = logarithmic embedding of x (arch. component)
3334 : * NB: [J,z] = idealred(I) --> I = y J, with {y} = - z
3335 : * G = g Uir - {Ga}, Uir = Ui - WX
3336 : * g = G Ur - {ga}, Ur = U - DY */
3337 64001 : G = cgetg(l,t_VEC);
3338 64001 : Ga= cgetg(l,t_MAT);
3339 64001 : Ge= cgetg(l,t_COL);
3340 64001 : z = init_famat(NULL);
3341 92784 : for (j = 1; j < l; j++)
3342 : {
3343 28784 : GEN I = genback(z, nf, Vbase, gel(Uir,j));
3344 28784 : gel(G,j) = gel(I,1); /* generator, order cyc[j] */
3345 28784 : gel(Ge,j)= gel(I,2);
3346 28784 : gel(Ga,j)= nf_cxlog(nf, gel(I,2), prec);
3347 28784 : if (!gel(Ga,j)) pari_err_PREC("class_group_gen");
3348 : }
3349 : /* {ga} = - {GD}Y + G U - g = - {GD}Y - {Ga} U - gW X U
3350 : = - gW (X Ur + V Y) - {Ga}Ur */
3351 64000 : M2 = ZM_neg(ZM_add(ZM_mul(X,Ur), ZM_mul(V,Y)));
3352 64000 : setlg(cyc,l); setlg(V,l); setlg(D,l);
3353 : /* G D =: {GD} = g (Ui - W X) D - {Ga}D = g W (V - X D) - {Ga}D
3354 : * NB: Ui D = W V. gW is given by (first l-1 cols of) C */
3355 64000 : M1 = ZM_sub(V, ZM_mul(X,D));
3356 63999 : *pclg2 = get_clg2(cyc, Ga, C, Ur, Ge, M1, M2);
3357 64001 : return mkvec3(ZV_prod(cyc), cyc, G);
3358 : }
3359 :
3360 : /* compute principal ideals corresponding to (gen[i]^cyc[i]) */
3361 : static GEN
3362 4956 : makecycgen(GEN bnf)
3363 : {
3364 4956 : GEN cyc = bnf_get_cyc(bnf), gen = bnf_get_gen(bnf), nf = bnf_get_nf(bnf);
3365 4956 : GEN h, y, GD = bnf_get_GD(bnf), W = bnf_get_W(bnf); /* HNF */
3366 4956 : GEN SUnits = bnf_get_sunits(bnf);
3367 4956 : GEN X = SUnits? gel(SUnits,1): NULL, C = SUnits? gel(SUnits,3): NULL;
3368 : long e, i, l;
3369 :
3370 4956 : if (DEBUGLEVEL) pari_warn(warner,"completing bnf (building cycgen)");
3371 4956 : h = cgetg_copy(gen, &l);
3372 11613 : for (i = 1; i < l; i++)
3373 : {
3374 6657 : GEN gi = gel(gen,i), ci = gel(cyc,i);
3375 6657 : if (X && equalii(ci, gcoeff(W,i,i)))
3376 : {
3377 : long j;
3378 8589 : for (j = i+1; j < l; j++)
3379 3213 : if (signe(gcoeff(W,i,j))) break;
3380 5550 : if (j == i) { gel(h,i) = mkmat2(X, gel(C,i)); continue; }
3381 : }
3382 6657 : if (abscmpiu(ci, 5) < 0)
3383 : {
3384 5544 : GEN N = ZM_det_triangular(gi);
3385 5544 : y = isprincipalarch(bnf,gel(GD,i), N, ci, gen_1, &e);
3386 5544 : if (y && fact_ok(nf,y,NULL,mkvec(gi),mkvec(ci)))
3387 : {
3388 4556 : gel(h,i) = to_famat_shallow(y,gen_1);
3389 4556 : continue;
3390 : }
3391 : }
3392 2101 : y = isprincipalfact(bnf, NULL, mkvec(gi), mkvec(ci), nf_GENMAT|nf_FORCE);
3393 2101 : gel(h,i) = gel(y,2);
3394 : }
3395 4956 : return h;
3396 : }
3397 :
3398 : static GEN
3399 69 : get_y(GEN bnf, GEN W, GEN B, GEN C, GEN pFB, long j)
3400 : {
3401 69 : GEN y, nf = bnf_get_nf(bnf);
3402 69 : long e, lW = lg(W)-1;
3403 69 : GEN ex = (j<=lW)? gel(W,j): gel(B,j-lW);
3404 69 : GEN P = (j<=lW)? NULL: gel(pFB,j);
3405 69 : if (C)
3406 : { /* archimedean embeddings known: cheap trial */
3407 69 : GEN Nx = get_norm_fact_primes(pFB, ex, P);
3408 69 : y = isprincipalarch(bnf,gel(C,j), Nx,gen_1, gen_1, &e);
3409 69 : if (y && fact_ok(nf,y,P,pFB,ex)) return y;
3410 : }
3411 0 : y = isprincipalfact_or_fail(bnf, P, pFB, ex);
3412 0 : return typ(y) == t_INT? y: gel(y,2);
3413 : }
3414 : /* compute principal ideals corresponding to bnf relations */
3415 : static GEN
3416 20 : makematal(GEN bnf)
3417 : {
3418 20 : GEN W = bnf_get_W(bnf), B = bnf_get_B(bnf), C = bnf_get_C(bnf);
3419 : GEN pFB, ma, retry;
3420 20 : long lma, j, prec = 0;
3421 :
3422 20 : if (DEBUGLEVEL) pari_warn(warner,"completing bnf (building matal)");
3423 20 : lma=lg(W)+lg(B)-1;
3424 20 : pFB = bnf_get_vbase(bnf);
3425 20 : ma = cgetg(lma,t_VEC);
3426 20 : retry = vecsmalltrunc_init(lma);
3427 89 : for (j=lma-1; j>0; j--)
3428 : {
3429 69 : pari_sp av = avma;
3430 69 : GEN y = get_y(bnf, W, B, C, pFB, j);
3431 69 : if (typ(y) == t_INT)
3432 : {
3433 0 : long E = itos(y);
3434 0 : if (DEBUGLEVEL>1) err_printf("\n%ld done later at prec %ld\n",j,E);
3435 0 : set_avma(av);
3436 0 : vecsmalltrunc_append(retry, j);
3437 0 : if (E > prec) prec = E;
3438 : }
3439 : else
3440 : {
3441 69 : if (DEBUGLEVEL>1) err_printf("%ld ",j);
3442 69 : gel(ma,j) = gc_upto(av,y);
3443 : }
3444 : }
3445 20 : if (prec)
3446 : {
3447 0 : long k, l = lg(retry);
3448 0 : GEN y, nf = bnf_get_nf(bnf);
3449 0 : if (DEBUGLEVEL) pari_warn(warnprec,"makematal",prec);
3450 0 : nf = nfnewprec_shallow(nf,prec);
3451 0 : bnf = Buchall(nf, nf_FORCE, prec);
3452 0 : if (DEBUGLEVEL) err_printf("makematal, adding missing entries:");
3453 0 : for (k=1; k<l; k++)
3454 : {
3455 0 : pari_sp av = avma;
3456 0 : long j = retry[k];
3457 0 : y = get_y(bnf,W,B,NULL, pFB, j);
3458 0 : if (typ(y) == t_INT) pari_err_PREC("makematal");
3459 0 : if (DEBUGLEVEL>1) err_printf("%ld ",j);
3460 0 : gel(ma,j) = gc_upto(av,y);
3461 : }
3462 : }
3463 20 : if (DEBUGLEVEL>1) err_printf("\n");
3464 20 : return ma;
3465 : }
3466 :
3467 : enum { MATAL = 1, CYCGEN, UNITS };
3468 : GEN
3469 26726 : bnf_build_cycgen(GEN bnf)
3470 26726 : { return obj_checkbuild(bnf, CYCGEN, &makecycgen); }
3471 : GEN
3472 20 : bnf_build_matalpha(GEN bnf)
3473 20 : { return obj_checkbuild(bnf, MATAL, &makematal); }
3474 : GEN
3475 50738 : bnf_build_units(GEN bnf)
3476 50738 : { return obj_checkbuild(bnf, UNITS, &makeunits); }
3477 :
3478 : /* return fu in compact form if available; in terms of a fixed basis
3479 : * of S-units */
3480 : GEN
3481 70 : bnf_compactfu_mat(GEN bnf)
3482 : {
3483 70 : GEN X, U, SUnits = bnf_get_sunits(bnf);
3484 70 : if (!SUnits) return NULL;
3485 70 : X = gel(SUnits,1);
3486 70 : U = gel(SUnits,2); ZM_remove_unused(&U, &X);
3487 70 : return mkvec2(X, U);
3488 : }
3489 : /* return fu in compact form if available; individually as famat */
3490 : GEN
3491 37415 : bnf_compactfu(GEN bnf)
3492 : {
3493 37415 : GEN fu, X, U, SUnits = bnf_get_sunits(bnf);
3494 : long i, l;
3495 37415 : if (!SUnits) return NULL;
3496 37177 : X = gel(SUnits,1);
3497 37177 : U = gel(SUnits,2); l = lg(U); fu = cgetg(l, t_VEC);
3498 61201 : for (i = 1; i < l; i++)
3499 24024 : gel(fu,i) = famat_remove_trivial(mkmat2(X, gel(U,i)));
3500 37177 : return fu;
3501 : }
3502 : /* return expanded fu if available */
3503 : GEN
3504 285769 : bnf_has_fu(GEN bnf)
3505 : {
3506 285769 : GEN fu = obj_check(bnf, UNITS);
3507 285767 : if (fu) return vecsplice(fu, 1);
3508 264112 : fu = bnf_get_fu_nocheck(bnf);
3509 264111 : return (typ(fu) == t_MAT)? NULL: fu;
3510 : }
3511 : /* return expanded fu if available; build if cheap */
3512 : GEN
3513 285482 : bnf_build_cheapfu(GEN bnf)
3514 : {
3515 : GEN fu, SUnits;
3516 285482 : if ((fu = bnf_has_fu(bnf))) return fu;
3517 142 : if ((SUnits = bnf_get_sunits(bnf)))
3518 : {
3519 142 : pari_sp av = avma;
3520 142 : long e = gexpo(real_i(bnf_get_logfu(bnf)));
3521 142 : set_avma(av); if (e < 13) return vecsplice(bnf_build_units(bnf), 1);
3522 : }
3523 77 : return NULL;
3524 : }
3525 :
3526 : static GEN
3527 48887 : get_regulator(GEN A)
3528 : {
3529 48887 : pari_sp av = avma;
3530 : GEN R;
3531 :
3532 48887 : if (lg(A) == 1) return gen_1;
3533 48880 : R = det( rowslice(real_i(A), 1, lgcols(A)-2) );
3534 48880 : setabssign(R); return gc_leaf(av, R);
3535 : }
3536 :
3537 : /* return corrected archimedian components for elts of x (vector)
3538 : * (= log(sigma_i(x)) - log(|Nx|) / [K:Q]) */
3539 : static GEN
3540 40 : get_archclean(GEN nf, GEN x, long prec, int units)
3541 : {
3542 40 : long k, N, l = lg(x);
3543 40 : GEN M = cgetg(l, t_MAT);
3544 :
3545 40 : if (l == 1) return M;
3546 26 : N = nf_get_degree(nf);
3547 114 : for (k = 1; k < l; k++)
3548 : {
3549 88 : pari_sp av = avma;
3550 88 : GEN c = nf_cxlog(nf, gel(x,k), prec);
3551 88 : if (!c || (!units && !(c = cleanarch(c, N, NULL,prec)))) return NULL;
3552 88 : gel(M,k) = gc_GEN(av, c);
3553 : }
3554 26 : return M;
3555 : }
3556 : static void
3557 77 : SUnits_archclean(GEN nf, GEN SUnits, GEN *pmun, GEN *pC, long prec)
3558 : {
3559 77 : GEN ipi, M, X = gel(SUnits,1), U = gel(SUnits,2), G = gel(SUnits,3);
3560 77 : long k, N = nf_get_degree(nf), l = lg(X);
3561 :
3562 77 : M = cgetg(l, t_MAT);
3563 3640 : for (k = 1; k < l; k++)
3564 3563 : if (!(gel(M,k) = nf_cxlog(nf, gel(X,k), prec))) return;
3565 77 : ipi = invr(mppi(prec));
3566 77 : *pmun = cleanarch(RgM_ZM_mul(M, U), N, ipi, prec); /* not cleanarchunit ! */
3567 77 : if (*pmun) *pC = cleanarch(RgM_ZM_mul(M, G), N, ipi, prec);
3568 : }
3569 :
3570 : GEN
3571 97 : bnfnewprec_shallow(GEN bnf, long prec)
3572 : {
3573 97 : GEN nf0 = bnf_get_nf(bnf), nf, v, fu, matal, y, A, C;
3574 97 : GEN SUnits = bnf_get_sunits(bnf), Ur, Ga, Ge, M1, M2;
3575 97 : long r1, r2, prec0 = prec;
3576 :
3577 97 : nf_get_sign(nf0, &r1, &r2);
3578 97 : if (SUnits)
3579 : {
3580 77 : fu = matal = NULL;
3581 77 : prec += nbits2extraprec(gexpo(SUnits));
3582 : }
3583 : else
3584 : {
3585 20 : fu = bnf_build_units(bnf);
3586 20 : fu = vecslice(fu, 2, lg(fu)-1);
3587 20 : if (r1 + r2 > 1) {
3588 13 : long e = gexpo(bnf_get_logfu(bnf)) + 1 - TWOPOTBITS_IN_LONG;
3589 13 : if (e >= 0) prec += nbits2extraprec(e);
3590 : }
3591 20 : matal = bnf_build_matalpha(bnf);
3592 : }
3593 :
3594 97 : if (DEBUGLEVEL && prec0 != prec) pari_warn(warnprec,"bnfnewprec",prec);
3595 97 : for(C = NULL;;)
3596 0 : {
3597 97 : pari_sp av = avma;
3598 97 : nf = nfnewprec_shallow(nf0,prec);
3599 97 : if (SUnits)
3600 77 : SUnits_archclean(nf, SUnits, &A, &C, prec);
3601 : else
3602 : {
3603 20 : A = get_archclean(nf, fu, prec, 1);
3604 20 : if (A) C = get_archclean(nf, matal, prec, 0);
3605 : }
3606 97 : if (C) break;
3607 0 : set_avma(av); prec = precdbl(prec);
3608 0 : if (DEBUGLEVEL) pari_warn(warnprec,"bnfnewprec(extra)",prec);
3609 : }
3610 97 : y = leafcopy(bnf);
3611 97 : gel(y,3) = A;
3612 97 : gel(y,4) = C;
3613 97 : gel(y,7) = nf;
3614 97 : gel(y,8) = v = leafcopy(gel(bnf,8));
3615 97 : gel(v,2) = get_regulator(A);
3616 97 : v = gel(bnf,9);
3617 97 : if (lg(v) < 7) pari_err_TYPE("bnfnewprec [obsolete bnf format]", bnf);
3618 97 : Ur = gel(v,1);
3619 97 : Ge = gel(v,4);
3620 97 : Ga = nfV_cxlog(nf, Ge, prec);
3621 97 : M1 = gel(v,5);
3622 97 : M2 = gel(v,6);
3623 97 : gel(y,9) = get_clg2(bnf_get_cyc(bnf), Ga, C, Ur, Ge, M1, M2);
3624 97 : return y;
3625 : }
3626 : GEN
3627 7 : bnfnewprec(GEN bnf, long prec)
3628 : {
3629 7 : pari_sp av = avma;
3630 7 : return gc_GEN(av, bnfnewprec_shallow(checkbnf(bnf), prec));
3631 : }
3632 :
3633 : GEN
3634 0 : bnrnewprec_shallow(GEN bnr, long prec)
3635 : {
3636 0 : GEN y = cgetg(7,t_VEC);
3637 : long i;
3638 0 : gel(y,1) = bnfnewprec_shallow(bnr_get_bnf(bnr), prec);
3639 0 : for (i=2; i<7; i++) gel(y,i) = gel(bnr,i);
3640 0 : return y;
3641 : }
3642 : GEN
3643 7 : bnrnewprec(GEN bnr, long prec)
3644 : {
3645 7 : GEN y = cgetg(7,t_VEC);
3646 : long i;
3647 7 : checkbnr(bnr);
3648 7 : gel(y,1) = bnfnewprec(bnr_get_bnf(bnr), prec);
3649 42 : for (i=2; i<7; i++) gel(y,i) = gcopy(gel(bnr,i));
3650 7 : return y;
3651 : }
3652 :
3653 : static GEN
3654 65205 : buchall_end(GEN nf,GEN res, GEN clg2, GEN W, GEN B, GEN A, GEN C,GEN Vbase)
3655 : {
3656 65205 : GEN z = obj_init(9, 3);
3657 65205 : gel(z,1) = W;
3658 65205 : gel(z,2) = B;
3659 65205 : gel(z,3) = A;
3660 65205 : gel(z,4) = C;
3661 65205 : gel(z,5) = Vbase;
3662 65205 : gel(z,6) = gen_0;
3663 65205 : gel(z,7) = nf;
3664 65205 : gel(z,8) = res;
3665 65205 : gel(z,9) = clg2;
3666 65205 : return z;
3667 : }
3668 :
3669 : GEN
3670 2632 : bnfinit0(GEN P, long flag, GEN data, long prec)
3671 : {
3672 2632 : double c1 = 0., c2 = 0.;
3673 2632 : long fl, relpid = BNF_RELPID, max_fact = MAXTRY_FACT, idex = 0, nbthr = 1, s;
3674 :
3675 2632 : if (data)
3676 : {
3677 21 : long lx = lg(data);
3678 21 : if (typ(data) != t_VEC || lx > 7) pari_err_TYPE("bnfinit",data);
3679 21 : switch(lx)
3680 : {
3681 0 : case 7: nbthr = itou(gel(data,6)); if (nbthr <= 1) nbthr = 1-nbthr;
3682 0 : case 6: idex = itou(gel(data,5));
3683 0 : case 5: s = itou(gel(data,4)); if (s) max_fact = s;
3684 0 : case 4: s = itos(gel(data,3)); if (s) relpid = s < 0 ? 0 : s;
3685 14 : case 3: c2 = gtodouble(gel(data,2));
3686 21 : case 2: c1 = gtodouble(gel(data,1));
3687 : }
3688 : }
3689 2632 : switch(flag)
3690 : {
3691 1778 : case 2:
3692 1778 : case 0: fl = 0; break;
3693 854 : case 1: fl = nf_FORCE; break;
3694 0 : default: pari_err_FLAG("bnfinit");
3695 : return NULL; /* LCOV_EXCL_LINE */
3696 : }
3697 2632 : return Buchall_param(P, c1, c2, relpid, max_fact, idex, nbthr, fl, prec);
3698 : }
3699 : GEN
3700 62578 : Buchall(GEN P, long flag, long prec)
3701 62578 : { return Buchall_param(P, 0., 0., BNF_RELPID, MAXTRY_FACT, 0, 1, flag & nf_FORCE, prec); }
3702 :
3703 : static GEN
3704 1204 : Buchall_deg1(GEN nf)
3705 : {
3706 1204 : GEN v = cgetg(1,t_VEC), m = cgetg(1,t_MAT);
3707 1204 : GEN res, W, A, B, C, Vbase = cgetg(1,t_COL);
3708 1204 : GEN fu = v, R = gen_1, zu = mkvec2(gen_2, gen_m1);
3709 1204 : GEN clg1 = mkvec3(gen_1,v,v), clg2 = mkvecn(6, m,m,m,v,m,m);
3710 :
3711 1204 : W = A = B = C = m; res = mkvec5(clg1, R, gen_1, zu, fu);
3712 1204 : return buchall_end(nf,res,clg2,W,B,A,C,Vbase);
3713 : }
3714 :
3715 : /* return (small set of) indices of columns generating the same lattice as x.
3716 : * Assume HNF(x) is inexpensive (few rows, many columns).
3717 : * Dichotomy approach since interesting columns may be at the very end */
3718 : GEN
3719 64007 : extract_full_lattice(GEN x)
3720 : {
3721 64007 : long dj, j, k, l = lg(x);
3722 : GEN h, h2, H, v;
3723 :
3724 64007 : if (l < 200) return NULL; /* not worth it */
3725 :
3726 1 : v = vecsmalltrunc_init(l);
3727 1 : H = ZM_hnf(x);
3728 1 : h = cgetg(1, t_MAT);
3729 1 : dj = 1;
3730 43 : for (j = 1; j < l; )
3731 : {
3732 43 : pari_sp av = avma;
3733 43 : long lv = lg(v);
3734 :
3735 145 : for (k = 0; k < dj; k++) v[lv+k] = j+k;
3736 43 : setlg(v, lv + dj);
3737 43 : h2 = ZM_hnf(vecpermute(x, v));
3738 43 : if (ZM_equal(h, h2))
3739 : { /* these dj columns can be eliminated */
3740 17 : set_avma(av); setlg(v, lv);
3741 17 : j += dj;
3742 17 : if (j >= l) break;
3743 17 : dj <<= 1;
3744 17 : if (j + dj >= l) { dj = (l - j) >> 1; if (!dj) dj = 1; }
3745 : }
3746 26 : else if (dj > 1)
3747 : { /* at least one interesting column, try with first half of this set */
3748 17 : set_avma(av); setlg(v, lv);
3749 17 : dj >>= 1; /* > 0 */
3750 : }
3751 : else
3752 : { /* this column should be kept */
3753 9 : if (ZM_equal(h2, H)) break;
3754 8 : h = h2; j++;
3755 : }
3756 : }
3757 1 : return v;
3758 : }
3759 :
3760 : static void
3761 64073 : init_rel(RELCACHE_t *cache, FB_t *F, long add_need)
3762 : {
3763 64073 : const long n = F->KC + add_need; /* expected # of needed relations */
3764 : long i, j, k, p;
3765 : GEN c, P;
3766 : GEN R;
3767 :
3768 64073 : if (DEBUGLEVEL) err_printf("KCZ = %ld, KC = %ld, n = %ld\n", F->KCZ,F->KC,n);
3769 64073 : reallocate(cache, 10*n + 50); /* make room for lots of relations */
3770 64074 : cache->chk = cache->base;
3771 64074 : cache->end = cache->base + n;
3772 64074 : cache->relsup = add_need;
3773 64074 : cache->last = cache->base;
3774 64074 : cache->missing = lg(cache->basis) - 1;
3775 306267 : for (i = 1; i <= F->KCZ; i++)
3776 : { /* trivial relations (p) = prod P^e */
3777 242193 : p = F->FB[i]; P = gel(F->LV,p);
3778 242193 : if (!isclone(P)) continue;
3779 :
3780 : /* all prime divisors in FB */
3781 169012 : c = zero_Flv(F->KC); k = F->iLP[p];
3782 169010 : R = c; c += k;
3783 539697 : for (j = lg(P)-1; j; j--) c[j] = pr_get_e(gel(P,j));
3784 169009 : add_rel(cache, F, R, k+1, pr_get_p(gel(P,1)), 0);
3785 : }
3786 64074 : }
3787 :
3788 : /* Let z = \zeta_n in nf. List of not-obviously-dependent generators for
3789 : * cyclotomic units modulo torsion in Q(z) [independent when n a prime power]:
3790 : * - z^a - 1, n/(a,n) not a prime power, a \nmid n unless a=1, 1 <= a < n/2
3791 : * - (Z^a - 1)/(Z - 1), p^k || n, Z = z^{n/p^k}, (p,a) = 1, 1 < a <= (p^k-1)/2
3792 : */
3793 : GEN
3794 64074 : nfcyclotomicunits(GEN nf, GEN zu)
3795 : {
3796 64074 : long n = itos(gel(zu, 1)), n2, lP, i, a;
3797 : GEN z, fa, P, E, L, mz, powz;
3798 64074 : if (n <= 6) return cgetg(1, t_VEC);
3799 :
3800 1911 : z = algtobasis(nf,gel(zu, 2));
3801 1911 : if ((n & 3) == 2) { n = n >> 1; z = ZC_neg(z); } /* ensure n != 2 (mod 4) */
3802 1911 : n2 = n/2;
3803 1911 : mz = zk_multable(nf, z); /* multiplication by z */
3804 1911 : powz = cgetg(n2, t_VEC); gel(powz,1) = z;
3805 6286 : for (i = 2; i < n2; i++) gel(powz,i) = ZM_ZC_mul(mz, gel(powz,i-1));
3806 : /* powz[i] = z^i */
3807 :
3808 1911 : L = vectrunc_init(n);
3809 1911 : fa = factoru(n);
3810 1911 : P = gel(fa,1); lP = lg(P);
3811 1911 : E = gel(fa,2);
3812 4613 : for (i = 1; i < lP; i++)
3813 : { /* second kind */
3814 2702 : long p = P[i], k = E[i], pk = upowuu(p,k), pk2 = (pk-1) / 2;
3815 2702 : GEN u = gen_1;
3816 4970 : for (a = 2; a <= pk2; a++)
3817 : {
3818 2268 : u = nfadd(nf, u, gel(powz, (n/pk) * (a-1))); /* = (Z^a-1)/(Z-1) */
3819 2268 : if (a % p) vectrunc_append(L, u);
3820 : }
3821 : }
3822 6160 : if (lP > 2) for (a = 1; a < n2; a++)
3823 : { /* first kind, when n not a prime power */
3824 : ulong p;
3825 4249 : if (a > 1 && (n % a == 0 || uisprimepower(n/ugcd(a,n), &p))) continue;
3826 1869 : vectrunc_append(L, nfadd(nf, gel(powz, a), gen_m1));
3827 : }
3828 1911 : return L;
3829 : }
3830 : static void
3831 64074 : add_cyclotomic_units(GEN nf, GEN zu, RELCACHE_t *cache, FB_t *F)
3832 : {
3833 64074 : pari_sp av = avma;
3834 64074 : GEN L = nfcyclotomicunits(nf, zu);
3835 64073 : long i, l = lg(L);
3836 64073 : if (l > 1)
3837 : {
3838 1911 : GEN R = zero_Flv(F->KC);
3839 5950 : for(i = 1; i < l; i++) add_rel(cache, F, R, F->KC+1, gel(L,i), 0);
3840 : }
3841 64073 : set_avma(av);
3842 64073 : }
3843 :
3844 : static GEN
3845 122448 : trim_list(FB_t *F)
3846 : {
3847 122448 : pari_sp av = avma;
3848 122448 : GEN v, L_jid = F->L_jid, minidx = F->minidx, present = zero_Flv(F->KC);
3849 122448 : long i, j, imax = minss(lg(L_jid), F->KC + 1);
3850 :
3851 122448 : v = cgetg(imax, t_VECSMALL);
3852 1332157 : for (i = j = 1; i < imax; i++)
3853 : {
3854 1209709 : long k = minidx[ L_jid[i] ];
3855 1209709 : if (!present[k]) { v[j++] = L_jid[i]; present[k] = 1; }
3856 : }
3857 122448 : setlg(v, j); return gc_leaf(av, v);
3858 : }
3859 :
3860 : /* x t_INT or primitive ZC */
3861 : static void
3862 1649 : try_elt(RELCACHE_t *cache, FB_t *F, GEN nf, GEN x, FACT *fact)
3863 : {
3864 1649 : pari_sp av = avma;
3865 : long nz;
3866 : GEN R;
3867 :
3868 1649 : if (typ(x) == t_INT /* 2nd path can't fail */
3869 1649 : || !can_factor(F, nf, NULL, x, nfnorm(nf, x), fact)) return;
3870 : /* smooth element */
3871 1428 : R = set_fact(F, fact, NULL, &nz);
3872 : /* make sure we get maximal rank first, then allow all relations */
3873 1428 : (void)add_rel(cache, F, R, nz, x, 0);
3874 1428 : set_avma(av);
3875 : }
3876 :
3877 : static void
3878 55507 : matenlarge(GEN C, long h)
3879 : {
3880 55507 : GEN _0 = zerocol(h);
3881 : long i;
3882 1264920 : for (i = lg(C); --i; ) gel(C,i) = shallowconcat(gel(C,i), _0);
3883 55508 : }
3884 :
3885 : /* E = floating point embeddings */
3886 : static GEN
3887 55507 : matbotidembs(RELCACHE_t *cache, GEN E)
3888 : {
3889 55507 : long w = cache->last - cache->chk, h = cache->last - cache->base;
3890 55507 : long j, d = h - w, hE = nbrows(E);
3891 55507 : GEN y = cgetg(w+1,t_MAT), _0 = zerocol(h);
3892 217579 : for (j = 1; j <= w; j++)
3893 : {
3894 162072 : GEN c = shallowconcat(gel(E,j), _0);
3895 162072 : if (d + j >= 1) gel(c, d + j + hE) = gen_1;
3896 162072 : gel(y,j) = c;
3897 : }
3898 55507 : return y;
3899 : }
3900 : static GEN
3901 62485 : matbotid(RELCACHE_t *cache)
3902 : {
3903 62485 : long w = cache->last - cache->chk, h = cache->last - cache->base;
3904 62485 : long j, d = h - w;
3905 62485 : GEN y = cgetg(w+1,t_MAT);
3906 851548 : for (j = 1; j <= w; j++)
3907 : {
3908 789064 : GEN c = zerocol(h);
3909 789063 : if (d + j >= 1) gel(c, d + j) = gen_1;
3910 789063 : gel(y,j) = c;
3911 : }
3912 62484 : return y;
3913 : }
3914 :
3915 : static long
3916 73 : myprecdbl(long prec, GEN C)
3917 : {
3918 73 : long p = prec < 1280? precdbl(prec): (long)(prec * 1.5);
3919 73 : if (C) p = maxss(p, minss(3*p, prec + nbits2extraprec(gexpo(C))));
3920 73 : return p;
3921 : }
3922 :
3923 : static GEN
3924 57721 : _nfnewprec(GEN nf, long prec, long *isclone)
3925 : {
3926 57721 : GEN NF = gclone(nfnewprec_shallow(nf, prec));
3927 57721 : if (*isclone) gunclone(nf);
3928 57721 : *isclone = 1; return NF;
3929 : }
3930 :
3931 : /* In small_norm, LLL reduction produces v0 in I such that
3932 : * T2(v0) <= (4/3)^((n-1)/2) (NI sqrt(disc(K)))^(2/n)
3933 : * NI <= LIMCMAX^2. We consider v with T2(v) ~ T2(v0), hence
3934 : * Nv <= ((4/3)^((n-1)/2) / n)^(n/2) LIMCMAX^2 sqrt(disc(K)) */
3935 : static long
3936 63993 : small_norm_prec(long N, double LOGD, long LIMCMAX)
3937 : {
3938 63993 : double a = N/2. * ((N-1)/2.*log(4./3) - log((double)N));
3939 63993 : double b = 2*log((double)LIMCMAX) + LOGD/2;
3940 63993 : return nbits2prec(BITS_IN_LONG + (a + b) / M_LN2);
3941 : }
3942 :
3943 : /* Nrelid = nb relations per ideal, possibly 0. If flag is set, keep data in
3944 : * algebraic form. */
3945 : GEN
3946 65213 : Buchall_param(GEN P, double cbach, double cbach2, long Nrelid, long max_fact, long idex, long nbthr, long flag, long prec)
3947 : {
3948 : pari_timer T;
3949 65213 : pari_sp av0 = avma, av, av2;
3950 : long PREC, N, R1, R2, RU, low, high, LIMC0, LIMC, LIMC2, LIMCMAX, zc, i;
3951 65213 : long LIMres, bit = 0, flag_nfinit = 0, nfisclone = 0;
3952 65213 : long nreldep, sfb_trials, need, old_need, precdouble = 0, TRIES = 0;
3953 : long done_small, small_fail, fail_limit, squash_index;
3954 : double LOGD, LOGD2, lim;
3955 65213 : GEN computed = NULL, fu = NULL, zu, nf, D, A, W, R, h, Ce, PERM;
3956 : GEN small_multiplier, auts, cyclic, embs, SUnits;
3957 : GEN res, L, invhr, B, C, lambda, dep, clg1, clg2, Vbase;
3958 65213 : const char *precpb = NULL;
3959 65213 : REL_t *old_cache = NULL;
3960 : nfmaxord_t nfT;
3961 : RELCACHE_t cache;
3962 : FB_t F;
3963 : GRHcheck_t GRHcheck;
3964 : FACT *fact;
3965 :
3966 65213 : if (DEBUGLEVEL) timer_start(&T);
3967 65213 : P = get_nfpol(P, &nf);
3968 65193 : if (degpol(P)==2) Nrelid = 0;
3969 65193 : if (nf)
3970 3871 : D = nf_get_disc(nf);
3971 : else
3972 : {
3973 61322 : nfinit_basic(&nfT, P);
3974 61327 : D = nfT.dK;
3975 61327 : if (!ZX_is_monic(nfT.T0))
3976 : {
3977 14 : pari_warn(warner,"nonmonic polynomial in bnfinit, using polredbest");
3978 14 : flag_nfinit = nf_RED;
3979 : }
3980 : }
3981 65197 : PREC = maxss(DEFAULTPREC, prec);
3982 65199 : N = degpol(P);
3983 65198 : if (N <= 1)
3984 : {
3985 1204 : if (!nf) nf = nfinit_complete(&nfT, flag_nfinit, PREC);
3986 1204 : return gc_GEN(av0, Buchall_deg1(nf));
3987 : }
3988 63994 : D = absi_shallow(D);
3989 63994 : LOGD = dbllog2(D) * M_LN2;
3990 63994 : LOGD2 = LOGD*LOGD;
3991 63994 : LIMCMAX = (long)(4.*LOGD2);
3992 63994 : if (nf) PREC = maxss(PREC, nf_get_prec(nf));
3993 63994 : PREC = maxss(PREC, nbits2prec((long)(LOGD2 * 0.02) + N*N));
3994 :
3995 63994 : if (nf) PREC = maxss(PREC, nf_get_prec(nf));
3996 63994 : PREC = maxss(PREC, nbits2prec((long)(LOGD2 * 0.02) + N*N));
3997 63993 : PREC = maxss(PREC, small_norm_prec(N, LOGD, LIMCMAX));
3998 63993 : if (DEBUGLEVEL) err_printf("PREC = %ld\n", PREC);
3999 :
4000 63993 : if (!nf)
4001 60332 : nf = nfinit_complete(&nfT, flag_nfinit, PREC);
4002 3661 : else if (nf_get_prec(nf) < PREC)
4003 161 : nf = nfnewprec_shallow(nf, PREC);
4004 63999 : zu = nfrootsof1(nf);
4005 63996 : gel(zu,2) = nf_to_scalar_or_alg(nf, gel(zu,2));
4006 :
4007 63996 : nf_get_sign(nf, &R1, &R2); RU = R1+R2;
4008 63996 : auts = automorphism_matrices(nf, &cyclic);
4009 64000 : F.embperm = automorphism_perms(nf_get_M(nf), auts, cyclic, R1, R2, N);
4010 64000 : if (DEBUGLEVEL)
4011 : {
4012 0 : timer_printf(&T, "nfinit & nfrootsof1");
4013 0 : err_printf("%s bnf: R1 = %ld, R2 = %ld\nD = %Ps\n",
4014 : flag? "Algebraic": "Floating point", R1,R2, D);
4015 : }
4016 64000 : if (LOGD < 20.)
4017 : { /* tiny disc, Minkowski may be smaller than Bach */
4018 62537 : lim = exp(-N + R2 * log(4/M_PI) + LOGD/2) * sqrt(2*M_PI*N);
4019 62537 : if (lim < 3) lim = 3;
4020 : }
4021 : else /* to be ignored */
4022 1463 : lim = -1;
4023 64000 : if (cbach > 12.) {
4024 0 : if (cbach2 < cbach) cbach2 = cbach;
4025 0 : cbach = 12.;
4026 : }
4027 64000 : if (cbach < 0.)
4028 0 : pari_err_DOMAIN("Buchall","Bach constant","<",gen_0,dbltor(cbach));
4029 :
4030 64000 : cache.base = NULL; F.subFB = NULL; F.LP = NULL; SUnits = Ce = NULL;
4031 64000 : init_GRHcheck(&GRHcheck, N, R1, LOGD);
4032 64001 : high = low = LIMC0 = maxss((long)(cbach2*LOGD2), 1);
4033 312174 : while (!GRHchk(nf, &GRHcheck, high)) { low = high; high *= 2; }
4034 248216 : while (high - low > 1)
4035 : {
4036 184217 : long test = (low+high)/2;
4037 184217 : if (GRHchk(nf, &GRHcheck, test)) high = test; else low = test;
4038 : }
4039 63999 : LIMC2 = (high == LIMC0+1 && GRHchk(nf, &GRHcheck, LIMC0))? LIMC0: high;
4040 63999 : if (LIMC2 > LIMCMAX) LIMC2 = LIMCMAX;
4041 : /* Assuming GRH, {P, NP <= LIMC2} generate Cl(K) */
4042 63999 : if (DEBUGLEVEL) err_printf("LIMC2 = %ld\n", LIMC2);
4043 63999 : LIMC0 = (long)(cbach*LOGD2); /* initial value for LIMC */
4044 63999 : LIMC = cbach? LIMC0: LIMC2; /* use {P, NP <= LIMC} as a factorbase */
4045 63999 : LIMC = maxss(LIMC, nthideal(&GRHcheck, nf, N));
4046 64000 : if (DEBUGLEVEL) timer_printf(&T, "computing Bach constant");
4047 64000 : LIMres = primeneeded(N, R1, R2, LOGD);
4048 64001 : cache_prime_dec(&GRHcheck, LIMres, nf);
4049 : /* invhr ~ 2^r1 (2pi)^r2 / sqrt(D) w * Res(zeta_K, s=1) = 1 / hR */
4050 128000 : invhr = gmul(gdiv(gmul2n(powru(mppi(DEFAULTPREC), R2), RU),
4051 64000 : mulri(gsqrt(D,DEFAULTPREC),gel(zu,1))),
4052 : compute_invres(&GRHcheck, LIMres));
4053 64000 : if (DEBUGLEVEL) timer_printf(&T, "computing inverse of hR");
4054 64000 : av = avma;
4055 :
4056 66250 : START:
4057 66250 : if (DEBUGLEVEL) timer_start(&T);
4058 66250 : if (TRIES) LIMC = bnf_increase_LIMC(LIMC,LIMCMAX);
4059 66250 : if (DEBUGLEVEL && LIMC > LIMC0)
4060 0 : err_printf("%s*** Bach constant: %f\n", TRIES?"\n":"", LIMC/LOGD2);
4061 66250 : if (cache.base)
4062 : {
4063 : REL_t *rel;
4064 3693 : for (i = 1, rel = cache.base + 1; rel < cache.last; rel++)
4065 3620 : if (rel->m) i++;
4066 73 : computed = cgetg(i, t_VEC);
4067 3693 : for (i = 1, rel = cache.base + 1; rel < cache.last; rel++)
4068 3620 : if (rel->m) gel(computed, i++) = rel->m;
4069 73 : computed = gclone(computed); delete_cache(&cache);
4070 : }
4071 66250 : TRIES++; set_avma(av);
4072 66250 : if (F.LP) delete_FB(&F);
4073 66250 : if (LIMC2 < LIMC) LIMC2 = LIMC;
4074 66250 : if (DEBUGLEVEL) { err_printf("LIMC = %ld, LIMC2 = %ld\n",LIMC,LIMC2); }
4075 :
4076 66250 : FBgen(&F, nf, N, LIMC, LIMC2, &GRHcheck);
4077 66247 : if (!F.KC) goto START;
4078 66247 : av = avma;
4079 66247 : subFBgen(&F,auts,cyclic,lim < 0? LIMC2: mindd(lim,LIMC2),MINSFB);
4080 66251 : if (lg(F.subFB) == 1) goto START;
4081 64074 : if (DEBUGLEVEL)
4082 0 : timer_printf(&T, "factorbase (#subFB = %ld) and ideal permutations",
4083 0 : lg(F.subFB)-1);
4084 :
4085 64074 : fact = (FACT*)stack_malloc((F.KC+1)*sizeof(FACT));
4086 64074 : PERM = leafcopy(F.perm); /* to be restored in case of precision increase */
4087 64072 : cache.basis = zero_Flm_copy(F.KC,F.KC);
4088 64073 : small_multiplier = zero_Flv(F.KC);
4089 64073 : done_small = small_fail = squash_index = zc = sfb_trials = nreldep = 0;
4090 64073 : fail_limit = F.KC + 1;
4091 64073 : W = A = R = NULL;
4092 64073 : av2 = avma;
4093 64073 : init_rel(&cache, &F, RELSUP + RU-1);
4094 64074 : old_need = need = cache.end - cache.last;
4095 64074 : add_cyclotomic_units(nf, zu, &cache, &F);
4096 64073 : if (DEBUGLEVEL) err_printf("\n");
4097 64073 : cache.end = cache.last + need;
4098 :
4099 64073 : if (computed)
4100 : {
4101 1722 : for (i = 1; i < lg(computed); i++)
4102 1649 : try_elt(&cache, &F, nf, gel(computed, i), fact);
4103 73 : gunclone(computed);
4104 73 : if (DEBUGLEVEL && i > 1)
4105 0 : timer_printf(&T, "including already computed relations");
4106 73 : need = 0;
4107 : }
4108 :
4109 : do
4110 : {
4111 : GEN Ar, C0;
4112 : do
4113 : {
4114 122601 : pari_sp av4 = avma;
4115 122601 : if (need > 0)
4116 : {
4117 122448 : long oneed = cache.end - cache.last;
4118 : /* Test below can be true if small_norm did not find enough linearly
4119 : * dependent relations */
4120 122448 : if (need < oneed) need = oneed;
4121 122448 : pre_allocate(&cache, need+lg(auts)-1+(R ? lg(W)-1 : 0));
4122 122448 : cache.end = cache.last + need;
4123 122448 : F.L_jid = trim_list(&F);
4124 : }
4125 122600 : if (need > 0 && Nrelid > 0 && (done_small <= F.KC+1 || A) &&
4126 70353 : small_fail <= fail_limit &&
4127 70353 : cache.last < cache.base + 2*F.KC+2*RU+RELSUP /* heuristic */)
4128 : {
4129 66762 : long j, k, LIE = (R && lg(W) > 1 && (done_small % 2));
4130 66762 : REL_t *last = cache.last;
4131 66762 : pari_sp av3 = avma;
4132 66762 : if (LIE)
4133 : { /* We have full rank for class group and unit. The following tries to
4134 : * improve the prime group lattice by looking for relations involving
4135 : * the primes generating the class group. */
4136 3393 : long n = lg(W)-1; /* need n relations to squash the class group */
4137 3393 : F.L_jid = vecslice(F.perm, 1, n);
4138 3393 : cache.end = cache.last + n;
4139 : /* Lie to the add_rel subsystem: pretend we miss relations involving
4140 : * the primes generating the class group (and only those). */
4141 3393 : cache.missing = n;
4142 10620 : for ( ; n > 0; n--) mael(cache.basis, F.perm[n], F.perm[n]) = 0;
4143 : }
4144 66762 : j = done_small % (F.KC+1);
4145 66762 : if (j && !A)
4146 : { /* Prevent considering both P_iP_j and P_jP_i in small_norm */
4147 : /* Not all elements end up in F.L_jid (eliminated by hnfspec/add or
4148 : * by trim_list): keep track of which ideals are being considered
4149 : * at each run. */
4150 414 : long mj = small_multiplier[j];
4151 6569 : for (i = k = 1; i < lg(F.L_jid); i++)
4152 6155 : if (F.L_jid[i] > mj)
4153 : {
4154 6155 : small_multiplier[F.L_jid[i]] = j;
4155 6155 : F.L_jid[k++] = F.L_jid[i];
4156 : }
4157 414 : setlg(F.L_jid, k);
4158 : }
4159 66762 : if (lg(F.L_jid) > 1) small_norm(&cache, &F, nf, Nrelid, max_fact, idex, nbthr, fact, j);
4160 66764 : F.L_jid = F.perm; set_avma(av3);
4161 66764 : if (!A && cache.last != last) small_fail = 0; else small_fail++;
4162 66764 : if (LIE)
4163 : { /* restore add_rel subsystem: undo above lie */
4164 3393 : long n = lg(W) - 1;
4165 10620 : for ( ; n > 0; n--) mael(cache.basis, F.perm[n], F.perm[n]) = 1;
4166 3393 : cache.missing = 0;
4167 : }
4168 66764 : cache.end = cache.last;
4169 66764 : done_small++;
4170 66764 : need = F.sfb_chg = 0;
4171 : }
4172 122602 : if (need > 0)
4173 : { /* Random relations */
4174 55685 : if (++nreldep > F.MAXDEPSIZESFB) {
4175 14 : if (++sfb_trials > SFB_MAX && LIMC < LIMCMAX/2) goto START;
4176 14 : F.sfb_chg = sfb_INCREASE;
4177 14 : nreldep = 0;
4178 : }
4179 55671 : else if (!(nreldep % F.MAXDEPSFB))
4180 26405 : F.sfb_chg = sfb_CHANGE;
4181 55685 : if (F.sfb_chg && !subFB_change(&F)) goto START;
4182 55612 : rnd_rel(&cache, &F, nf, max_fact, nbthr, fact);
4183 55613 : F.L_jid = F.perm;
4184 : }
4185 122530 : if (DEBUGLEVEL) timer_start(&T);
4186 122530 : if (precpb)
4187 : {
4188 : REL_t *rel;
4189 80 : if (DEBUGLEVEL)
4190 : {
4191 0 : char str[64]; sprintf(str,"Buchall_param (%s)",precpb);
4192 0 : pari_warn(warnprec,str,PREC);
4193 : }
4194 80 : nf = _nfnewprec(nf, PREC, &nfisclone);
4195 80 : precdouble++; precpb = NULL;
4196 :
4197 80 : if (flag)
4198 : { /* recompute embs only, no need to redo HNF */
4199 38 : long j, le = lg(embs), lC = lg(C);
4200 38 : GEN E, M = nf_get_M(nf);
4201 38 : set_avma(av4);
4202 12611 : for (rel = cache.base+1, i = 1; i < le; i++,rel++)
4203 12573 : gel(embs,i) = rel_embed(rel, &F, embs, i, M, RU, R1, PREC);
4204 38 : E = RgM_ZM_mul(embs, rowslice(C, RU+1, nbrows(C)));
4205 12611 : for (j = 1; j < lC; j++)
4206 65595 : for (i = 1; i <= RU; i++) gcoeff(C,i,j) = gcoeff(E,i,j);
4207 38 : av4 = avma;
4208 : }
4209 : else
4210 : { /* recompute embs + HNF */
4211 10318 : for(i = 1; i < lg(PERM); i++) F.perm[i] = PERM[i];
4212 42 : cache.chk = cache.base;
4213 42 : W = NULL;
4214 : }
4215 80 : if (DEBUGLEVEL) timer_printf(&T, "increasing accuracy");
4216 : }
4217 122530 : set_avma(av4);
4218 122530 : if (cache.chk != cache.last)
4219 : { /* Reduce relation matrices */
4220 122410 : long l = cache.last - cache.chk + 1, j;
4221 122410 : GEN mat = cgetg(l, t_MAT);
4222 : REL_t *rel;
4223 :
4224 1126924 : for (j=1,rel = cache.chk + 1; j < l; rel++,j++) gel(mat,j) = rel->R;
4225 122408 : if (!flag || W)
4226 : {
4227 59924 : embs = get_embs(&F, &cache, nf, embs, PREC);
4228 59924 : if (DEBUGLEVEL && timer_get(&T) > 1)
4229 0 : timer_printf(&T, "floating point embeddings");
4230 : }
4231 122409 : if (!W)
4232 : { /* never reduced before */
4233 64116 : C = flag? matbotid(&cache): embs;
4234 64115 : W = hnfspec_i(mat, F.perm, &dep, &B, &C, F.subFB ? lg(F.subFB)-1:0);
4235 64116 : if (DEBUGLEVEL)
4236 0 : timer_printf(&T, "hnfspec [%ld x %ld]", lg(F.perm)-1, l-1);
4237 64116 : if (flag)
4238 : {
4239 62485 : PREC += nbits2extraprec(gexpo(C));
4240 62485 : if (nf_get_prec(nf) < PREC) nf = _nfnewprec(nf, PREC, &nfisclone);
4241 62485 : embs = get_embs(&F, &cache, nf, embs, PREC);
4242 62485 : C = vconcat(RgM_ZM_mul(embs, C), C);
4243 : }
4244 64116 : if (DEBUGLEVEL)
4245 0 : timer_printf(&T, "hnfspec floating points");
4246 : }
4247 : else
4248 : {
4249 58293 : long k = lg(embs);
4250 58293 : GEN E = vecslice(embs, k-l+1,k-1);
4251 58294 : if (flag)
4252 : {
4253 55507 : E = matbotidembs(&cache, E);
4254 55507 : matenlarge(C, cache.last - cache.chk);
4255 : }
4256 58294 : W = hnfadd_i(W, F.perm, &dep, &B, &C, mat, E);
4257 58293 : if (DEBUGLEVEL)
4258 0 : timer_printf(&T, "hnfadd (%ld + %ld)", l-1, lg(dep)-1);
4259 : }
4260 122409 : (void)gc_all(av2, 5, &W,&C,&B,&dep,&embs);
4261 122410 : cache.chk = cache.last;
4262 : }
4263 120 : else if (!W)
4264 : {
4265 0 : need = old_need;
4266 0 : F.L_jid = vecslice(F.perm, 1, need);
4267 0 : continue;
4268 : }
4269 122530 : need = F.KC - (lg(W)-1) - (lg(B)-1);
4270 122530 : if (!need && cache.missing)
4271 : { /* The test above will never be true except if 27449|class number.
4272 : * Ensure that if we have maximal rank for the ideal lattice, then
4273 : * cache.missing == 0. */
4274 14 : for (i = 1; cache.missing; i++)
4275 7 : if (!mael(cache.basis, i, i))
4276 : {
4277 : long j;
4278 7 : cache.missing--; mael(cache.basis, i, i) = 1;
4279 427 : for (j = i+1; j <= F.KC; j++) mael(cache.basis, j, i) = 0;
4280 : }
4281 : }
4282 122530 : zc = (lg(C)-1) - (lg(B)-1) - (lg(W)-1);
4283 122530 : if (RU-1-zc > 0) need = minss(need + RU-1-zc, F.KC); /* for units */
4284 122530 : if (need)
4285 : { /* dependent rows */
4286 811 : F.L_jid = vecslice(F.perm, 1, need);
4287 811 : vecsmall_sort(F.L_jid);
4288 811 : if (need != old_need) { nreldep = 0; old_need = need; }
4289 : }
4290 : else
4291 : { /* If the relation lattice is too small, check will be > 1 and we will
4292 : * do a new run of small_norm/rnd_rel asking for 1 relation. This often
4293 : * gives a relation involving L_jid[1]. We rotate the first element of
4294 : * L_jid in order to increase the probability of finding relations that
4295 : * increases the lattice. */
4296 121719 : long j, n = lg(W) - 1;
4297 121719 : if (n > 1 && squash_index % n)
4298 : {
4299 9025 : F.L_jid = leafcopy(F.perm);
4300 36564 : for (j = 1; j <= n; j++)
4301 27539 : F.L_jid[j] = F.perm[1 + (j + squash_index - 1) % n];
4302 : }
4303 : else
4304 112694 : F.L_jid = F.perm;
4305 121719 : squash_index++;
4306 : }
4307 : }
4308 122530 : while (need);
4309 :
4310 121719 : if (!A)
4311 : {
4312 64081 : small_fail = old_need = 0;
4313 64081 : fail_limit = maxss(F.KC / FAIL_DIVISOR, MINFAIL);
4314 : }
4315 121719 : A = vecslice(C, 1, zc); /* cols corresponding to units */
4316 121718 : if (flag) A = rowslice(A, 1, RU);
4317 121718 : Ar = real_i(A);
4318 121719 : R = compute_multiple_of_R(Ar, RU, N, &need, &bit, &lambda);
4319 121719 : if (need < old_need) small_fail = 0;
4320 : #if 0 /* A good idea if we are indeed stuck but needs tuning */
4321 : /* we have computed way more relations than should be necessary */
4322 : if (TRIES < 3 && LIMC < LIMCMAX / 8 &&
4323 : cache.last - cache.base > 10 * F.KC) goto START;
4324 : #endif
4325 121719 : old_need = need;
4326 121719 : if (!lambda)
4327 11 : { precpb = "bestappr"; PREC = myprecdbl(PREC, flag? C: NULL); continue; }
4328 121708 : if (!R)
4329 : { /* not full rank for units */
4330 25361 : if (!need)
4331 0 : { precpb = "regulator"; PREC = myprecdbl(PREC, flag? C: NULL); }
4332 25361 : continue;
4333 : }
4334 96347 : if (cache.last==old_cache) { need=1; continue; }
4335 96263 : old_cache = cache.last;
4336 96263 : h = ZM_det_triangular(W);
4337 96263 : if (DEBUGLEVEL) err_printf("\n#### Tentative class number: %Ps\n", h);
4338 96263 : i = compute_R(lambda, mulir(h,invhr), &L, &R);
4339 96260 : if (DEBUGLEVEL)
4340 : {
4341 0 : err_printf("\n");
4342 0 : timer_printf(&T, "computing regulator and check");
4343 : }
4344 96260 : switch(i)
4345 : {
4346 32192 : case fupb_RELAT:
4347 32192 : need = 1; /* not enough relations */
4348 32192 : continue;
4349 62 : case fupb_PRECI: /* prec problem unless we cheat on Bach constant */
4350 62 : if ((precdouble&7) == 7 && LIMC <= LIMCMAX/2) goto START;
4351 62 : precpb = "compute_R"; PREC = myprecdbl(PREC, flag? C: NULL);
4352 62 : continue;
4353 : }
4354 : /* DONE */
4355 :
4356 64006 : if (F.KCZ2 > F.KCZ)
4357 : {
4358 7 : if (F.sfb_chg && !subFB_change(&F)) goto START;
4359 7 : if (!be_honest(&F, nf, auts, fact)) goto START;
4360 7 : if (DEBUGLEVEL) timer_printf(&T, "to be honest");
4361 : }
4362 64006 : F.KCZ2 = 0; /* be honest only once */
4363 :
4364 : /* fundamental units */
4365 : {
4366 64006 : GEN AU, CU, U, v = extract_full_lattice(L); /* L may be large */
4367 64007 : CU = NULL;
4368 64007 : if (v) { A = vecpermute(A, v); L = vecpermute(L, v); }
4369 : /* arch. components of fund. units */
4370 64007 : U = ZM_lll(L, 0.99, LLL_IM);
4371 64008 : U = ZM_mul(U, lll(RgM_ZM_mul(real_i(A), U)));
4372 64008 : if (DEBUGLEVEL) timer_printf(&T, "units LLL");
4373 64008 : AU = RgM_ZM_mul(A, U);
4374 64008 : A = cleanarchunit(AU, N, NULL, PREC);
4375 64008 : if (RU > 1 /* if there are fund units, test we have correct regulator */
4376 48797 : && (!A || lg(A) < RU || expo(subrr(get_regulator(A), R)) > -1))
4377 7 : {
4378 7 : long add = nbits2extraprec( gexpo(AU) + 64 ) - gprecision(AU);
4379 7 : long t = maxss(PREC * 0.15, add);
4380 7 : if (!A && DEBUGLEVEL) err_printf("### Incorrect units lognorm");
4381 7 : precpb = "cleanarch"; PREC += maxss(t, EXTRAPREC64); continue;
4382 : }
4383 64000 : if (flag)
4384 : {
4385 62425 : long l = lgcols(C) - RU;
4386 : REL_t *rel;
4387 62425 : SUnits = cgetg(l, t_COL);
4388 1011070 : for (rel = cache.base+1, i = 1; i < l; i++,rel++)
4389 948644 : set_rel_alpha(rel, auts, SUnits, i);
4390 62426 : if (RU > 1)
4391 : {
4392 47712 : GEN c = v? vecpermute(C,v): vecslice(C,1,zc);
4393 47712 : CU = ZM_mul(rowslice(c, RU+1, nbrows(c)), U);
4394 : }
4395 : }
4396 64001 : if (DEBUGLEVEL) err_printf("\n#### Computing fundamental units\n");
4397 64001 : fu = getfu(nf, &A, CU? &U: NULL, PREC);
4398 64001 : CU = CU? ZM_mul(CU, U): cgetg(1, t_MAT);
4399 64001 : if (DEBUGLEVEL) timer_printf(&T, "getfu");
4400 64001 : Ce = vecslice(C, zc+1, lg(C)-1);
4401 64001 : if (flag) SUnits = mkvec4(SUnits, CU, rowslice(Ce, RU+1, nbrows(Ce)),
4402 : utoipos(LIMC));
4403 : }
4404 : /* class group generators */
4405 64001 : if (flag) Ce = rowslice(Ce, 1, RU);
4406 64001 : C0 = Ce; Ce = cleanarch(Ce, N, NULL, PREC);
4407 64000 : if (!Ce) {
4408 0 : long add = nbits2extraprec( gexpo(C0) + 64 ) - gprecision(C0);
4409 0 : precpb = "cleanarch"; PREC += maxss(add, 1);
4410 : }
4411 64000 : if (DEBUGLEVEL) timer_printf(&T, "cleanarch");
4412 121717 : } while (need || precpb);
4413 :
4414 64000 : Vbase = vecpermute(F.LP, F.perm);
4415 64001 : if (!fu) fu = cgetg(1, t_MAT);
4416 64001 : if (!SUnits) SUnits = gen_1;
4417 64001 : clg1 = class_group_gen(nf,W,Ce,Vbase,PREC, &clg2);
4418 64001 : res = mkvec5(clg1, R, SUnits, zu, fu);
4419 64001 : res = buchall_end(nf,res,clg2,W,B,A,Ce,Vbase);
4420 64001 : delete_FB(&F);
4421 64001 : res = gc_GEN(av0, res);
4422 64001 : if (flag) obj_insert_shallow(res, MATAL, cgetg(1,t_VEC));
4423 64001 : if (nfisclone) gunclone(nf);
4424 64001 : delete_cache(&cache);
4425 64001 : free_GRHcheck(&GRHcheck);
4426 64001 : return res;
4427 : }
|