Book a Demo!
CoCalc Logo Icon
StoreFeaturesDocsShareSupportNewsAboutPoliciesSign UpSign In
Download

Testing latest pari + WASM + node.js... and it works?! Wow.

28479 views
License: GPL3
ubuntu2004
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_quadclassunit
18
19
/*******************************************************************/
20
/* */
21
/* CLASS GROUP AND REGULATOR (McCURLEY, BUCHMANN) */
22
/* QUADRATIC FIELDS */
23
/* */
24
/*******************************************************************/
25
/* For largeprime() hashtable. Note that hashed pseudoprimes are odd (unless
26
* 2 | index), hence the low order bit is not useful. So we hash
27
* HASHBITS bits starting at bit 1, not bit 0 */
28
#define HASHBITS 10
29
static const long HASHT = 1L << HASHBITS;
30
31
static long
32
hash(long q) { return (q & ((1L << (HASHBITS+1)) - 1)) >> 1; }
33
#undef HASHBITS
34
35
/* See buch2.c:
36
* B->subFB contains split p such that \prod p > sqrt(B->Disc)
37
* B->powsubFB contains powers of forms in B->subFB */
38
#define RANDOM_BITS 4
39
static const long CBUCH = (1L<<RANDOM_BITS)-1;
40
41
struct buch_quad
42
{
43
ulong limhash;
44
long KC, KC2, PRECREG;
45
long *primfact, *exprimfact, **hashtab;
46
GEN FB, numFB;
47
GEN powsubFB, vperm, subFB, badprim;
48
struct qfr_data *q;
49
};
50
51
/*******************************************************************/
52
/* */
53
/* Routines related to binary quadratic forms (for internal use) */
54
/* */
55
/*******************************************************************/
56
/* output canonical representative wrt projection Cl^+ --> Cl (a > 0) */
57
static GEN
58
qfr3_canon(GEN x, struct qfr_data *S)
59
{
60
GEN a = gel(x,1), c = gel(x,3);
61
if (signe(a) < 0) {
62
if (absequalii(a,c)) return qfr3_rho(x, S);
63
setsigne(a, 1);
64
setsigne(c,-1);
65
}
66
return x;
67
}
68
static GEN
69
qfr3_canon_safe(GEN x, struct qfr_data *S)
70
{
71
GEN a = gel(x,1), c = gel(x,3);
72
if (signe(a) < 0) {
73
if (absequalii(a,c)) return qfr3_rho(x, S);
74
gel(x,1) = negi(a);
75
gel(x,3) = negi(c);
76
}
77
return x;
78
}
79
static GEN
80
qfr5_canon(GEN x, struct qfr_data *S)
81
{
82
GEN a = gel(x,1), c = gel(x,3);
83
if (signe(a) < 0) {
84
if (absequalii(a,c)) return qfr5_rho(x, S);
85
setsigne(a, 1);
86
setsigne(c,-1);
87
}
88
return x;
89
}
90
static GEN
91
QFR5_comp(GEN x,GEN y, struct qfr_data *S)
92
{ return qfr5_canon(qfr5_comp(x,y,S), S); }
93
static GEN
94
QFR3_comp(GEN x, GEN y, struct qfr_data *S)
95
{ return qfr3_canon(qfr3_comp(x,y,S), S); }
96
97
/* compute rho^n(x) */
98
static GEN
99
qfr5_rho_pow(GEN x, long n, struct qfr_data *S)
100
{
101
long i;
102
pari_sp av = avma;
103
for (i=1; i<=n; i++)
104
{
105
x = qfr5_rho(x,S);
106
if (gc_needed(av,1))
107
{
108
if(DEBUGMEM>1) pari_warn(warnmem,"qfr5_rho_pow");
109
x = gerepilecopy(av, x);
110
}
111
}
112
return gerepilecopy(av, x);
113
}
114
115
static GEN
116
qfr5_pf(struct qfr_data *S, long p, long prec)
117
{
118
GEN y = primeform_u(S->D,p);
119
return qfr5_canon(qfr5_red(qfr_to_qfr5(y,prec), S), S);
120
}
121
122
static GEN
123
qfr3_pf(struct qfr_data *S, long p)
124
{
125
GEN y = primeform_u(S->D,p);
126
return qfr3_canon(qfr3_red(y, S), S);
127
}
128
129
#define qfi_pf primeform_u
130
131
/* Warning: ex[0] not set in general */
132
static GEN
133
init_form(struct buch_quad *B, GEN ex,
134
GEN (*comp)(GEN,GEN,struct qfr_data *S))
135
{
136
long i, l = lg(B->powsubFB);
137
GEN F = NULL;
138
for (i=1; i<l; i++)
139
if (ex[i])
140
{
141
GEN t = gmael(B->powsubFB,i,ex[i]);
142
F = F? comp(F, t, B->q): t;
143
}
144
return F;
145
}
146
static GEN
147
qfr5_factorback(struct buch_quad *B, GEN ex) { return init_form(B, ex, &QFR5_comp); }
148
149
static GEN
150
QFI_comp(GEN x, GEN y, struct qfr_data *S) { (void)S; return qfbcomp_i(x,y); }
151
152
static GEN
153
qfi_factorback(struct buch_quad *B, GEN ex) { return init_form(B, ex, &QFI_comp); }
154
155
static GEN
156
random_form(struct buch_quad *B, GEN ex,
157
GEN (*comp)(GEN,GEN, struct qfr_data *S))
158
{
159
long i, l = lg(ex);
160
pari_sp av = avma;
161
GEN F;
162
for(;;)
163
{
164
for (i=1; i<l; i++) ex[i] = random_bits(RANDOM_BITS);
165
if ((F = init_form(B, ex, comp))) return F;
166
set_avma(av);
167
}
168
}
169
static GEN
170
qfr3_random(struct buch_quad *B,GEN ex){ return random_form(B, ex, &QFR3_comp); }
171
static GEN
172
qfi_random(struct buch_quad *B,GEN ex) { return random_form(B, ex, &QFI_comp); }
173
174
/*******************************************************************/
175
/* */
176
/* Common subroutines */
177
/* */
178
/*******************************************************************/
179
long
180
bnf_increase_LIMC(long LIMC, long LIMCMAX)
181
{
182
if (LIMC >= LIMCMAX) pari_err_BUG("Buchmann's algorithm");
183
if (LIMC <= LIMCMAX/40) /* cbach <= 0.3 */
184
LIMC *= 2;
185
else if (LIMCMAX < 60) /* \Delta_K <= 9 */
186
LIMC++;
187
else
188
LIMC += LIMCMAX / 60; /* cbach += 0.2 */
189
if (LIMC > LIMCMAX) LIMC = LIMCMAX;
190
return LIMC;
191
}
192
193
/* Is |q| <= p ? */
194
static int
195
isless_iu(GEN q, ulong p) {
196
long l = lgefint(q);
197
return l==2 || (l == 3 && uel(q,2) <= p);
198
}
199
200
static long
201
factorquad(struct buch_quad *B, GEN f, long nFB, ulong limp)
202
{
203
ulong X;
204
long i, lo = 0;
205
GEN x = gel(f,1), FB = B->FB, P = B->primfact, E = B->exprimfact;
206
207
for (i=1; lgefint(x) > 3; i++)
208
{
209
ulong p = uel(FB,i), r;
210
GEN q = absdiviu_rem(x, p, &r);
211
if (!r)
212
{
213
long k = 0;
214
do { k++; x = q; q = absdiviu_rem(x, p, &r); } while (!r);
215
lo++; P[lo] = p; E[lo] = k;
216
}
217
if (isless_iu(q,p)) {
218
if (lgefint(x) == 3) { X = uel(x,2); goto END; }
219
return 0;
220
}
221
if (i == nFB) return 0;
222
}
223
X = uel(x,2);
224
if (X == 1) { P[0] = 0; return 1; }
225
for (;; i++)
226
{ /* single precision affair, split for efficiency */
227
ulong p = uel(FB,i);
228
ulong q = X / p, r = X % p; /* gcc makes a single div */
229
if (!r)
230
{
231
long k = 0;
232
do { k++; X = q; q = X / p; r = X % p; } while (!r);
233
lo++; P[lo] = p; E[lo] = k;
234
}
235
if (q <= p) break;
236
if (i == nFB) return 0;
237
}
238
END:
239
if (X > B->limhash) return 0;
240
if (X != 1 && X <= limp)
241
{
242
if (B->badprim && ugcdui(X, B->badprim) > 1) return 0;
243
lo++; P[lo] = X; E[lo] = 1;
244
X = 1;
245
}
246
P[0] = lo; return X;
247
}
248
249
/* Check for a "large prime relation" involving q; q may not be prime */
250
static long *
251
largeprime(struct buch_quad *B, long q, GEN ex, long np, long nrho)
252
{
253
const long hashv = hash(q);
254
long *pt, i, l = lg(B->subFB);
255
256
for (pt = B->hashtab[hashv]; ; pt = (long*) pt[0])
257
{
258
if (!pt)
259
{
260
pt = (long*) pari_malloc((l+3) * sizeof(long));
261
*pt++ = nrho; /* nrho = pt[-3] */
262
*pt++ = np; /* np = pt[-2] */
263
*pt++ = q; /* q = pt[-1] */
264
pt[0] = (long)B->hashtab[hashv];
265
for (i=1; i<l; i++) pt[i]=ex[i];
266
B->hashtab[hashv]=pt; return NULL;
267
}
268
if (pt[-1] == q) break;
269
}
270
for(i=1; i<l; i++)
271
if (pt[i] != ex[i]) return pt;
272
return (pt[-2]==np)? NULL: pt;
273
}
274
275
static void
276
clearhash(long **hash)
277
{
278
long *pt;
279
long i;
280
for (i=0; i<HASHT; i++) {
281
for (pt = hash[i]; pt; ) {
282
void *z = (void*)(pt - 3);
283
pt = (long*) pt[0]; pari_free(z);
284
}
285
hash[i] = NULL;
286
}
287
}
288
289
/* last prime stored */
290
ulong
291
GRH_last_prime(GRHcheck_t *S) { return (S->primes + S->nprimes-1)->p; }
292
/* ensure that S->primes can hold at least nb primes */
293
void
294
GRH_ensure(GRHcheck_t *S, long nb)
295
{
296
if (S->maxprimes <= nb)
297
{
298
do S->maxprimes *= 2; while (S->maxprimes <= nb);
299
pari_realloc_ip((void**)&S->primes, S->maxprimes*sizeof(*S->primes));
300
}
301
}
302
/* cache data for all primes up to the LIM */
303
static void
304
cache_prime_quad(GRHcheck_t *S, ulong LIM, GEN D)
305
{
306
GRHprime_t *pr;
307
long nb;
308
309
if (S->limp >= LIM) return;
310
nb = (long)primepi_upper_bound((double)LIM); /* #{p <= LIM} <= nb */
311
GRH_ensure(S, nb+1); /* room for one extra prime */
312
for (pr = S->primes + S->nprimes;;)
313
{
314
ulong p = u_forprime_next(&(S->P));
315
pr->p = p;
316
pr->logp = log((double)p);
317
pr->dec = (GEN)kroiu(D,p);
318
S->nprimes++;
319
pr++;
320
/* store up to nextprime(LIM) included */
321
if (p >= LIM) { S->limp = p; break; }
322
}
323
}
324
325
static GEN
326
compute_invresquad(GRHcheck_t *S, long LIMC)
327
{
328
pari_sp av = avma;
329
GEN invres = real_1(DEFAULTPREC);
330
double limp = log((double)LIMC) / 2;
331
GRHprime_t *pr;
332
long i;
333
for (pr = S->primes, i = S->nprimes; i > 0; pr++, i--)
334
{
335
long s = (long)pr->dec;
336
if (s)
337
{
338
ulong p = pr->p;
339
if (s > 0 || pr->logp <= limp)
340
/* Both p and P contribute */
341
invres = mulur(p - s, divru(invres, p));
342
else if (s<0)
343
/* Only p contributes */
344
invres = mulur(p, divru(invres, p - 1));
345
}
346
}
347
return gerepileuptoleaf(av, invres);
348
}
349
350
/* p | conductor of order of disc D ? */
351
static int
352
is_bad(GEN D, ulong p)
353
{
354
pari_sp av = avma;
355
if (p == 2)
356
{
357
long r = mod16(D) >> 1;
358
if (r && signe(D) < 0) r = 8-r;
359
return (r < 4);
360
}
361
return gc_bool(av, dvdii(D, sqru(p))); /* p^2 | D ? */
362
}
363
364
/* returns the n-th suitable ideal for the factorbase */
365
static long
366
nthidealquad(GEN D, long n)
367
{
368
pari_sp av = avma;
369
forprime_t S;
370
ulong p;
371
(void)u_forprime_init(&S, 2, ULONG_MAX);
372
while ((p = u_forprime_next(&S)) && n > 0)
373
if (!is_bad(D, p) && kroiu(D, p) >= 0) n--;
374
return gc_long(av, p);
375
}
376
377
static int
378
quadGRHchk(GEN D, GRHcheck_t *S, ulong LIMC)
379
{
380
double logC = log((double)LIMC), SA = 0, SB = 0;
381
long i;
382
383
cache_prime_quad(S, LIMC, D);
384
for (i = 0;; i++)
385
{
386
GRHprime_t *pr = S->primes+i;
387
ulong p = pr->p;
388
long M;
389
double logNP, q, A, B;
390
if (p > LIMC) break;
391
if ((long)pr->dec < 0)
392
{
393
logNP = 2 * pr->logp;
394
q = 1/(double)p;
395
}
396
else
397
{
398
logNP = pr->logp;
399
q = 1/sqrt((double)p);
400
}
401
A = logNP * q; B = logNP * A; M = (long)(logC/logNP);
402
if (M > 1)
403
{
404
double inv1_q = 1 / (1-q);
405
A *= (1 - pow(q, (double) M)) * inv1_q;
406
B *= (1 - pow(q, (double) M)*(M+1 - M*q)) * inv1_q * inv1_q;
407
}
408
if ((long)pr->dec>0) { SA += 2*A;SB += 2*B; } else { SA += A; SB += B; }
409
if (p == LIMC) break;
410
}
411
return GRHok(S, logC, SA, SB);
412
}
413
414
/* C2 >= C1; create B->FB, B->numFB; set B->badprim. Return L(kro_D, 1) */
415
static void
416
FBquad(struct buch_quad *B, ulong C2, ulong C1, GRHcheck_t *S)
417
{
418
GEN D = B->q->D;
419
long i;
420
pari_sp av;
421
GRHprime_t *pr;
422
423
cache_prime_quad(S, C2, D);
424
pr = S->primes;
425
B->numFB = cgetg(C2+1, t_VECSMALL);
426
B->FB = cgetg(C2+1, t_VECSMALL);
427
av = avma;
428
B->KC = 0; i = 0;
429
B->badprim = gen_1;
430
for (;; pr++) /* p <= C2 */
431
{
432
ulong p = pr->p;
433
if (!B->KC && p > C1) B->KC = i;
434
if (p > C2) break;
435
switch ((long)pr->dec)
436
{
437
case -1: break; /* inert */
438
case 0: /* ramified */
439
if (is_bad(D, p)) { B->badprim = muliu(B->badprim, p); break; }
440
/* fall through */
441
default: /* split */
442
i++; B->numFB[p] = i; B->FB[i] = p; break;
443
}
444
if (p == C2)
445
{
446
if (!B->KC) B->KC = i;
447
break;
448
}
449
}
450
B->KC2 = i;
451
setlg(B->FB, B->KC2+1);
452
if (B->badprim != gen_1)
453
B->badprim = gerepileuptoint(av, B->badprim);
454
else
455
{
456
B->badprim = NULL; set_avma(av);
457
}
458
}
459
460
/* create B->vperm, return B->subFB */
461
static GEN
462
subFBquad(struct buch_quad *B, GEN D, double PROD, long minSFB)
463
{
464
long i, j, lgsub = 1, ino = 1, lv = B->KC+1;
465
double prod = 1.;
466
pari_sp av;
467
GEN no;
468
469
B->vperm = cgetg(lv, t_VECSMALL);
470
av = avma;
471
no = cgetg(lv, t_VECSMALL);
472
for (j = 1; j < lv; j++)
473
{
474
ulong p = uel(B->FB,j);
475
if (!umodiu(D, p)) no[ino++] = j; /* ramified */
476
else
477
{
478
B->vperm[lgsub++] = j;
479
prod *= p;
480
if (lgsub > minSFB && prod > PROD) break;
481
}
482
}
483
/* lgsub >= 1 otherwise quadGRHchk is false */
484
i = lgsub;
485
for (j = 1; j < ino;i++,j++) B->vperm[i] = no[j];
486
for ( ; i < lv; i++) B->vperm[i] = i;
487
no = gclone(vecslice(B->vperm, 1, lgsub-1));
488
set_avma(av); return no;
489
}
490
491
/* assume n >= 1, x[i][j] = B->subFB[i]^j, for j = 1..n */
492
static GEN
493
powsubFBquad(struct buch_quad *B, long n)
494
{
495
pari_sp av = avma;
496
long i,j, l = lg(B->subFB);
497
GEN F, y, x = cgetg(l, t_VEC), D = B->q->D;
498
499
if (B->PRECREG) /* real */
500
{
501
for (i=1; i<l; i++)
502
{
503
F = qfr5_pf(B->q, B->FB[B->subFB[i]], B->PRECREG);
504
y = cgetg(n+1, t_VEC); gel(x,i) = y;
505
gel(y,1) = F;
506
for (j=2; j<=n; j++) gel(y,j) = QFR5_comp(gel(y,j-1), F, B->q);
507
}
508
}
509
else /* imaginary */
510
{
511
for (i=1; i<l; i++)
512
{
513
F = qfi_pf(D, B->FB[B->subFB[i]]);
514
y = cgetg(n+1, t_VEC); gel(x,i) = y;
515
gel(y,1) = F;
516
for (j=2; j<=n; j++) gel(y,j) = qfbcomp_i(gel(y,j-1), F);
517
}
518
}
519
x = gclone(x); set_avma(av); return x;
520
}
521
522
static void
523
sub_fact(struct buch_quad *B, GEN col, GEN F)
524
{
525
GEN b = gel(F,2);
526
long i;
527
for (i=1; i<=B->primfact[0]; i++)
528
{
529
ulong p = B->primfact[i], k = B->numFB[p];
530
long e = B->exprimfact[i];
531
if (umodiu(b, p<<1) > p) e = -e;
532
col[k] -= e;
533
}
534
}
535
static void
536
add_fact(struct buch_quad *B, GEN col, GEN F)
537
{
538
GEN b = gel(F,2);
539
long i;
540
for (i=1; i<=B->primfact[0]; i++)
541
{
542
ulong p = B->primfact[i], k = B->numFB[p];
543
long e = B->exprimfact[i];
544
if (umodiu(b, p<<1) > p) e = -e;
545
col[k] += e;
546
}
547
}
548
549
static GEN
550
get_clgp(struct buch_quad *B, GEN W, GEN *ptD)
551
{
552
GEN res, init, u1, D = ZM_snf_group(W,NULL,&u1);
553
long i, j, l = lg(W), c = lg(D);
554
555
res=cgetg(c,t_VEC); init = cgetg(l,t_VEC);
556
for (i=1; i<l; i++) gel(init,i) = primeform_u(B->q->D, B->FB[B->vperm[i]]);
557
for (j=1; j<c; j++)
558
{
559
GEN g = NULL;
560
if (signe(B->q->D) > 0)
561
{
562
for (i=1; i<l; i++)
563
{
564
GEN t, u = gcoeff(u1,i,j);
565
if (!signe(u)) continue;
566
t = qfr3_pow(gel(init,i), u, B->q);
567
g = g? qfr3_comp(g, t, B->q): t;
568
}
569
g = qfr3_to_qfr(qfr3_canon_safe(qfr3_red(g, B->q), B->q), B->q->D);
570
}
571
else
572
{
573
for (i=1; i<l; i++)
574
{
575
GEN t, u = gcoeff(u1,i,j);
576
if (!signe(u)) continue;
577
t = powgi(gel(init,i), u);
578
g = g? qfbcomp_i(g, t): t;
579
}
580
}
581
gel(res,j) = g;
582
}
583
*ptD = D; return res;
584
}
585
586
static long
587
trivial_relations(struct buch_quad *B, GEN mat, GEN C)
588
{
589
long i, j = 0;
590
GEN col, D = B->q->D;
591
for (i = 1; i <= B->KC; i++)
592
{ /* ramified prime ==> trivial relation */
593
if (umodiu(D, B->FB[i])) continue;
594
col = zero_zv(B->KC);
595
col[i] = 2; j++;
596
gel(mat,j) = col;
597
gel(C,j) = gen_0;
598
}
599
return j;
600
}
601
602
static void
603
dbg_all(pari_timer *T, const char *phase, long s, long n)
604
{
605
err_printf("\n");
606
timer_printf(T, "%s rel [#rel/#test = %ld/%ld]", phase,s,n);
607
}
608
609
/* Imaginary Quadratic fields */
610
611
static void
612
imag_relations(struct buch_quad *B, long need, long *pc, ulong LIMC, GEN mat)
613
{
614
pari_timer T;
615
long lgsub = lg(B->subFB), current = *pc, nbtest = 0, s = 0;
616
long i, fpc;
617
pari_sp av;
618
GEN col, form, ex = cgetg(lgsub, t_VECSMALL);
619
620
if (!current) current = 1;
621
if (DEBUGLEVEL>2) timer_start(&T);
622
av = avma;
623
for(;;)
624
{
625
if (s >= need) break;
626
set_avma(av);
627
form = qfi_random(B,ex);
628
form = qfbcomp_i(form, qfi_pf(B->q->D, B->FB[current]));
629
nbtest++; fpc = factorquad(B,form,B->KC,LIMC);
630
if (!fpc)
631
{
632
if (DEBUGLEVEL>3) err_printf(".");
633
if ((nbtest & 0xff) == 0 && ++current > B->KC) current = 1;
634
continue;
635
}
636
if (fpc > 1)
637
{
638
long *fpd = largeprime(B,fpc,ex,current,0);
639
ulong b1, b2, p;
640
GEN form2;
641
if (!fpd)
642
{
643
if (DEBUGLEVEL>3) err_printf(".");
644
continue;
645
}
646
form2 = qfbcomp_i(qfi_factorback(B,fpd), qfi_pf(B->q->D, B->FB[fpd[-2]]));
647
p = fpc << 1;
648
b1 = umodiu(gel(form2,2), p);
649
b2 = umodiu(gel(form,2), p);
650
if (b1 != b2 && b1+b2 != p) continue;
651
652
col = gel(mat,++s);
653
add_fact(B,col, form);
654
(void)factorquad(B,form2,B->KC,LIMC);
655
if (b1==b2)
656
{
657
for (i=1; i<lgsub; i++) col[B->subFB[i]] += fpd[i]-ex[i];
658
sub_fact(B, col, form2); col[fpd[-2]]++;
659
}
660
else
661
{
662
for (i=1; i<lgsub; i++) col[B->subFB[i]] += -fpd[i]-ex[i];
663
add_fact(B, col, form2); col[fpd[-2]]--;
664
}
665
if (DEBUGLEVEL>2) err_printf(" %ldP",s);
666
}
667
else
668
{
669
col = gel(mat,++s);
670
for (i=1; i<lgsub; i++) col[B->subFB[i]] = -ex[i];
671
add_fact(B, col, form);
672
if (DEBUGLEVEL>2) err_printf(" %ld",s);
673
}
674
col[current]--;
675
if (++current > B->KC) current = 1;
676
}
677
if (DEBUGLEVEL>2) dbg_all(&T, "random", s, nbtest);
678
*pc = current;
679
}
680
681
static int
682
imag_be_honest(struct buch_quad *B)
683
{
684
long p, fpc, s = B->KC, nbtest = 0;
685
GEN F, ex = cgetg(lg(B->subFB), t_VECSMALL);
686
pari_sp av = avma;
687
688
while (s<B->KC2)
689
{
690
p = B->FB[s+1]; if (DEBUGLEVEL>2) err_printf(" %ld",p);
691
F = qfbcomp_i(qfi_pf(B->q->D, p), qfi_random(B, ex));
692
fpc = factorquad(B,F,s,p-1);
693
if (fpc == 1) { nbtest=0; s++; }
694
else
695
if (++nbtest > 40) return 0;
696
set_avma(av);
697
}
698
return 1;
699
}
700
701
/* Real Quadratic fields */
702
703
static void
704
real_relations(struct buch_quad *B, long need, long *pc, long lim, ulong LIMC, GEN mat, GEN C)
705
{
706
pari_timer T;
707
long lgsub = lg(B->subFB), prec = B->PRECREG, current = *pc, nbtest=0, s=0;
708
long i, fpc, endcycle, rhoacc, rho;
709
/* in a 2nd phase, don't include FB[current] but run along the cyle
710
* ==> get more units */
711
int first = (current == 0);
712
pari_sp av, av1;
713
GEN d, col, form, form0, form1, ex = cgetg(lgsub, t_VECSMALL);
714
715
if (DEBUGLEVEL>2) timer_start(&T);
716
if (!current) current = 1;
717
if (lim > need) lim = need;
718
av = avma;
719
for(;;)
720
{
721
if (s >= need) break;
722
if (first && s >= lim) {
723
first = 0;
724
if (DEBUGLEVEL>2) dbg_all(&T, "initial", s, nbtest);
725
}
726
set_avma(av); form = qfr3_random(B, ex);
727
if (!first)
728
form = QFR3_comp(form, qfr3_pf(B->q, B->FB[current]), B->q);
729
av1 = avma;
730
form0 = form; form1 = NULL;
731
endcycle = rhoacc = 0;
732
rho = -1;
733
734
CYCLE:
735
if (endcycle || rho > 5000)
736
{
737
if (++current > B->KC) current = 1;
738
continue;
739
}
740
if (gc_needed(av,1))
741
{
742
if(DEBUGMEM>1) pari_warn(warnmem,"real_relations");
743
gerepileall(av1, form1? 2: 1, &form, &form1);
744
}
745
if (rho < 0) rho = 0; /* first time in */
746
else
747
{
748
form = qfr3_rho(form, B->q); rho++;
749
rhoacc++;
750
if (first)
751
endcycle = (absequalii(gel(form,1),gel(form0,1))
752
&& equalii(gel(form,2),gel(form0,2)));
753
else
754
{
755
if (absequalii(gel(form,1), gel(form,3))) /* a = -c */
756
{
757
if (absequalii(gel(form,1),gel(form0,1)) &&
758
equalii(gel(form,2),gel(form0,2))) continue;
759
form = qfr3_rho(form, B->q); rho++;
760
rhoacc++;
761
}
762
else
763
{ setsigne(form[1],1); setsigne(form[3],-1); }
764
if (equalii(gel(form,1),gel(form0,1)) &&
765
equalii(gel(form,2),gel(form0,2))) continue;
766
}
767
}
768
nbtest++; fpc = factorquad(B,form,B->KC,LIMC);
769
if (!fpc)
770
{
771
if (DEBUGLEVEL>3) err_printf(".");
772
goto CYCLE;
773
}
774
if (fpc > 1)
775
{ /* look for Large Prime relation */
776
long *fpd = largeprime(B,fpc,ex,first? 0: current,rhoacc);
777
ulong b1, b2, p;
778
GEN form2;
779
if (!fpd)
780
{
781
if (DEBUGLEVEL>3) err_printf(".");
782
goto CYCLE;
783
}
784
if (!form1)
785
{
786
form1 = qfr5_factorback(B,ex);
787
if (!first)
788
form1 = QFR5_comp(form1, qfr5_pf(B->q, B->FB[current], prec), B->q);
789
}
790
form1 = qfr5_rho_pow(form1, rho, B->q);
791
rho = 0;
792
793
form2 = qfr5_factorback(B,fpd);
794
if (fpd[-2])
795
form2 = QFR5_comp(form2, qfr5_pf(B->q, B->FB[fpd[-2]], prec), B->q);
796
form2 = qfr5_rho_pow(form2, fpd[-3], B->q);
797
if (!absequalii(gel(form2,1),gel(form2,3)))
798
{
799
setsigne(form2[1], 1);
800
setsigne(form2[3],-1);
801
}
802
p = fpc << 1;
803
b1 = umodiu(gel(form2,2), p);
804
b2 = umodiu(gel(form1,2), p);
805
if (b1 != b2 && b1+b2 != p) goto CYCLE;
806
807
col = gel(mat,++s);
808
add_fact(B, col, form1);
809
(void)factorquad(B,form2,B->KC,LIMC);
810
if (b1==b2)
811
{
812
for (i=1; i<lgsub; i++) col[B->subFB[i]] += fpd[i]-ex[i];
813
sub_fact(B,col, form2);
814
if (fpd[-2]) col[fpd[-2]]++;
815
d = qfr5_dist(subii(gel(form1,4),gel(form2,4)),
816
divrr(gel(form1,5),gel(form2,5)), prec);
817
}
818
else
819
{
820
for (i=1; i<lgsub; i++) col[B->subFB[i]] += -fpd[i]-ex[i];
821
add_fact(B, col, form2);
822
if (fpd[-2]) col[fpd[-2]]--;
823
d = qfr5_dist(addii(gel(form1,4),gel(form2,4)),
824
mulrr(gel(form1,5),gel(form2,5)), prec);
825
}
826
if (DEBUGLEVEL>2) err_printf(" %ldP",s);
827
}
828
else
829
{ /* standard relation */
830
if (!form1)
831
{
832
form1 = qfr5_factorback(B, ex);
833
if (!first)
834
form1 = QFR5_comp(form1, qfr5_pf(B->q, B->FB[current], prec), B->q);
835
}
836
form1 = qfr5_rho_pow(form1, rho, B->q);
837
rho = 0;
838
839
col = gel(mat,++s);
840
for (i=1; i<lgsub; i++) col[B->subFB[i]] = -ex[i];
841
add_fact(B, col, form1);
842
d = qfr5_dist(gel(form1,4), gel(form1,5), prec);
843
if (DEBUGLEVEL>2) err_printf(" %ld",s);
844
}
845
affrr(d, gel(C,s));
846
if (first)
847
{
848
if (s >= lim) continue;
849
goto CYCLE;
850
}
851
else
852
{
853
col[current]--;
854
if (++current > B->KC) current = 1;
855
}
856
}
857
if (DEBUGLEVEL>2) dbg_all(&T, "random", s, nbtest);
858
*pc = current;
859
}
860
861
static int
862
real_be_honest(struct buch_quad *B)
863
{
864
long p, fpc, s = B->KC, nbtest = 0;
865
GEN F,F0, ex = cgetg(lg(B->subFB), t_VECSMALL);
866
pari_sp av = avma;
867
868
while (s<B->KC2)
869
{
870
p = B->FB[s+1]; if (DEBUGLEVEL>2) err_printf(" %ld",p);
871
F = QFR3_comp(qfr3_random(B, ex), qfr3_pf(B->q, p), B->q);
872
for (F0 = F;;)
873
{
874
fpc = factorquad(B,F,s,p-1);
875
if (fpc == 1) { nbtest=0; s++; break; }
876
if (++nbtest > 40) return 0;
877
F = qfr3_canon(qfr3_rho(F, B->q), B->q);
878
if (equalii(gel(F,1),gel(F0,1))
879
&& equalii(gel(F,2),gel(F0,2))) break;
880
}
881
set_avma(av);
882
}
883
return 1;
884
}
885
886
static GEN
887
gcdreal(GEN a,GEN b)
888
{
889
if (!signe(a)) return mpabs_shallow(b);
890
if (!signe(b)) return mpabs_shallow(a);
891
if (typ(a)==t_INT)
892
{
893
if (typ(b)==t_INT) return gcdii(a,b);
894
a = itor(a, lg(b));
895
}
896
else if (typ(b)==t_INT)
897
b = itor(b, lg(a));
898
if (expo(a)<-5) return absr(b);
899
if (expo(b)<-5) return absr(a);
900
a = absr(a); b = absr(b);
901
while (expo(b) >= -5 && signe(b))
902
{
903
long e;
904
GEN r, q = gcvtoi(divrr(a,b),&e);
905
if (e > 0) return NULL;
906
r = subrr(a, mulir(q,b)); a = b; b = r;
907
}
908
return mpabs_shallow(a);
909
}
910
911
static int
912
get_R(struct buch_quad *B, GEN C, long sreg, GEN z, GEN *ptR)
913
{
914
GEN R = gen_1;
915
double c;
916
long i;
917
918
if (B->PRECREG)
919
{
920
R = mpabs_shallow(gel(C,1));
921
for (i=2; i<=sreg; i++)
922
{
923
R = gcdreal(gel(C,i), R);
924
if (!R) return fupb_PRECI;
925
}
926
if (gexpo(R) <= -3)
927
{
928
if (DEBUGLEVEL>2) err_printf("regulator is zero.\n");
929
return fupb_RELAT;
930
}
931
if (DEBUGLEVEL>2) err_printf("#### Tentative regulator: %Ps\n",R);
932
}
933
c = gtodouble(gmul(z, R));
934
if (c < 0.8 || c > 1.3) return fupb_RELAT;
935
*ptR = R; return fupb_NONE;
936
}
937
938
static int
939
quad_be_honest(struct buch_quad *B)
940
{
941
int r;
942
if (B->KC2 <= B->KC) return 1;
943
if (DEBUGLEVEL>2)
944
err_printf("be honest for primes from %ld to %ld\n", B->FB[B->KC+1],B->FB[B->KC2]);
945
r = B->PRECREG? real_be_honest(B): imag_be_honest(B);
946
if (DEBUGLEVEL>2) err_printf("\n");
947
return r;
948
}
949
950
GEN
951
Buchquad(GEN D, double cbach, double cbach2, long prec)
952
{
953
const long MAXRELSUP = 7, SFB_MAX = 3;
954
pari_timer T;
955
pari_sp av0 = avma, av, av2;
956
const long RELSUP = 5;
957
long i, s, current, triv, sfb_trials, nrelsup, nreldep, need, nsubFB, minSFB;
958
ulong low, high, LIMC0, LIMC, LIMC2, LIMCMAX, cp;
959
GEN W, cyc, res, gen, dep, mat, C, extraC, B, R, invhr, h = NULL; /*-Wall*/
960
double drc, sdrc, lim, LOGD, LOGD2;
961
GRHcheck_t GRHcheck;
962
struct qfr_data q;
963
struct buch_quad BQ;
964
int FIRST = 1;
965
966
check_quaddisc(D, &s, /*junk*/&i, "Buchquad");
967
R = NULL; /* -Wall */
968
BQ.q = &q;
969
q.D = D;
970
if (s < 0)
971
{
972
if (abscmpiu(q.D,4) <= 0)
973
{
974
GEN z = cgetg(5,t_VEC);
975
gel(z,1) = gel(z,4) = gen_1; gel(z,2) = gel(z,3) = cgetg(1,t_VEC);
976
return z;
977
}
978
prec = BQ.PRECREG = 0;
979
} else {
980
BQ.PRECREG = maxss(prec+EXTRAPRECWORD, nbits2prec(2*expi(q.D) + 128));
981
}
982
if (DEBUGLEVEL>2) timer_start(&T);
983
BQ.primfact = new_chunk(100);
984
BQ.exprimfact = new_chunk(100);
985
BQ.hashtab = (long**) new_chunk(HASHT);
986
for (i=0; i<HASHT; i++) BQ.hashtab[i] = NULL;
987
988
drc = fabs(gtodouble(q.D));
989
LOGD = log(drc);
990
LOGD2 = LOGD * LOGD;
991
992
sdrc = lim = sqrt(drc);
993
if (!BQ.PRECREG) lim /= sqrt(3.);
994
cp = (ulong)exp(sqrt(LOGD*log(LOGD)/8.0));
995
if (cp < 20) cp = 20;
996
if (cbach > 6.) {
997
if (cbach2 < cbach) cbach2 = cbach;
998
cbach = 6.;
999
}
1000
if (cbach < 0.)
1001
pari_err_DOMAIN("Buchquad","Bach constant","<",gen_0,dbltor(cbach));
1002
av = avma;
1003
BQ.powsubFB = BQ.subFB = NULL;
1004
minSFB = (expi(D) > 15)? 3: 2;
1005
init_GRHcheck(&GRHcheck, 2, BQ.PRECREG? 2: 0, LOGD);
1006
high = low = LIMC0 = maxss((long)(cbach2*LOGD2), 1);
1007
LIMCMAX = (long)(6.*LOGD2);
1008
/* 97/1223 below to ensure a good enough approximation of residue */
1009
cache_prime_quad(&GRHcheck, expi(D) < 16 ? 97: 1223, D);
1010
while (!quadGRHchk(D, &GRHcheck, high))
1011
{
1012
low = high;
1013
high *= 2;
1014
}
1015
while (high - low > 1)
1016
{
1017
long test = (low+high)/2;
1018
if (quadGRHchk(D, &GRHcheck, test))
1019
high = test;
1020
else
1021
low = test;
1022
}
1023
if (high == LIMC0+1 && quadGRHchk(D, &GRHcheck, LIMC0))
1024
LIMC2 = LIMC0;
1025
else
1026
LIMC2 = high;
1027
if (LIMC2 > LIMCMAX) LIMC2 = LIMCMAX;
1028
LIMC0 = (long)(cbach*LOGD2);
1029
LIMC = cbach ? LIMC0 : LIMC2;
1030
LIMC = maxss(LIMC, nthidealquad(D, 2));
1031
1032
/* LIMC = Max(cbach*(log D)^2, exp(sqrt(log D loglog D) / 8)) */
1033
START:
1034
do
1035
{
1036
if (!FIRST) LIMC = bnf_increase_LIMC(LIMC,LIMCMAX);
1037
if (DEBUGLEVEL>2 && LIMC > LIMC0)
1038
err_printf("%s*** Bach constant: %f\n", FIRST?"":"\n", LIMC/LOGD2);
1039
FIRST = 0; set_avma(av);
1040
guncloneNULL(BQ.subFB);
1041
guncloneNULL(BQ.powsubFB);
1042
clearhash(BQ.hashtab);
1043
if (LIMC < cp) LIMC = cp;
1044
if (LIMC2 < LIMC) LIMC2 = LIMC;
1045
if (BQ.PRECREG) qfr_data_init(q.D, BQ.PRECREG, &q);
1046
1047
FBquad(&BQ, LIMC2, LIMC, &GRHcheck);
1048
if (DEBUGLEVEL>2) timer_printf(&T, "factor base");
1049
BQ.subFB = subFBquad(&BQ, q.D, lim + 0.5, minSFB);
1050
if (DEBUGLEVEL>2) timer_printf(&T, "subFBquad = %Ps",
1051
vecpermute(BQ.FB, BQ.subFB));
1052
nsubFB = lg(BQ.subFB) - 1;
1053
}
1054
while (nsubFB < (expi(D) > 15 ? 3 : 2));
1055
/* invhr = 2^r1 (2pi)^r2 / sqrt(D) w ~ L(chi,1) / hR */
1056
invhr = gmul(dbltor((BQ.PRECREG?2.:M_PI)/sdrc),
1057
compute_invresquad(&GRHcheck, LIMC));
1058
BQ.powsubFB = powsubFBquad(&BQ,CBUCH+1);
1059
if (DEBUGLEVEL>2) timer_printf(&T, "powsubFBquad");
1060
BQ.limhash = (LIMC & HIGHMASK)? (HIGHBIT>>1): LIMC*LIMC;
1061
1062
need = BQ.KC + RELSUP - 2;
1063
current = 0;
1064
W = NULL;
1065
sfb_trials = nreldep = nrelsup = 0;
1066
s = nsubFB + RELSUP;
1067
av2 = avma;
1068
1069
do
1070
{
1071
if ((nreldep & 3) == 1 || (nrelsup & 7) == 1) {
1072
if (DEBUGLEVEL>2) err_printf("*** Changing sub factor base\n");
1073
gunclone(BQ.subFB);
1074
gunclone(BQ.powsubFB);
1075
BQ.subFB = gclone(vecslice(BQ.vperm, 1, nsubFB));
1076
BQ.powsubFB = powsubFBquad(&BQ,CBUCH+1);
1077
if (DEBUGLEVEL>2) timer_printf(&T, "powsubFBquad");
1078
clearhash(BQ.hashtab);
1079
}
1080
need += 2;
1081
mat = cgetg(need+1, t_MAT);
1082
extraC = cgetg(need+1, t_VEC);
1083
if (!W) { /* first time */
1084
C = extraC;
1085
triv = trivial_relations(&BQ, mat, C);
1086
if (DEBUGLEVEL>2) err_printf("KC = %ld, need %ld relations\n", BQ.KC, need);
1087
} else {
1088
triv = 0;
1089
if (DEBUGLEVEL>2) err_printf("...need %ld more relations\n", need);
1090
}
1091
if (BQ.PRECREG) {
1092
for (i = triv+1; i<=need; i++) {
1093
gel(mat,i) = zero_zv(BQ.KC);
1094
gel(extraC,i) = cgetr(BQ.PRECREG);
1095
}
1096
real_relations(&BQ, need - triv, &current, s,LIMC,mat + triv,extraC + triv);
1097
} else {
1098
for (i = triv+1; i<=need; i++) {
1099
gel(mat,i) = zero_zv(BQ.KC);
1100
gel(extraC,i) = gen_0;
1101
}
1102
imag_relations(&BQ, need - triv, &current, LIMC,mat + triv);
1103
}
1104
1105
if (!W)
1106
W = hnfspec_i(mat,BQ.vperm,&dep,&B,&C,nsubFB);
1107
else
1108
W = hnfadd_i(W,BQ.vperm,&dep,&B,&C, mat,extraC);
1109
gerepileall(av2, 4, &W,&C,&B,&dep);
1110
need = BQ.KC - (lg(W)-1) - (lg(B)-1);
1111
if (need)
1112
{
1113
if (++nreldep > 15 && cbach < 1) goto START;
1114
continue;
1115
}
1116
1117
h = ZM_det_triangular(W);
1118
if (DEBUGLEVEL>2) err_printf("\n#### Tentative class number: %Ps\n", h);
1119
1120
switch(get_R(&BQ, C, (lg(C)-1) - (lg(B)-1) - (lg(W)-1), mulir(h,invhr), &R))
1121
{
1122
case fupb_PRECI:
1123
BQ.PRECREG = precdbl(BQ.PRECREG);
1124
FIRST = 1; goto START;
1125
1126
case fupb_RELAT:
1127
if (++nrelsup > MAXRELSUP)
1128
{
1129
if (++sfb_trials > SFB_MAX && cbach <= 1) goto START;
1130
if (nsubFB < minss(10,BQ.KC)) nsubFB++;
1131
}
1132
need = minss(BQ.KC, nrelsup);
1133
}
1134
}
1135
while (need);
1136
/* DONE */
1137
if (!quad_be_honest(&BQ)) goto START;
1138
if (DEBUGLEVEL>2) timer_printf(&T, "be honest");
1139
clearhash(BQ.hashtab);
1140
free_GRHcheck(&GRHcheck);
1141
1142
gen = get_clgp(&BQ,W,&cyc);
1143
gunclone(BQ.subFB);
1144
gunclone(BQ.powsubFB);
1145
res = cgetg(5,t_VEC);
1146
gel(res,1) = h;
1147
gel(res,2) = cyc;
1148
gel(res,3) = gen;
1149
gel(res,4) = R; return gerepilecopy(av0,res);
1150
}
1151
1152
GEN
1153
buchimag(GEN D, GEN c, GEN c2, GEN REL)
1154
{ (void)REL; return Buchquad(D,gtodouble(c),gtodouble(c2),0); }
1155
1156
GEN
1157
buchreal(GEN D, GEN flag, GEN c, GEN c2, GEN REL, long prec) {
1158
if (signe(flag)) pari_err_IMPL("narrow class group");
1159
(void)REL; return Buchquad(D,gtodouble(c),gtodouble(c2),prec);
1160
}
1161
1162
GEN
1163
quadclassunit0(GEN x, long flag, GEN data, long prec)
1164
{
1165
long lx;
1166
double c1 = 0.0, c2 = 0.0;
1167
1168
if (!data) lx=1;
1169
else
1170
{
1171
lx = lg(data);
1172
if (typ(data)!=t_VEC) pari_err_TYPE("quadclassunit", data);
1173
if (lx > 7) pari_err_DIM("quadclassunit [tech vector]");
1174
if (lx > 3) lx = 3;
1175
}
1176
switch(lx)
1177
{
1178
case 3: c2 = gtodouble(gel(data,2));
1179
case 2: c1 = gtodouble(gel(data,1));
1180
}
1181
if (flag) pari_err_IMPL("narrow class group");
1182
return Buchquad(x,c1,c2,prec);
1183
}
1184
1185