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
Check the License for details. You should have received a copy of it, along
11
with the package; see the file 'COPYING'. If not, write to the Free Software
12
Foundation, Inc., 51 Franklin Street, Fifth Floor, Boston, MA 02110-1301 USA. */
13
14
/*******************************************************************/
15
/* */
16
/* MAXIMAL ORDERS */
17
/* */
18
/*******************************************************************/
19
#include "pari.h"
20
#include "paripriv.h"
21
22
#define DEBUGLEVEL DEBUGLEVEL_nf
23
24
/* allow p = -1 from factorizations, avoid oo loop on p = 1 */
25
static long
26
safe_Z_pvalrem(GEN x, GEN p, GEN *z)
27
{
28
if (is_pm1(p))
29
{
30
if (signe(p) > 0) return gvaluation(x,p); /*error*/
31
*z = absi(x); return 1;
32
}
33
return Z_pvalrem(x, p, z);
34
}
35
/* D an integer, P a ZV, return a factorization matrix for D over P, removing
36
* entries with 0 exponent. */
37
static GEN
38
fact_from_factors(GEN D, GEN P, long flag)
39
{
40
long i, l = lg(P), iq = 1;
41
GEN Q = cgetg(l+1,t_COL);
42
GEN E = cgetg(l+1,t_COL);
43
for (i=1; i<l; i++)
44
{
45
GEN p = gel(P,i);
46
long k;
47
if (flag && !equalim1(p))
48
{
49
p = gcdii(p, D);
50
if (is_pm1(p)) continue;
51
}
52
k = safe_Z_pvalrem(D, p, &D);
53
if (k) { gel(Q,iq) = p; gel(E,iq) = utoipos(k); iq++; }
54
}
55
D = absi_shallow(D);
56
if (!equali1(D))
57
{
58
long k = Z_isanypower(D, &D);
59
if (!k) k = 1;
60
gel(Q,iq) = D; gel(E,iq) = utoipos(k); iq++;
61
}
62
setlg(Q,iq);
63
setlg(E,iq); return mkmat2(Q,E);
64
}
65
66
/* d a t_INT; f a t_MAT factorisation of some t_INT sharing some divisors
67
* with d, or a prime (t_INT). Return a factorization F of d: "primes"
68
* entries in f _may_ be composite, and are included as is in d. */
69
static GEN
70
update_fact(GEN d, GEN f)
71
{
72
GEN P;
73
switch (typ(f))
74
{
75
case t_INT: case t_VEC: case t_COL: return f;
76
case t_MAT:
77
if (lg(f) == 3) { P = gel(f,1); break; }
78
/*fall through*/
79
default:
80
pari_err_TYPE("nfbasis [factorization expected]",f);
81
return NULL;/*LCOV_EXCL_LINE*/
82
}
83
return fact_from_factors(d, P, 1);
84
}
85
86
/* T = C T0(X/L); C = L^d / lt(T0), d = deg(T)
87
* disc T = C^2(d - 1) L^-(d(d-1)) disc T0 = (L^d / lt(T0)^2)^(d-1) disc T0 */
88
static GEN
89
set_disc(nfmaxord_t *S)
90
{
91
GEN L, dT;
92
long d;
93
if (S->T0 == S->T) return ZX_disc(S->T);
94
d = degpol(S->T0);
95
L = S->unscale;
96
if (typ(L) == t_FRAC && abscmpii(gel(L,1), gel(L,2)) < 0)
97
dT = ZX_disc(S->T); /* more efficient */
98
else
99
{
100
GEN l0 = leading_coeff(S->T0);
101
GEN a = gpowgs(gdiv(gpowgs(L, d), sqri(l0)), d-1);
102
dT = gmul(a, ZX_disc(S->T0)); /* more efficient */
103
}
104
return S->dT = dT;
105
}
106
107
/* dT != 0 */
108
static GEN
109
poldiscfactors_i(GEN T, GEN dT, long flag)
110
{
111
GEN U, fa, Z, E, P, Tp;
112
long i, l;
113
114
fa = absZ_factor_limit_strict(dT, minuu(tridiv_bound(dT), maxprime()), &U);
115
if (!U) return fa;
116
Z = mkcol(gel(U,1)); P = gel(fa,1); Tp = NULL;
117
while (lg(Z) != 1)
118
{ /* pop and handle last element of Z */
119
GEN p = gel(Z, lg(Z)-1), r;
120
setlg(Z, lg(Z)-1);
121
if (!Tp) /* first time: p is composite and not a power */
122
Tp = ZX_deriv(T);
123
else
124
{
125
(void)Z_isanypower(p, &p);
126
if ((flag || lgefint(p)==3) && BPSW_psp(p))
127
{ P = vec_append(P, p); continue; }
128
}
129
r = FpX_gcd_check(T, Tp, p);
130
if (r)
131
Z = shallowconcat(Z, Z_cba(r, diviiexact(p,r)));
132
else if (flag)
133
P = shallowconcat(P, gel(Z_factor(p),1));
134
else
135
P = vec_append(P, p);
136
}
137
ZV_sort_inplace(P); l = lg(P); E = cgetg(l, t_COL);
138
for (i = 1; i < l; i++) gel(E,i) = utoipos(Z_pvalrem(dT, gel(P,i), &dT));
139
return mkmat2(P,E);
140
}
141
142
GEN
143
poldiscfactors(GEN T, long flag)
144
{
145
pari_sp av = avma;
146
GEN dT;
147
if (typ(T) != t_POL || !RgX_is_ZX(T)) pari_err_TYPE("poldiscfactors",T);
148
if (flag < 0 || flag > 1) pari_err_FLAG("poldiscfactors");
149
dT = ZX_disc(T);
150
if (!signe(dT)) retmkvec2(gen_0, Z_factor(gen_0));
151
return gerepilecopy(av, mkvec2(dT, poldiscfactors_i(T, dT, flag)));
152
}
153
154
static void
155
nfmaxord_check_args(nfmaxord_t *S, GEN T, long flag)
156
{
157
GEN dT, L, E, P, fa = NULL;
158
pari_timer t;
159
long l, ty = typ(T);
160
161
if (DEBUGLEVEL) timer_start(&t);
162
if (ty == t_VEC) {
163
if (lg(T) != 3) pari_err_TYPE("nfmaxord",T);
164
fa = gel(T,2); T = gel(T,1); ty = typ(T);
165
}
166
if (ty != t_POL) pari_err_TYPE("nfmaxord",T);
167
T = Q_primpart(T);
168
if (degpol(T) <= 0) pari_err_CONSTPOL("nfmaxord");
169
RgX_check_ZX(T, "nfmaxord");
170
S->T0 = T;
171
S->T = T = ZX_Q_normalize(T, &L);
172
S->unscale = L;
173
S->dT = dT = set_disc(S);
174
if (!signe(dT)) pari_err_IRREDPOL("nfmaxord",T);
175
if (fa)
176
{
177
const long MIN = 100; /* include at least all p < 101 */
178
GEN P0 = NULL, U;
179
if (!isint1(L)) fa = update_fact(dT, fa);
180
switch(typ(fa))
181
{
182
case t_MAT:
183
if (!is_Z_factornon0(fa)) pari_err_TYPE("nfmaxord",fa);
184
P0 = gel(fa,1); /* fall through */
185
case t_VEC: case t_COL:
186
if (!P0)
187
{
188
if (!RgV_is_ZV(fa)) pari_err_TYPE("nfmaxord",fa);
189
P0 = fa;
190
}
191
P = gel(absZ_factor_limit_strict(dT, MIN, &U), 1);
192
if (lg(P) != 0) { settyp(P, typ(P0)); P0 = shallowconcat(P0,P); }
193
P0 = ZV_sort_uniq(P0);
194
fa = fact_from_factors(dT, P0, 0);
195
break;
196
case t_INT:
197
fa = absZ_factor_limit(dT, (signe(fa) <= 0)? 1: maxuu(itou(fa), MIN));
198
break;
199
default:
200
pari_err_TYPE("nfmaxord",fa);
201
}
202
}
203
else
204
fa = poldiscfactors_i(T, dT, !(flag & nf_PARTIALFACT));
205
P = gel(fa,1); l = lg(P);
206
E = gel(fa,2);
207
if (l > 1 && is_pm1(gel(P,1)))
208
{
209
l--;
210
P = vecslice(P, 2, l);
211
E = vecslice(E, 2, l);
212
}
213
S->dTP = P;
214
S->dTE = vec_to_vecsmall(E);
215
if (DEBUGLEVEL>2) timer_printf(&t, "disc. factorisation");
216
}
217
218
static int
219
fnz(GEN x,long j)
220
{
221
long i;
222
for (i=1; i<j; i++)
223
if (signe(gel(x,i))) return 0;
224
return 1;
225
}
226
/* return list u[i], 2 by 2 coprime with the same prime divisors as ab */
227
static GEN
228
get_coprimes(GEN a, GEN b)
229
{
230
long i, k = 1;
231
GEN u = cgetg(3, t_COL);
232
gel(u,1) = a;
233
gel(u,2) = b;
234
/* u1,..., uk 2 by 2 coprime */
235
while (k+1 < lg(u))
236
{
237
GEN d, c = gel(u,k+1);
238
if (is_pm1(c)) { k++; continue; }
239
for (i=1; i<=k; i++)
240
{
241
GEN ui = gel(u,i);
242
if (is_pm1(ui)) continue;
243
d = gcdii(c, ui);
244
if (d == gen_1) continue;
245
c = diviiexact(c, d);
246
gel(u,i) = diviiexact(ui, d);
247
u = vec_append(u, d);
248
}
249
gel(u,++k) = c;
250
}
251
for (i = k = 1; i < lg(u); i++)
252
if (!is_pm1(gel(u,i))) gel(u,k++) = gel(u,i);
253
setlg(u, k); return u;
254
}
255
256
/*******************************************************************/
257
/* */
258
/* ROUND 4 */
259
/* */
260
/*******************************************************************/
261
typedef struct {
262
/* constants */
263
long pisprime; /* -1: unknown, 1: prime, 0: composite */
264
GEN p, f; /* goal: factor f p-adically */
265
long df;
266
GEN pdf; /* p^df = reduced discriminant of f */
267
long mf; /* */
268
GEN psf, pmf; /* stability precision for f, wanted precision for f */
269
long vpsf; /* v_p(p_f) */
270
/* these are updated along the way */
271
GEN phi; /* a p-integer, in Q[X] */
272
GEN phi0; /* a p-integer, in Q[X] from testb2 / testc2, to be composed with
273
* phi when correct precision is known */
274
GEN chi; /* characteristic polynomial of phi (mod psc) in Z[X] */
275
GEN nu; /* irreducible divisor of chi mod p, in Z[X] */
276
GEN invnu; /* numerator ( 1/ Mod(nu, chi) mod pmr ) */
277
GEN Dinvnu;/* denominator ( ... ) */
278
long vDinvnu; /* v_p(Dinvnu) */
279
GEN prc, psc; /* reduced discriminant of chi, stability precision for chi */
280
long vpsc; /* v_p(p_c) */
281
GEN ns, nsf, precns; /* cached Newton sums for nsf and their precision */
282
} decomp_t;
283
static GEN maxord_i(decomp_t *S, GEN p, GEN f, long mf, GEN w, long flag);
284
static GEN dbasis(GEN p, GEN f, long mf, GEN alpha, GEN U);
285
static GEN maxord(GEN p,GEN f,long mf);
286
static GEN ZX_Dedekind(GEN F, GEN *pg, GEN p);
287
288
/* Warning: data computed for T = ZX_Q_normalize(T0). If S.unscale !=
289
* gen_1, caller must take steps to correct the components if it wishes
290
* to stick to the original T0. Return a vector of p-maximal orders, for
291
* those p s.t p^2 | disc(T) [ = S->dTP ]*/
292
static GEN
293
get_maxord(nfmaxord_t *S, GEN T0, long flag)
294
{
295
VOLATILE GEN P, E, O;
296
VOLATILE long lP, i, k;
297
298
nfmaxord_check_args(S, T0, flag);
299
P = S->dTP; lP = lg(P);
300
E = S->dTE;
301
O = cgetg(1, t_VEC);
302
for (i=1; i<lP; i++)
303
{
304
VOLATILE pari_sp av;
305
/* includes the silly case where P[i] = -1 */
306
if (E[i] <= 1) { O = vec_append(O, gen_1); continue; }
307
av = avma;
308
pari_CATCH(CATCH_ALL) {
309
GEN N, u, err = pari_err_last();
310
long l;
311
switch(err_get_num(err))
312
{
313
case e_INV:
314
{
315
GEN p, x = err_get_compo(err, 2);
316
if (typ(x) == t_INTMOD)
317
{ /* caught false prime, update factorization */
318
p = gcdii(gel(x,1), gel(x,2));
319
u = diviiexact(gel(x,1),p);
320
if (DEBUGLEVEL) pari_warn(warner,"impossible inverse: %Ps", x);
321
gerepileall(av, 2, &p, &u);
322
323
u = get_coprimes(p, u); l = lg(u);
324
/* no small factors, but often a prime power */
325
for (k = 1; k < l; k++) (void)Z_isanypower(gel(u,k), &gel(u,k));
326
break;
327
}
328
/* fall through */
329
}
330
case e_PRIME: case e_IRREDPOL:
331
{ /* we're here because we failed BPSW_isprime(), no point in
332
* reporting a possible counter-example to the BPSW test */
333
GEN p = gel(P,i);
334
set_avma(av);
335
if (DEBUGLEVEL)
336
pari_warn(warner,"large composite in nfmaxord:loop(), %Ps", p);
337
if (expi(p) < 100) /* factor should require ~20ms for this */
338
u = gel(Z_factor(p), 1);
339
else
340
{ /* give up, probably not maximal */
341
GEN B, g, k = ZX_Dedekind(S->T, &g, p);
342
k = FpX_normalize(k, p);
343
B = dbasis(p, S->T, E[i], NULL, FpX_div(S->T,k,p));
344
O = vec_append(O, B);
345
pari_CATCH_reset(); continue;
346
}
347
break;
348
}
349
default: pari_err(0, err);
350
return NULL;/*LCOV_EXCL_LINE*/
351
}
352
l = lg(u);
353
gel(P,i) = gel(u,1);
354
P = shallowconcat(P, vecslice(u, 2, l-1));
355
av = avma;
356
N = S->dT; E[i] = Z_pvalrem(N, gel(P,i), &N);
357
for (k=lP, lP=lg(P); k < lP; k++) E[k] = Z_pvalrem(N, gel(P,k), &N);
358
} pari_RETRY {
359
if (DEBUGLEVEL>2) err_printf("Treating p^k = %Ps^%ld\n",P[i],E[i]);
360
O = vec_append(O, maxord(gel(P,i),S->T,E[i]));
361
} pari_ENDCATCH;
362
}
363
S->dTP = P; return O;
364
}
365
366
/* M a QM, return denominator of diagonal. All denominators are powers of
367
* a given integer */
368
static GEN
369
diag_denom(GEN M)
370
{
371
GEN d = gen_1;
372
long j, l = lg(M);
373
for (j=1; j<l; j++)
374
{
375
GEN t = gcoeff(M,j,j);
376
if (typ(t) == t_INT) continue;
377
t = gel(t,2);
378
if (abscmpii(t,d) > 0) d = t;
379
}
380
return d;
381
}
382
static void
383
setPE(GEN D, GEN P, GEN *pP, GEN *pE)
384
{
385
long k, j, l = lg(P);
386
GEN P2, E2;
387
*pP = P2 = cgetg(l, t_COL);
388
*pE = E2 = cgetg(l, t_VECSMALL);
389
for (k = j = 1; j < l; j++)
390
{
391
long v = Z_pvalrem(D, gel(P,j), &D);
392
if (v) { gel(P2,k) = gel(P,j); E2[k] = v; k++; }
393
}
394
setlg(P2, k);
395
setlg(E2, k);
396
}
397
void
398
nfmaxord(nfmaxord_t *S, GEN T0, long flag)
399
{
400
GEN O = get_maxord(S, T0, flag);
401
GEN f = S->T, P = S->dTP, a = NULL, da = NULL;
402
long n = degpol(f), lP = lg(P), i, j, k;
403
int centered = 0;
404
pari_sp av = avma;
405
/* r1 & basden not initialized here */
406
S->r1 = -1;
407
S->basden = NULL;
408
for (i=1; i<lP; i++)
409
{
410
GEN M, db, b = gel(O,i);
411
if (b == gen_1) continue;
412
db = diag_denom(b);
413
if (db == gen_1) continue;
414
415
/* db = denom(b), (da,db) = 1. Compute da Im(b) + db Im(a) */
416
b = Q_muli_to_int(b,db);
417
if (!da) { da = db; a = b; }
418
else
419
{ /* optimization: easy as long as both matrix are diagonal */
420
j=2; while (j<=n && fnz(gel(a,j),j) && fnz(gel(b,j),j)) j++;
421
k = j-1; M = cgetg(2*n-k+1,t_MAT);
422
for (j=1; j<=k; j++)
423
{
424
gel(M,j) = gel(a,j);
425
gcoeff(M,j,j) = mulii(gcoeff(a,j,j),gcoeff(b,j,j));
426
}
427
/* could reduce mod M(j,j) but not worth it: usually close to da*db */
428
for ( ; j<=n; j++) gel(M,j) = ZC_Z_mul(gel(a,j), db);
429
for ( ; j<=2*n-k; j++) gel(M,j) = ZC_Z_mul(gel(b,j+k-n), da);
430
da = mulii(da,db);
431
a = ZM_hnfmodall_i(M, da, hnf_MODID|hnf_CENTER);
432
gerepileall(av, 2, &a, &da);
433
centered = 1;
434
}
435
}
436
if (da)
437
{
438
GEN index = diviiexact(da, gcoeff(a,1,1));
439
for (j=2; j<=n; j++) index = mulii(index, diviiexact(da, gcoeff(a,j,j)));
440
if (!centered) a = ZM_hnfcenter(a);
441
a = RgM_Rg_div(a, da);
442
S->index = index;
443
S->dK = diviiexact(S->dT, sqri(index));
444
}
445
else
446
{
447
S->index = gen_1;
448
S->dK = S->dT;
449
a = matid(n);
450
}
451
setPE(S->dK, P, &S->dKP, &S->dKE);
452
S->basis = RgM_to_RgXV(a, varn(f));
453
}
454
GEN
455
nfbasis(GEN x, GEN *pdK)
456
{
457
pari_sp av = avma;
458
nfmaxord_t S;
459
GEN B;
460
nfmaxord(&S, x, 0);
461
B = RgXV_unscale(S.basis, S.unscale);
462
if (pdK) *pdK = S.dK;
463
gerepileall(av, pdK? 2: 1, &B, pdK); return B;
464
}
465
/* field discriminant: faster than nfmaxord, use local data only */
466
static GEN
467
maxord_disc(nfmaxord_t *S, GEN x)
468
{
469
GEN O = get_maxord(S, x, 0), I = gen_1;
470
long n = degpol(S->T), lP = lg(O), i, j;
471
for (i = 1; i < lP; i++)
472
{
473
GEN b = gel(O,i);
474
if (b == gen_1) continue;
475
for (j = 1; j <= n; j++)
476
{
477
GEN c = gcoeff(b,j,j);
478
if (typ(c) == t_FRAC) I = mulii(I, gel(c,2)) ;
479
}
480
}
481
return diviiexact(S->dT, sqri(I));
482
}
483
GEN
484
nfdisc(GEN x)
485
{
486
pari_sp av = avma;
487
nfmaxord_t S;
488
return gerepileuptoint(av, maxord_disc(&S, x));
489
}
490
GEN
491
nfdiscfactors(GEN x)
492
{
493
pari_sp av = avma;
494
GEN E, P, D, nf = checknf_i(x);
495
if (nf)
496
{
497
D = nf_get_disc(nf);
498
P = nf_get_ramified_primes(nf);
499
}
500
else
501
{
502
nfmaxord_t S;
503
D = maxord_disc(&S, x);
504
P = S.dTP;
505
}
506
setPE(D, P, &P, &E); settyp(P, t_COL);
507
return gerepilecopy(av, mkvec2(D, mkmat2(P, zc_to_ZC(E))));
508
}
509
510
static ulong
511
Flx_checkdeflate(GEN x)
512
{
513
ulong d = 0, i, lx = (ulong)lg(x);
514
for (i=3; i<lx; i++)
515
if (x[i]) { d = ugcd(d,i-2); if (d == 1) break; }
516
return d;
517
}
518
519
/* product of (monic) irreducible factors of f over Fp[X]
520
* Assume f reduced mod p, otherwise valuation at x may be wrong */
521
static GEN
522
Flx_radical(GEN f, ulong p)
523
{
524
long v0 = Flx_valrem(f, &f);
525
ulong du, d, e;
526
GEN u;
527
528
d = Flx_checkdeflate(f);
529
if (!d) return v0? polx_Flx(f[1]): pol1_Flx(f[1]);
530
if (u_lvalrem(d,p, &e)) f = Flx_deflate(f, d/e); /* f(x^p^i) -> f(x) */
531
u = Flx_gcd(f, Flx_deriv(f, p), p); /* (f,f') */
532
du = degpol(u);
533
if (du)
534
{
535
if (du == (ulong)degpol(f))
536
f = Flx_radical(Flx_deflate(f,p), p);
537
else
538
{
539
u = Flx_normalize(u, p);
540
f = Flx_div(f, u, p);
541
if (p <= du)
542
{
543
GEN w = (degpol(f) >= degpol(u))? Flx_rem(f, u, p): f;
544
w = Flxq_powu(w, du, u, p);
545
w = Flx_div(u, Flx_gcd(w,u,p), p); /* u / gcd(u, v^(deg u-1)) */
546
f = Flx_mul(f, Flx_radical(Flx_deflate(w,p), p), p);
547
}
548
}
549
}
550
if (v0) f = Flx_shift(f, 1);
551
return f;
552
}
553
/* Assume f reduced mod p, otherwise valuation at x may be wrong */
554
static GEN
555
FpX_radical(GEN f, GEN p)
556
{
557
GEN u;
558
long v0;
559
if (lgefint(p) == 3)
560
{
561
ulong q = p[2];
562
return Flx_to_ZX( Flx_radical(ZX_to_Flx(f, q), q) );
563
}
564
v0 = ZX_valrem(f, &f);
565
u = FpX_gcd(f,FpX_deriv(f, p), p);
566
if (degpol(u)) f = FpX_div(f, u, p);
567
if (v0) f = RgX_shift(f, 1);
568
return f;
569
}
570
/* f / a */
571
static GEN
572
zx_z_div(GEN f, ulong a)
573
{
574
long i, l = lg(f);
575
GEN g = cgetg(l, t_VECSMALL);
576
g[1] = f[1];
577
for (i = 2; i < l; i++) g[i] = f[i] / a;
578
return g;
579
}
580
/* Dedekind criterion; return k = gcd(g,h, (f-gh)/p), where
581
* f = \prod f_i^e_i, g = \prod f_i, h = \prod f_i^{e_i-1}
582
* k = 1 iff Z[X]/(f) is p-maximal */
583
static GEN
584
ZX_Dedekind(GEN F, GEN *pg, GEN p)
585
{
586
GEN k, h, g, f, f2;
587
ulong q = p[2];
588
if (lgefint(p) == 3 && q < (1UL << BITS_IN_HALFULONG))
589
{
590
ulong q2 = q*q;
591
f2 = ZX_to_Flx(F, q2);
592
f = Flx_red(f2, q);
593
g = Flx_radical(f, q);
594
h = Flx_div(f, g, q);
595
k = zx_z_div(Flx_sub(f2, Flx_mul(g,h,q2), q2), q);
596
k = Flx_gcd(k, Flx_gcd(g,h,q), q);
597
k = Flx_to_ZX(k);
598
g = Flx_to_ZX(g);
599
}
600
else
601
{
602
f2 = FpX_red(F, sqri(p));
603
f = FpX_red(f2, p);
604
g = FpX_radical(f, p);
605
h = FpX_div(f, g, p);
606
k = ZX_Z_divexact(ZX_sub(f2, ZX_mul(g,h)), p);
607
k = FpX_gcd(FpX_red(k, p), FpX_gcd(g,h,p), p);
608
}
609
*pg = g; return k;
610
}
611
612
/* p-maximal order of Z[x]/f; mf = v_p(Disc(f)) or < 0 [unknown].
613
* Return gen_1 if p-maximal */
614
static GEN
615
maxord(GEN p, GEN f, long mf)
616
{
617
const pari_sp av = avma;
618
GEN res, g, k = ZX_Dedekind(f, &g, p);
619
long dk = degpol(k);
620
if (DEBUGLEVEL>2) err_printf(" ZX_Dedekind: gcd has degree %ld\n", dk);
621
if (!dk) { set_avma(av); return gen_1; }
622
if (mf < 0) mf = ZpX_disc_val(f, p);
623
k = FpX_normalize(k, p);
624
if (2*dk >= mf-1)
625
res = dbasis(p, f, mf, NULL, FpX_div(f,k,p));
626
else
627
{
628
GEN w, F1, F2;
629
decomp_t S;
630
F1 = FpX_factor(k,p);
631
F2 = FpX_factor(FpX_div(g,k,p),p);
632
w = merge_sort_uniq(gel(F1,1),gel(F2,1),(void*)cmpii,&gen_cmp_RgX);
633
res = maxord_i(&S, p, f, mf, w, 0);
634
}
635
return gerepilecopy(av,res);
636
}
637
/* T monic separable ZX, p prime */
638
GEN
639
ZpX_primedec(GEN T, GEN p)
640
{
641
const pari_sp av = avma;
642
GEN w, F1, F2, res, g, k = ZX_Dedekind(T, &g, p);
643
decomp_t S;
644
if (!degpol(k)) return zm_to_ZM(FpX_degfact(T, p));
645
k = FpX_normalize(k, p);
646
F1 = FpX_factor(k,p);
647
F2 = FpX_factor(FpX_div(g,k,p),p);
648
w = merge_sort_uniq(gel(F1,1),gel(F2,1),(void*)cmpii,&gen_cmp_RgX);
649
res = maxord_i(&S, p, T, ZpX_disc_val(T, p), w, -1);
650
if (!res)
651
{
652
long f = degpol(S.nu), e = degpol(T) / f;
653
set_avma(av); retmkmat2(mkcols(f), mkcols(e));
654
}
655
return gerepilecopy(av,res);
656
}
657
658
static GEN
659
Zlx_sylvester_echelon(GEN f1, GEN f2, long early_abort, ulong p, ulong pm)
660
{
661
long j, n = degpol(f1);
662
GEN h, a = cgetg(n+1,t_MAT);
663
f1 = Flx_get_red(f1, pm);
664
h = Flx_rem(f2,f1,pm);
665
for (j=1;; j++)
666
{
667
gel(a,j) = Flx_to_Flv(h, n);
668
if (j == n) break;
669
h = Flx_rem(Flx_shift(h, 1), f1, pm);
670
}
671
return zlm_echelon(a, early_abort, p, pm);
672
}
673
/* Sylvester's matrix, mod p^m (assumes f1 monic). If early_abort
674
* is set, return NULL if one pivot is 0 mod p^m */
675
static GEN
676
ZpX_sylvester_echelon(GEN f1, GEN f2, long early_abort, GEN p, GEN pm)
677
{
678
long j, n = degpol(f1);
679
GEN h, a = cgetg(n+1,t_MAT);
680
h = FpXQ_red(f2,f1,pm);
681
for (j=1;; j++)
682
{
683
gel(a,j) = RgX_to_RgC(h, n);
684
if (j == n) break;
685
h = FpX_rem(RgX_shift_shallow(h, 1), f1, pm);
686
}
687
return ZpM_echelon(a, early_abort, p, pm);
688
}
689
690
/* polynomial gcd mod p^m (assumes f1 monic). Return a QpX ! */
691
static GEN
692
Zlx_gcd(GEN f1, GEN f2, ulong p, ulong pm)
693
{
694
pari_sp av = avma;
695
GEN a = Zlx_sylvester_echelon(f1,f2,0,p,pm);
696
long c, l = lg(a), sv = f1[1];
697
for (c = 1; c < l; c++)
698
{
699
ulong t = ucoeff(a,c,c);
700
if (t)
701
{
702
a = Flx_to_ZX(Flv_to_Flx(gel(a,c), sv));
703
if (t == 1) return gerepilecopy(av, a);
704
return gerepileupto(av, RgX_Rg_div(a, utoipos(t)));
705
}
706
}
707
set_avma(av);
708
a = cgetg(2,t_POL); a[1] = sv; return a;
709
}
710
GEN
711
ZpX_gcd(GEN f1, GEN f2, GEN p, GEN pm)
712
{
713
pari_sp av = avma;
714
GEN a;
715
long c, l, v;
716
if (lgefint(pm) == 3)
717
{
718
ulong q = pm[2];
719
return Zlx_gcd(ZX_to_Flx(f1, q), ZX_to_Flx(f2,q), p[2], q);
720
}
721
a = ZpX_sylvester_echelon(f1,f2,0,p,pm);
722
l = lg(a); v = varn(f1);
723
for (c = 1; c < l; c++)
724
{
725
GEN t = gcoeff(a,c,c);
726
if (signe(t))
727
{
728
a = RgV_to_RgX(gel(a,c), v);
729
if (equali1(t)) return gerepilecopy(av, a);
730
return gerepileupto(av, RgX_Rg_div(a, t));
731
}
732
}
733
set_avma(av); return pol_0(v);
734
}
735
736
/* Return m > 0, such that p^m ~ 2^16 for initial value of m; p > 1 */
737
static long
738
init_m(GEN p)
739
{
740
if (lgefint(p) > 3) return 1;
741
return (long)(16 / log2(p[2]));
742
}
743
744
/* reduced resultant mod p^m (assumes x monic) */
745
GEN
746
ZpX_reduced_resultant(GEN x, GEN y, GEN p, GEN pm)
747
{
748
pari_sp av = avma;
749
GEN z;
750
if (lgefint(pm) == 3)
751
{
752
ulong q = pm[2];
753
z = Zlx_sylvester_echelon(ZX_to_Flx(x,q), ZX_to_Flx(y,q),0,p[2],q);
754
if (lg(z) > 1)
755
{
756
ulong c = ucoeff(z,1,1);
757
if (c) { set_avma(av); return utoipos(c); }
758
}
759
}
760
else
761
{
762
z = ZpX_sylvester_echelon(x,y,0,p,pm);
763
if (lg(z) > 1)
764
{
765
GEN c = gcoeff(z,1,1);
766
if (signe(c)) return gerepileuptoint(av, c);
767
}
768
}
769
set_avma(av); return gen_0;
770
}
771
/* Assume Res(f,g) divides p^M. Return Res(f, g), using dynamic p-adic
772
* precision (until result is nonzero or p^M). */
773
GEN
774
ZpX_reduced_resultant_fast(GEN f, GEN g, GEN p, long M)
775
{
776
GEN R, q = NULL;
777
long m;
778
m = init_m(p); if (m < 1) m = 1;
779
for(;; m <<= 1) {
780
if (M < 2*m) break;
781
q = q? sqri(q): powiu(p, m); /* p^m */
782
R = ZpX_reduced_resultant(f,g, p, q); if (signe(R)) return R;
783
}
784
q = powiu(p, M);
785
R = ZpX_reduced_resultant(f,g, p, q); return signe(R)? R: q;
786
}
787
788
/* v_p(Res(x,y) mod p^m), assumes (lc(x),p) = 1 */
789
static long
790
ZpX_resultant_val_i(GEN x, GEN y, GEN p, GEN pm)
791
{
792
pari_sp av = avma;
793
GEN z;
794
long i, l, v;
795
if (lgefint(pm) == 3)
796
{
797
ulong q = pm[2], pp = p[2];
798
z = Zlx_sylvester_echelon(ZX_to_Flx(x,q), ZX_to_Flx(y,q), 1, pp, q);
799
if (!z) return gc_long(av,-1); /* failure */
800
v = 0; l = lg(z);
801
for (i = 1; i < l; i++) v += u_lval(ucoeff(z,i,i), pp);
802
}
803
else
804
{
805
z = ZpX_sylvester_echelon(x, y, 1, p, pm);
806
if (!z) return gc_long(av,-1); /* failure */
807
v = 0; l = lg(z);
808
for (i = 1; i < l; i++) v += Z_pval(gcoeff(z,i,i), p);
809
}
810
return v;
811
}
812
813
/* assume (lc(f),p) = 1; no assumption on g */
814
long
815
ZpX_resultant_val(GEN f, GEN g, GEN p, long M)
816
{
817
pari_sp av = avma;
818
GEN q = NULL;
819
long v, m;
820
m = init_m(p); if (m < 2) m = 2;
821
for(;; m <<= 1) {
822
if (m > M) m = M;
823
q = q? sqri(q): powiu(p, m); /* p^m */
824
v = ZpX_resultant_val_i(f,g, p, q); if (v >= 0) return gc_long(av,v);
825
if (m == M) return gc_long(av,M);
826
}
827
}
828
829
/* assume f separable and (lc(f),p) = 1 */
830
long
831
ZpX_disc_val(GEN f, GEN p)
832
{
833
pari_sp av = avma;
834
long v;
835
if (degpol(f) == 1) return 0;
836
v = ZpX_resultant_val(f, ZX_deriv(f), p, LONG_MAX);
837
return gc_long(av,v);
838
}
839
840
/* *e a ZX, *d, *z in Z, *d = p^(*vd). Simplify e / d by cancelling a
841
* common factor p^v; if z!=NULL, update it by cancelling the same power of p */
842
static void
843
update_den(GEN p, GEN *e, GEN *d, long *vd, GEN *z)
844
{
845
GEN newe;
846
long ve = ZX_pvalrem(*e, p, &newe);
847
if (ve) {
848
GEN newd;
849
long v = minss(*vd, ve);
850
if (v) {
851
if (v == *vd)
852
{ /* rare, denominator cancelled */
853
if (ve != v) newe = ZX_Z_mul(newe, powiu(p, ve - v));
854
newd = gen_1;
855
*vd = 0;
856
if (z) *z =diviiexact(*z, powiu(p, v));
857
}
858
else
859
{ /* v = ve < vd, generic case */
860
GEN q = powiu(p, v);
861
newd = diviiexact(*d, q);
862
*vd -= v;
863
if (z) *z = diviiexact(*z, q);
864
}
865
*e = newe;
866
*d = newd;
867
}
868
}
869
}
870
871
/* return denominator, a power of p */
872
static GEN
873
QpX_denom(GEN x)
874
{
875
long i, l = lg(x);
876
GEN maxd = gen_1;
877
for (i=2; i<l; i++)
878
{
879
GEN d = gel(x,i);
880
if (typ(d) == t_FRAC && cmpii(gel(d,2), maxd) > 0) maxd = gel(d,2);
881
}
882
return maxd;
883
}
884
static GEN
885
QpXV_denom(GEN x)
886
{
887
long l = lg(x), i;
888
GEN maxd = gen_1;
889
for (i = 1; i < l; i++)
890
{
891
GEN d = QpX_denom(gel(x,i));
892
if (cmpii(d, maxd) > 0) maxd = d;
893
}
894
return maxd;
895
}
896
897
static GEN
898
QpX_remove_denom(GEN x, GEN p, GEN *pdx, long *pv)
899
{
900
*pdx = QpX_denom(x);
901
if (*pdx == gen_1) { *pv = 0; *pdx = NULL; }
902
else {
903
x = Q_muli_to_int(x,*pdx);
904
*pv = Z_pval(*pdx, p);
905
}
906
return x;
907
}
908
909
/* p^v * f o g mod (T,q). q = p^vq */
910
static GEN
911
compmod(GEN p, GEN f, GEN g, GEN T, GEN q, long v)
912
{
913
GEN D = NULL, z, df, dg, qD;
914
long vD = 0, vdf, vdg;
915
916
f = QpX_remove_denom(f, p, &df, &vdf);
917
if (typ(g) == t_VEC) /* [num,den,v_p(den)] */
918
{ vdg = itos(gel(g,3)); dg = gel(g,2); g = gel(g,1); }
919
else
920
g = QpX_remove_denom(g, p, &dg, &vdg);
921
if (df) { D = df; vD = vdf; }
922
if (dg) {
923
long degf = degpol(f);
924
D = mul_content(D, powiu(dg, degf));
925
vD += degf * vdg;
926
}
927
qD = D ? mulii(q, D): q;
928
if (dg) f = FpX_rescale(f, dg, qD);
929
z = FpX_FpXQ_eval(f, g, T, qD);
930
if (!D) {
931
if (v) {
932
if (v > 0)
933
z = ZX_Z_mul(z, powiu(p, v));
934
else
935
z = RgX_Rg_div(z, powiu(p, -v));
936
}
937
return z;
938
}
939
update_den(p, &z, &D, &vD, NULL);
940
qD = mulii(D,q);
941
if (v) vD -= v;
942
z = FpX_center_i(z, qD, shifti(qD,-1));
943
if (vD > 0)
944
z = RgX_Rg_div(z, powiu(p, vD));
945
else if (vD < 0)
946
z = ZX_Z_mul(z, powiu(p, -vD));
947
return z;
948
}
949
950
/* fast implementation of ZM_hnfmodid(M, D) / D, D = p^k */
951
static GEN
952
ZpM_hnfmodid(GEN M, GEN p, GEN D)
953
{
954
long i, l = lg(M);
955
M = RgM_Rg_div(ZpM_echelon(M,0,p,D), D);
956
for (i = 1; i < l; i++)
957
if (gequal0(gcoeff(M,i,i))) gcoeff(M,i,i) = gen_1;
958
return M;
959
}
960
961
/* Return Z-basis for Z[a] + U(a)/p Z[a] in Z[t]/(f), mf = v_p(disc f), U
962
* a ZX. Special cases: a = t is coded as NULL, U = 0 is coded as NULL */
963
static GEN
964
dbasis(GEN p, GEN f, long mf, GEN a, GEN U)
965
{
966
long n = degpol(f), i, dU;
967
GEN b, h;
968
969
if (n == 1) return matid(1);
970
if (a && gequalX(a)) a = NULL;
971
if (DEBUGLEVEL>5)
972
{
973
err_printf(" entering Dedekind Basis with parameters p=%Ps\n",p);
974
err_printf(" f = %Ps,\n a = %Ps\n",f, a? a: pol_x(varn(f)));
975
}
976
if (a)
977
{
978
GEN pd = powiu(p, mf >> 1);
979
GEN da, pdp = mulii(pd,p), D = pdp;
980
long vda;
981
dU = U ? degpol(U): 0;
982
b = cgetg(n+1, t_MAT);
983
h = scalarpol(pd, varn(f));
984
a = QpX_remove_denom(a, p, &da, &vda);
985
if (da) D = mulii(D, da);
986
gel(b,1) = scalarcol_shallow(pd, n);
987
for (i=2; i<=n; i++)
988
{
989
if (i == dU+1)
990
h = compmod(p, U, mkvec3(a,da,stoi(vda)), f, pdp, (mf>>1) - 1);
991
else
992
{
993
h = FpXQ_mul(h, a, f, D);
994
if (da) h = ZX_Z_divexact(h, da);
995
}
996
gel(b,i) = RgX_to_RgC(h,n);
997
}
998
return ZpM_hnfmodid(b, p, pd);
999
}
1000
else
1001
{
1002
if (!U) return matid(n);
1003
dU = degpol(U);
1004
if (dU == n) return matid(n);
1005
U = FpX_normalize(U, p);
1006
b = cgetg(n+1, t_MAT);
1007
for (i = 1; i <= dU; i++) gel(b,i) = vec_ei(n, i);
1008
h = RgX_Rg_div(U, p);
1009
for ( ; i <= n; i++)
1010
{
1011
gel(b, i) = RgX_to_RgC(h,n);
1012
if (i == n) break;
1013
h = RgX_shift_shallow(h,1);
1014
}
1015
return b;
1016
}
1017
}
1018
1019
static GEN
1020
get_partial_order_as_pols(GEN p, GEN f)
1021
{
1022
GEN O = maxord(p, f, -1);
1023
long v = varn(f);
1024
return O == gen_1? pol_x_powers(degpol(f), v): RgM_to_RgXV(O, v);
1025
}
1026
1027
static long
1028
p_is_prime(decomp_t *S)
1029
{
1030
if (S->pisprime < 0) S->pisprime = BPSW_psp(S->p);
1031
return S->pisprime;
1032
}
1033
static GEN ZpX_monic_factor_squarefree(GEN f, GEN p, long prec);
1034
1035
/* if flag = 0, maximal order, else factorization to precision r = flag */
1036
static GEN
1037
Decomp(decomp_t *S, long flag)
1038
{
1039
pari_sp av = avma;
1040
GEN fred, pr2, pr, pk, ph2, ph, b1, b2, a, e, de, f1, f2, dt, th, chip;
1041
GEN p = S->p;
1042
long vde, vdt, k, r = maxss(flag, 2*S->df + 1);
1043
1044
if (DEBUGLEVEL>5) err_printf(" entering Decomp: %Ps^%ld\n f = %Ps\n",
1045
p, r, S->f);
1046
else if (DEBUGLEVEL>2) err_printf(" entering Decomp\n");
1047
chip = FpX_red(S->chi, p);
1048
if (!FpX_valrem(chip, S->nu, p, &b1))
1049
{
1050
if (!p_is_prime(S)) pari_err_PRIME("Decomp",p);
1051
pari_err_BUG("Decomp (not a factor)");
1052
}
1053
b2 = FpX_div(chip, b1, p);
1054
a = FpX_mul(FpXQ_inv(b2, b1, p), b2, p);
1055
/* E = e / de, e in Z[X], de in Z, E = a(phi) mod (f, p) */
1056
th = QpX_remove_denom(S->phi, p, &dt, &vdt);
1057
if (dt)
1058
{
1059
long dega = degpol(a);
1060
vde = dega * vdt;
1061
de = powiu(dt, dega);
1062
pr = mulii(p, de);
1063
a = FpX_rescale(a, dt, pr);
1064
}
1065
else
1066
{
1067
vde = 0;
1068
de = gen_1;
1069
pr = p;
1070
}
1071
e = FpX_FpXQ_eval(a, th, S->f, pr);
1072
update_den(p, &e, &de, &vde, NULL);
1073
1074
pk = p; k = 1;
1075
/* E, (1 - E) tend to orthogonal idempotents in Zp[X]/(f) */
1076
while (k < r + vde)
1077
{ /* E <-- E^2(3-2E) mod p^2k, with E = e/de */
1078
GEN D;
1079
pk = sqri(pk); k <<= 1;
1080
e = ZX_mul(ZX_sqr(e), Z_ZX_sub(mului(3,de), gmul2n(e,1)));
1081
de= mulii(de, sqri(de));
1082
vde *= 3;
1083
D = mulii(pk, de);
1084
e = FpX_rem(e, centermod(S->f, D), D); /* e/de defined mod pk */
1085
update_den(p, &e, &de, &vde, NULL);
1086
}
1087
/* required precision of the factors */
1088
pr = powiu(p, r); pr2 = shifti(pr, -1);
1089
ph = mulii(de,pr);ph2 = shifti(ph, -1);
1090
e = FpX_center_i(FpX_red(e, ph), ph, ph2);
1091
fred = FpX_red(S->f, ph);
1092
1093
f1 = ZpX_gcd(fred, Z_ZX_sub(de, e), p, ph); /* p-adic gcd(f, 1-e) */
1094
if (!is_pm1(de))
1095
{
1096
fred = FpX_red(fred, pr);
1097
f1 = FpX_red(f1, pr);
1098
}
1099
f2 = FpX_div(fred,f1, pr);
1100
f1 = FpX_center_i(f1, pr, pr2);
1101
f2 = FpX_center_i(f2, pr, pr2);
1102
1103
if (DEBUGLEVEL>5)
1104
err_printf(" leaving Decomp: f1 = %Ps\nf2 = %Ps\ne = %Ps\nde= %Ps\n", f1,f2,e,de);
1105
1106
if (flag < 0)
1107
{
1108
GEN m = vconcat(ZpX_primedec(f1, p), ZpX_primedec(f2, p));
1109
return sort_factor(m, (void*)&cmpii, &cmp_nodata);
1110
}
1111
else if (flag)
1112
{
1113
gerepileall(av, 2, &f1, &f2);
1114
return shallowconcat(ZpX_monic_factor_squarefree(f1, p, flag),
1115
ZpX_monic_factor_squarefree(f2, p, flag));
1116
} else {
1117
GEN D, d1, d2, B1, B2, M;
1118
long n, n1, n2, i;
1119
gerepileall(av, 4, &f1, &f2, &e, &de);
1120
D = de;
1121
B1 = get_partial_order_as_pols(p,f1); n1 = lg(B1)-1;
1122
B2 = get_partial_order_as_pols(p,f2); n2 = lg(B2)-1; n = n1+n2;
1123
d1 = QpXV_denom(B1);
1124
d2 = QpXV_denom(B2); if (cmpii(d1, d2) < 0) d1 = d2;
1125
if (d1 != gen_1) {
1126
B1 = Q_muli_to_int(B1, d1);
1127
B2 = Q_muli_to_int(B2, d1);
1128
D = mulii(d1, D);
1129
}
1130
fred = centermod_i(S->f, D, shifti(D,-1));
1131
M = cgetg(n+1, t_MAT);
1132
for (i=1; i<=n1; i++)
1133
gel(M,i) = RgX_to_RgC(FpX_rem(FpX_mul(gel(B1,i),e,D), fred, D), n);
1134
e = Z_ZX_sub(de, e); B2 -= n1;
1135
for ( ; i<=n; i++)
1136
gel(M,i) = RgX_to_RgC(FpX_rem(FpX_mul(gel(B2,i),e,D), fred, D), n);
1137
return ZpM_hnfmodid(M, p, D);
1138
}
1139
}
1140
1141
/* minimum extension valuation: L/E */
1142
static void
1143
vstar(GEN p,GEN h, long *L, long *E)
1144
{
1145
long first, j, k, v, w, m = degpol(h);
1146
1147
first = 1; k = 1; v = 0;
1148
for (j=1; j<=m; j++)
1149
{
1150
GEN c = gel(h, m-j+2);
1151
if (signe(c))
1152
{
1153
w = Z_pval(c,p);
1154
if (first || w*k < v*j) { v = w; k = j; }
1155
first = 0;
1156
}
1157
}
1158
/* v/k = min_j ( v_p(h_{m-j}) / j ) */
1159
w = (long)ugcd(v,k);
1160
*L = v/w;
1161
*E = k/w;
1162
}
1163
1164
static GEN
1165
redelt_i(GEN a, GEN N, GEN p, GEN *pda, long *pvda)
1166
{
1167
GEN z;
1168
a = Q_remove_denom(a, pda);
1169
*pvda = 0;
1170
if (*pda)
1171
{
1172
long v = Z_pvalrem(*pda, p, &z);
1173
if (v) {
1174
*pda = powiu(p, v);
1175
*pvda = v;
1176
N = mulii(*pda, N);
1177
}
1178
else
1179
*pda = NULL;
1180
if (!is_pm1(z)) a = ZX_Z_mul(a, Fp_inv(z, N));
1181
}
1182
return centermod(a, N);
1183
}
1184
/* reduce the element a modulo N [ a power of p ], taking first care of the
1185
* denominators */
1186
static GEN
1187
redelt(GEN a, GEN N, GEN p)
1188
{
1189
GEN da;
1190
long vda;
1191
a = redelt_i(a, N, p, &da, &vda);
1192
if (da) a = RgX_Rg_div(a, da);
1193
return a;
1194
}
1195
1196
/* compute the c first Newton sums modulo pp of the
1197
characteristic polynomial of a/d mod chi, d > 0 power of p (NULL = gen_1),
1198
a, chi in Zp[X], vda = v_p(da)
1199
ns = Newton sums of chi */
1200
static GEN
1201
newtonsums(GEN p, GEN a, GEN da, long vda, GEN chi, long c, GEN pp, GEN ns)
1202
{
1203
GEN va, pa, dpa, s;
1204
long j, k, vdpa, lns = lg(ns);
1205
pari_sp av;
1206
1207
a = centermod(a, pp); av = avma;
1208
dpa = pa = NULL; /* -Wall */
1209
vdpa = 0;
1210
va = zerovec(c);
1211
for (j = 1; j <= c; j++)
1212
{ /* pa/dpa = (a/d)^(j-1) mod (chi, pp), dpa = p^vdpa */
1213
long l;
1214
pa = j == 1? a: FpXQ_mul(pa, a, chi, pp);
1215
l = lg(pa); if (l == 2) break;
1216
if (lns < l) l = lns;
1217
1218
if (da) {
1219
dpa = j == 1? da: mulii(dpa, da);
1220
vdpa += vda;
1221
update_den(p, &pa, &dpa, &vdpa, &pp);
1222
}
1223
s = mulii(gel(pa,2), gel(ns,2)); /* k = 2 */
1224
for (k = 3; k < l; k++) s = addii(s, mulii(gel(pa,k), gel(ns,k)));
1225
if (da) {
1226
GEN r;
1227
s = dvmdii(s, dpa, &r);
1228
if (r != gen_0) return NULL;
1229
}
1230
gel(va,j) = centermodii(s, pp, shifti(pp,-1));
1231
1232
if (gc_needed(av, 1))
1233
{
1234
if(DEBUGMEM>1) pari_warn(warnmem, "newtonsums");
1235
gerepileall(av, dpa?4:3, &pa, &va, &pp, &dpa);
1236
}
1237
}
1238
for (; j <= c; j++) gel(va,j) = gen_0;
1239
return va;
1240
}
1241
1242
/* compute the characteristic polynomial of a/da mod chi (a in Z[X]), given
1243
* by its Newton sums to a precision of pp using Newton sums */
1244
static GEN
1245
newtoncharpoly(GEN pp, GEN p, GEN NS)
1246
{
1247
long n = lg(NS)-1, j, k;
1248
GEN c = cgetg(n + 2, t_VEC), pp2 = shifti(pp,-1);
1249
1250
gel(c,1) = (n & 1 ? gen_m1: gen_1);
1251
for (k = 2; k <= n+1; k++)
1252
{
1253
pari_sp av2 = avma;
1254
GEN s = gen_0;
1255
ulong z;
1256
long v = u_pvalrem(k - 1, p, &z);
1257
for (j = 1; j < k; j++)
1258
{
1259
GEN t = mulii(gel(NS,j), gel(c,k-j));
1260
if (!odd(j)) t = negi(t);
1261
s = addii(s, t);
1262
}
1263
if (v) {
1264
s = gdiv(s, powiu(p, v));
1265
if (typ(s) != t_INT) return NULL;
1266
}
1267
s = mulii(s, Fp_inv(utoipos(z), pp));
1268
gel(c,k) = gerepileuptoint(av2, Fp_center_i(s, pp, pp2));
1269
}
1270
for (k = odd(n)? 1: 2; k <= n+1; k += 2) gel(c,k) = negi(gel(c,k));
1271
return gtopoly(c, 0);
1272
}
1273
1274
static void
1275
manage_cache(decomp_t *S, GEN f, GEN pp)
1276
{
1277
GEN t = S->precns;
1278
1279
if (!t) t = mulii(S->pmf, powiu(S->p, S->df));
1280
if (cmpii(t, pp) < 0) t = pp;
1281
1282
if (!S->precns || !RgX_equal(f, S->nsf) || cmpii(S->precns, t) < 0)
1283
{
1284
if (DEBUGLEVEL>4)
1285
err_printf(" Precision for cached Newton sums for %Ps: %Ps -> %Ps\n",
1286
f, S->precns? S->precns: gen_0, t);
1287
S->nsf = f;
1288
S->ns = FpX_Newton(f, degpol(f), t);
1289
S->precns = t;
1290
}
1291
}
1292
1293
/* return NULL if a mod f is not an integer
1294
* The denominator of any integer in Zp[X]/(f) divides pdr */
1295
static GEN
1296
mycaract(decomp_t *S, GEN f, GEN a, GEN pp, GEN pdr)
1297
{
1298
pari_sp av;
1299
GEN d, chi, prec1, prec2, prec3, ns;
1300
long vd, n = degpol(f);
1301
1302
if (gequal0(a)) return pol_0(varn(f));
1303
1304
a = QpX_remove_denom(a, S->p, &d, &vd);
1305
prec1 = pp;
1306
if (lgefint(S->p) == 3)
1307
prec1 = mulii(prec1, powiu(S->p, factorial_lval(n, itou(S->p))));
1308
if (d)
1309
{
1310
GEN p1 = powiu(d, n);
1311
prec2 = mulii(prec1, p1);
1312
prec3 = mulii(prec1, gmin_shallow(mulii(p1, d), pdr));
1313
}
1314
else
1315
prec2 = prec3 = prec1;
1316
manage_cache(S, f, prec3);
1317
1318
av = avma;
1319
ns = newtonsums(S->p, a, d, vd, f, n, prec2, S->ns);
1320
if (!ns) return NULL;
1321
chi = newtoncharpoly(prec1, S->p, ns);
1322
if (!chi) return NULL;
1323
setvarn(chi, varn(f));
1324
return gerepileupto(av, centermod(chi, pp));
1325
}
1326
1327
static GEN
1328
get_nu(GEN chi, GEN p, long *ptl)
1329
{ /* split off powers of x first for efficiency */
1330
long v = ZX_valrem(FpX_red(chi,p), &chi), n;
1331
GEN P;
1332
if (!degpol(chi)) { *ptl = 1; return pol_x(varn(chi)); }
1333
P = gel(FpX_factor(chi,p), 1); n = lg(P)-1;
1334
*ptl = v? n+1: n; return gel(P,n);
1335
}
1336
1337
/* Factor characteristic polynomial chi of phi mod p. If it splits, update
1338
* S->{phi, chi, nu} and return 1. In any case, set *nu to an irreducible
1339
* factor mod p of chi */
1340
static int
1341
split_char(decomp_t *S, GEN chi, GEN phi, GEN phi0, GEN *nu)
1342
{
1343
long l;
1344
*nu = get_nu(chi, S->p, &l);
1345
if (l == 1) return 0; /* single irreducible factor: doesn't split */
1346
/* phi o phi0 mod (p, f) */
1347
S->phi = compmod(S->p, phi, phi0, S->f, S->p, 0);
1348
S->chi = chi;
1349
S->nu = *nu; return 1;
1350
}
1351
1352
/* Return the prime element in Zp[phi], a t_INT (iff *Ep = 1) or QX;
1353
* nup, chip are ZX. phi = NULL codes X
1354
* If *Ep < oE or Ep divides Ediv (!=0) return NULL (uninteresting) */
1355
static GEN
1356
getprime(decomp_t *S, GEN phi, GEN chip, GEN nup, long *Lp, long *Ep,
1357
long oE, long Ediv)
1358
{
1359
GEN z, chin, q, qp;
1360
long r, s;
1361
1362
if (phi && dvdii(constant_coeff(chip), S->psc))
1363
{
1364
chip = mycaract(S, S->chi, phi, S->pmf, S->prc);
1365
if (dvdii(constant_coeff(chip), S->pmf))
1366
chip = ZXQ_charpoly(phi, S->chi, varn(chip));
1367
}
1368
if (degpol(nup) == 1)
1369
{
1370
GEN c = gel(nup,2); /* nup = X + c */
1371
chin = signe(c)? RgX_translate(chip, negi(c)): chip;
1372
}
1373
else
1374
chin = ZXQ_charpoly(nup, chip, varn(chip));
1375
1376
vstar(S->p, chin, Lp, Ep);
1377
if (*Ep < oE || (Ediv && Ediv % *Ep == 0)) return NULL;
1378
1379
if (*Ep == 1) return S->p;
1380
(void)cbezout(*Lp, -*Ep, &r, &s); /* = 1 */
1381
if (r <= 0)
1382
{
1383
long t = 1 + ((-r) / *Ep);
1384
r += t * *Ep;
1385
s += t * *Lp;
1386
}
1387
/* r > 0 minimal such that r L/E - s = 1/E
1388
* pi = nu^r / p^s is an element of valuation 1/E,
1389
* so is pi + O(p) since 1/E < 1. May compute nu^r mod p^(s+1) */
1390
q = powiu(S->p, s); qp = mulii(q, S->p);
1391
nup = FpXQ_powu(nup, r, S->chi, qp);
1392
if (!phi) return RgX_Rg_div(nup, q); /* phi = X : no composition */
1393
z = compmod(S->p, nup, phi, S->chi, qp, -s);
1394
return signe(z)? z: NULL;
1395
}
1396
1397
static int
1398
update_phi(decomp_t *S)
1399
{
1400
GEN PHI = NULL, prc, psc, X = pol_x(varn(S->f));
1401
long k;
1402
for (k = 1;; k++)
1403
{
1404
prc = ZpX_reduced_resultant_fast(S->chi, ZX_deriv(S->chi), S->p, S->vpsc);
1405
if (!equalii(prc, S->psc)) break;
1406
1407
/* increase precision */
1408
S->vpsc = maxss(S->vpsf, S->vpsc + 1);
1409
S->psc = (S->vpsc == S->vpsf)? S->psf: mulii(S->psc, S->p);
1410
1411
PHI = S->phi;
1412
if (S->phi0) PHI = compmod(S->p, PHI, S->phi0, S->f, S->psc, 0);
1413
PHI = gadd(PHI, ZX_Z_mul(X, mului(k, S->p)));
1414
S->chi = mycaract(S, S->f, PHI, S->psc, S->pdf);
1415
}
1416
psc = mulii(sqri(prc), S->p);
1417
1418
if (!PHI) /* ok above for k = 1 */
1419
{
1420
PHI = S->phi;
1421
if (S->phi0) PHI = compmod(S->p, PHI, S->phi0, S->f, psc, 0);
1422
if (S->phi0 || cmpii(psc,S->psc) > 0)
1423
S->chi = mycaract(S, S->f, PHI, psc, S->pdf);
1424
}
1425
S->phi = PHI;
1426
S->chi = FpX_red(S->chi, psc);
1427
1428
/* may happen if p is unramified */
1429
if (is_pm1(prc)) return 0;
1430
S->psc = psc;
1431
S->vpsc = 2*Z_pval(prc, S->p) + 1;
1432
S->prc = mulii(prc, S->p); return 1;
1433
}
1434
1435
/* return 1 if at least 2 factors mod p ==> chi splits
1436
* Replace S->phi such that F increases (to D) */
1437
static int
1438
testb2(decomp_t *S, long D, GEN theta)
1439
{
1440
long v = varn(S->chi), dlim = degpol(S->chi)-1;
1441
GEN T0 = S->phi, chi, phi, nu;
1442
if (DEBUGLEVEL>4) err_printf(" Increasing Fa\n");
1443
for (;;)
1444
{
1445
phi = gadd(theta, random_FpX(dlim, v, S->p));
1446
chi = mycaract(S, S->chi, phi, S->psf, S->prc);
1447
/* phi nonprimary ? */
1448
if (split_char(S, chi, phi, T0, &nu)) return 1;
1449
if (degpol(nu) == D) break;
1450
}
1451
/* F_phi=lcm(F_alpha, F_theta)=D and E_phi=E_alpha */
1452
S->phi0 = T0;
1453
S->chi = chi;
1454
S->phi = phi;
1455
S->nu = nu; return 0;
1456
}
1457
1458
/* return 1 if at least 2 factors mod p ==> chi can be split.
1459
* compute a new S->phi such that E = lcm(Ea, Et);
1460
* A a ZX, T a t_INT (iff Et = 1, probably impossible ?) or QX */
1461
static int
1462
testc2(decomp_t *S, GEN A, long Ea, GEN T, long Et)
1463
{
1464
GEN c, chi, phi, nu, T0 = S->phi;
1465
1466
if (DEBUGLEVEL>4) err_printf(" Increasing Ea\n");
1467
if (Et == 1) /* same as other branch, split for efficiency */
1468
c = A; /* Et = 1 => s = 1, r = 0, t = 0 */
1469
else
1470
{
1471
long r, s, t;
1472
(void)cbezout(Ea, Et, &r, &s); t = 0;
1473
while (r < 0) { r = r + Et; t++; }
1474
while (s < 0) { s = s + Ea; t++; }
1475
1476
/* A^s T^r / p^t */
1477
c = RgXQ_mul(RgXQ_powu(A, s, S->chi), RgXQ_powu(T, r, S->chi), S->chi);
1478
c = RgX_Rg_div(c, powiu(S->p, t));
1479
c = redelt(c, S->psc, S->p);
1480
}
1481
phi = RgX_add(c, pol_x(varn(S->chi)));
1482
chi = mycaract(S, S->chi, phi, S->psf, S->prc);
1483
if (split_char(S, chi, phi, T0, &nu)) return 1;
1484
/* E_phi = lcm(E_alpha,E_theta) */
1485
S->phi0 = T0;
1486
S->chi = chi;
1487
S->phi = phi;
1488
S->nu = nu; return 0;
1489
}
1490
1491
/* Return h^(-degpol(P)) P(x * h) if result is integral, NULL otherwise */
1492
static GEN
1493
ZX_rescale_inv(GEN P, GEN h)
1494
{
1495
long i, l = lg(P);
1496
GEN Q = cgetg(l,t_POL), hi = h;
1497
gel(Q,l-1) = gel(P,l-1);
1498
for (i=l-2; i>=2; i--)
1499
{
1500
GEN r;
1501
gel(Q,i) = dvmdii(gel(P,i), hi, &r);
1502
if (signe(r)) return NULL;
1503
if (i == 2) break;
1504
hi = mulii(hi,h);
1505
}
1506
Q[1] = P[1]; return Q;
1507
}
1508
1509
/* x p^-eq nu^-er mod p */
1510
static GEN
1511
get_gamma(decomp_t *S, GEN x, long eq, long er)
1512
{
1513
GEN q, g = x, Dg = powiu(S->p, eq);
1514
long vDg = eq;
1515
if (er)
1516
{
1517
if (!S->invnu)
1518
{
1519
while (gdvd(S->chi, S->nu)) S->nu = RgX_Rg_add(S->nu, S->p);
1520
S->invnu = QXQ_inv(S->nu, S->chi);
1521
S->invnu = redelt_i(S->invnu, S->psc, S->p, &S->Dinvnu, &S->vDinvnu);
1522
}
1523
if (S->Dinvnu) {
1524
Dg = mulii(Dg, powiu(S->Dinvnu, er));
1525
vDg += er * S->vDinvnu;
1526
}
1527
q = mulii(S->p, Dg);
1528
g = ZX_mul(g, FpXQ_powu(S->invnu, er, S->chi, q));
1529
g = FpX_rem(g, S->chi, q);
1530
update_den(S->p, &g, &Dg, &vDg, NULL);
1531
g = centermod(g, mulii(S->p, Dg));
1532
}
1533
if (!is_pm1(Dg)) g = RgX_Rg_div(g, Dg);
1534
return g;
1535
}
1536
static GEN
1537
get_g(decomp_t *S, long Ea, long L, long E, GEN beta, GEN *pchig,
1538
long *peq, long *per)
1539
{
1540
long eq, er;
1541
GEN g, chig, chib = NULL;
1542
for(;;) /* at most twice */
1543
{
1544
if (L < 0)
1545
{
1546
chib = ZXQ_charpoly(beta, S->chi, varn(S->chi));
1547
vstar(S->p, chib, &L, &E);
1548
}
1549
eq = L / E; er = L*Ea / E - eq*Ea;
1550
/* floor(L Ea/E) = eq Ea + er */
1551
if (er || !chib)
1552
{ /* g might not be an integer ==> chig = NULL */
1553
g = get_gamma(S, beta, eq, er);
1554
chig = mycaract(S, S->chi, g, S->psc, S->prc);
1555
}
1556
else
1557
{ /* g = beta/p^eq, special case of the above */
1558
GEN h = powiu(S->p, eq);
1559
g = RgX_Rg_div(beta, h);
1560
chig = ZX_rescale_inv(chib, h); /* chib(x h) / h^N */
1561
if (chig) chig = FpX_red(chig, S->pmf);
1562
}
1563
/* either success or second consecutive failure */
1564
if (chig || chib) break;
1565
/* if g fails the v*-test, v(beta) was wrong. Retry once */
1566
L = -1;
1567
}
1568
*pchig = chig; *peq = eq; *per = er; return g;
1569
}
1570
1571
/* return 1 if at least 2 factors mod p ==> chi can be split */
1572
static int
1573
loop(decomp_t *S, long Ea)
1574
{
1575
pari_sp av = avma;
1576
GEN beta = FpXQ_powu(S->nu, Ea, S->chi, S->p);
1577
long N = degpol(S->f), v = varn(S->f);
1578
S->invnu = NULL;
1579
for (;;)
1580
{ /* beta tends to a factor of chi */
1581
long L, i, Fg, eq, er;
1582
GEN chig = NULL, d, g, nug;
1583
1584
if (DEBUGLEVEL>4) err_printf(" beta = %Ps\n", beta);
1585
L = ZpX_resultant_val(S->chi, beta, S->p, S->mf+1);
1586
if (L > S->mf) L = -1; /* from scratch */
1587
g = get_g(S, Ea, L, N, beta, &chig, &eq, &er);
1588
if (DEBUGLEVEL>4) err_printf(" (eq,er) = (%ld,%ld)\n", eq,er);
1589
/* g = beta p^-eq nu^-er (a unit), chig = charpoly(g) */
1590
if (split_char(S, chig, g,S->phi, &nug)) return 1;
1591
1592
Fg = degpol(nug);
1593
if (Fg == 1)
1594
{ /* frequent special case nug = x - d */
1595
long Le, Ee;
1596
GEN chie, nue, e, pie;
1597
d = negi(gel(nug,2));
1598
chie = RgX_translate(chig, d);
1599
nue = pol_x(v);
1600
e = RgX_Rg_sub(g, d);
1601
pie = getprime(S, e, chie, nue, &Le, &Ee, 0,Ea);
1602
if (pie) return testc2(S, S->nu, Ea, pie, Ee);
1603
}
1604
else
1605
{
1606
long Fa = degpol(S->nu), vdeng;
1607
GEN deng, numg, nume;
1608
if (Fa % Fg) return testb2(S, ulcm(Fa,Fg), g);
1609
/* nu & nug irreducible mod p, deg nug | deg nu. To improve beta, look
1610
* for a root d of nug in Fp[phi] such that v_p(g - d) > 0 */
1611
if (ZX_equal(nug, S->nu))
1612
d = pol_x(v);
1613
else
1614
{
1615
if (!p_is_prime(S)) pari_err_PRIME("FpX_ffisom",S->p);
1616
d = FpX_ffisom(nug, S->nu, S->p);
1617
}
1618
/* write g = numg / deng, e = nume / deng */
1619
numg = QpX_remove_denom(g, S->p, &deng, &vdeng);
1620
for (i = 1; i <= Fg; i++)
1621
{
1622
GEN chie, nue, e;
1623
if (i != 1) d = FpXQ_pow(d, S->p, S->nu, S->p); /* next root */
1624
nume = ZX_sub(numg, ZX_Z_mul(d, deng));
1625
/* test e = nume / deng */
1626
if (ZpX_resultant_val(S->chi, nume, S->p, vdeng*N+1) <= vdeng*N)
1627
continue;
1628
e = RgX_Rg_div(nume, deng);
1629
chie = mycaract(S, S->chi, e, S->psc, S->prc);
1630
if (split_char(S, chie, e,S->phi, &nue)) return 1;
1631
if (RgX_is_monomial(nue))
1632
{ /* v_p(e) = v_p(g - d) > 0 */
1633
long Le, Ee;
1634
GEN pie;
1635
pie = getprime(S, e, chie, nue, &Le, &Ee, 0,Ea);
1636
if (pie) return testc2(S, S->nu, Ea, pie, Ee);
1637
break;
1638
}
1639
}
1640
if (i > Fg)
1641
{
1642
if (!p_is_prime(S)) pari_err_PRIME("nilord",S->p);
1643
pari_err_BUG("nilord (no root)");
1644
}
1645
}
1646
if (eq) d = gmul(d, powiu(S->p, eq));
1647
if (er) d = gmul(d, gpowgs(S->nu, er));
1648
beta = gsub(beta, d);
1649
1650
if (gc_needed(av,1))
1651
{
1652
if (DEBUGMEM > 1) pari_warn(warnmem, "nilord");
1653
gerepileall(av, S->invnu? 6: 4, &beta, &(S->precns), &(S->ns), &(S->nsf), &(S->invnu), &(S->Dinvnu));
1654
}
1655
}
1656
}
1657
1658
/* E and F cannot decrease; return 1 if O = Zp[phi], 2 if we can get a
1659
* decomposition and 0 otherwise */
1660
static long
1661
progress(decomp_t *S, GEN *ppa, long *pE)
1662
{
1663
long E = *pE, F;
1664
GEN pa = *ppa;
1665
S->phi0 = NULL; /* no delayed composition */
1666
for(;;)
1667
{
1668
long l, La, Ea; /* N.B If E = 0, getprime cannot return NULL */
1669
GEN pia = getprime(S, NULL, S->chi, S->nu, &La, &Ea, E,0);
1670
if (pia) { /* success, we break out in THIS loop */
1671
pa = (typ(pia) == t_POL)? RgX_RgXQ_eval(pia, S->phi, S->f): pia;
1672
E = Ea;
1673
if (La == 1) break; /* no need to change phi so that nu = pia */
1674
}
1675
/* phi += prime elt */
1676
S->phi = typ(pa) == t_INT? RgX_Rg_add_shallow(S->phi, pa)
1677
: RgX_add(S->phi, pa);
1678
/* recompute char. poly. chi from scratch */
1679
S->chi = mycaract(S, S->f, S->phi, S->psf, S->pdf);
1680
S->nu = get_nu(S->chi, S->p, &l);
1681
if (l > 1) return 2;
1682
if (!update_phi(S)) return 1; /* unramified */
1683
if (pia) break;
1684
}
1685
*pE = E; *ppa = pa; F = degpol(S->nu);
1686
if (DEBUGLEVEL>4) err_printf(" (E, F) = (%ld,%ld)\n", E, F);
1687
if (E * F == degpol(S->f)) return 1;
1688
if (loop(S, E)) return 2;
1689
if (!update_phi(S)) return 1;
1690
return 0;
1691
}
1692
1693
/* flag != 0 iff we're looking for the p-adic factorization,
1694
in which case it is the p-adic precision we want */
1695
static GEN
1696
maxord_i(decomp_t *S, GEN p, GEN f, long mf, GEN w, long flag)
1697
{
1698
long oE, n = lg(w)-1; /* factor of largest degree */
1699
GEN opa, D = ZpX_reduced_resultant_fast(f, ZX_deriv(f), p, mf);
1700
S->pisprime = -1;
1701
S->p = p;
1702
S->mf = mf;
1703
S->nu = gel(w,n);
1704
S->df = Z_pval(D, p);
1705
S->pdf = powiu(p, S->df);
1706
S->phi = pol_x(varn(f));
1707
S->chi = S->f = f;
1708
if (n > 1) return Decomp(S, flag); /* FIXME: use bezout_lift_fact */
1709
1710
if (DEBUGLEVEL>4)
1711
err_printf(" entering Nilord: %Ps^%ld\n f = %Ps, nu = %Ps\n",
1712
p, S->df, S->f, S->nu);
1713
else if (DEBUGLEVEL>2) err_printf(" entering Nilord\n");
1714
S->psf = S->psc = mulii(sqri(D), p);
1715
S->vpsf = S->vpsc = 2*S->df + 1;
1716
S->prc = mulii(D, p);
1717
S->chi = FpX_red(S->f, S->psc);
1718
S->pmf = powiu(p, S->mf+1);
1719
S->precns = NULL;
1720
for(opa = NULL, oE = 0;;)
1721
{
1722
long n = progress(S, &opa, &oE);
1723
if (n == 1) return flag? NULL: dbasis(p, S->f, S->mf, S->phi, S->chi);
1724
if (n == 2) return Decomp(S, flag);
1725
}
1726
}
1727
1728
static int
1729
expo_is_squarefree(GEN e)
1730
{
1731
long i, l = lg(e);
1732
for (i=1; i<l; i++)
1733
if (e[i] != 1) return 0;
1734
return 1;
1735
}
1736
/* pure round 4 */
1737
static GEN
1738
ZpX_round4(GEN f, GEN p, GEN w, long prec)
1739
{
1740
decomp_t S;
1741
GEN L = maxord_i(&S, p, f, ZpX_disc_val(f,p), w, prec);
1742
return L? L: mkvec(f);
1743
}
1744
/* f a squarefree ZX with leading_coeff 1, degree > 0. Return list of
1745
* irreducible factors in Zp[X] (computed mod p^prec) */
1746
static GEN
1747
ZpX_monic_factor_squarefree(GEN f, GEN p, long prec)
1748
{
1749
pari_sp av = avma;
1750
GEN L, fa, w, e;
1751
long i, l;
1752
if (degpol(f) == 1) return mkvec(f);
1753
fa = FpX_factor(f,p); w = gel(fa,1); e = gel(fa,2);
1754
/* no repeated factors: Hensel lift */
1755
if (expo_is_squarefree(e)) return ZpX_liftfact(f, w, powiu(p,prec), p, prec);
1756
l = lg(w);
1757
if (l == 2)
1758
{
1759
L = ZpX_round4(f,p,w,prec);
1760
if (lg(L) == 2) { set_avma(av); return mkvec(f); }
1761
}
1762
else
1763
{ /* >= 2 factors mod p: partial Hensel lift */
1764
GEN D = ZpX_reduced_resultant_fast(f, ZX_deriv(f), p, ZpX_disc_val(f,p));
1765
long r = maxss(2*Z_pval(D,p)+1, prec);
1766
GEN W = cgetg(l, t_VEC);
1767
for (i = 1; i < l; i++)
1768
gel(W,i) = e[i] == 1? gel(w,i): FpX_powu(gel(w,i), e[i], p);
1769
L = ZpX_liftfact(f, W, powiu(p,r), p, r);
1770
for (i = 1; i < l; i++)
1771
gel(L,i) = e[i] == 1? mkvec(gel(L,i))
1772
: ZpX_round4(gel(L,i), p, mkvec(gel(w,i)), prec);
1773
L = shallowconcat1(L);
1774
}
1775
return gerepilecopy(av, L);
1776
}
1777
1778
/* assume f a ZX with leading_coeff 1, degree > 0 */
1779
GEN
1780
ZpX_monic_factor(GEN f, GEN p, long prec)
1781
{
1782
GEN poly, ex, P, E;
1783
long l, i;
1784
1785
if (degpol(f) == 1) return mkmat2(mkcol(f), mkcol(gen_1));
1786
poly = ZX_squff(f,&ex); l = lg(poly);
1787
P = cgetg(l, t_VEC);
1788
E = cgetg(l, t_VEC);
1789
for (i = 1; i < l; i++)
1790
{
1791
GEN L = ZpX_monic_factor_squarefree(gel(poly,i), p, prec);
1792
gel(P,i) = L; settyp(L, t_COL);
1793
gel(E,i) = const_col(lg(L)-1, utoipos(ex[i]));
1794
}
1795
return mkmat2(shallowconcat1(P), shallowconcat1(E));
1796
}
1797
1798
/* DT = multiple of disc(T) or NULL
1799
* Return a multiple of the denominator of an algebraic integer (in Q[X]/(T))
1800
* when expressed in terms of the power basis */
1801
GEN
1802
indexpartial(GEN T, GEN DT)
1803
{
1804
pari_sp av = avma;
1805
long i, nb;
1806
GEN fa, E, P, U, res = gen_1, dT = ZX_deriv(T);
1807
1808
if (!DT) DT = ZX_disc(T);
1809
fa = absZ_factor_limit_strict(DT, 0, &U);
1810
P = gel(fa,1);
1811
E = gel(fa,2); nb = lg(P)-1;
1812
for (i = 1; i <= nb; i++)
1813
{
1814
long e = itou(gel(E,i)), e2 = e >> 1;
1815
GEN p = gel(P,i), q = p;
1816
if (e2 >= 2) q = ZpX_reduced_resultant_fast(T, dT, p, e2);
1817
res = mulii(res, q);
1818
}
1819
if (U)
1820
{
1821
long e = itou(gel(U,2)), e2 = e >> 1;
1822
GEN p = gel(U,1), q = powiu(p, odd(e)? e2+1: e2);
1823
res = mulii(res, q);
1824
}
1825
return gerepileuptoint(av,res);
1826
}
1827
1828
/*******************************************************************/
1829
/* */
1830
/* 2-ELT REPRESENTATION FOR PRIME IDEALS (dividing index) */
1831
/* */
1832
/*******************************************************************/
1833
/* to compute norm of elt in basis form */
1834
typedef struct {
1835
long r1;
1836
GEN M; /* via embed_norm */
1837
1838
GEN D, w, T; /* via resultant if M = NULL */
1839
} norm_S;
1840
1841
static GEN
1842
get_norm(norm_S *S, GEN a)
1843
{
1844
if (S->M)
1845
{
1846
long e;
1847
GEN N = grndtoi( embed_norm(RgM_RgC_mul(S->M, a), S->r1), &e );
1848
if (e > -5) pari_err_PREC( "get_norm");
1849
return N;
1850
}
1851
if (S->w) a = RgV_RgC_mul(S->w, a);
1852
return ZX_resultant_all(S->T, a, S->D, 0);
1853
}
1854
static void
1855
init_norm(norm_S *S, GEN nf, GEN p)
1856
{
1857
GEN T = nf_get_pol(nf), M = nf_get_M(nf);
1858
long N = degpol(T), ex = gexpo(M) + gexpo(mului(8 * N, p));
1859
1860
S->r1 = nf_get_r1(nf);
1861
if (N * ex <= prec2nbits(gprecision(M)) - 20)
1862
{ /* enough prec to use embed_norm */
1863
S->M = M;
1864
S->D = NULL;
1865
S->w = NULL;
1866
S->T = NULL;
1867
}
1868
else
1869
{
1870
GEN w = leafcopy(nf_get_zkprimpart(nf)), D = nf_get_zkden(nf), Dp = sqri(p);
1871
long i;
1872
if (!equali1(D))
1873
{
1874
GEN w1 = D;
1875
long v = Z_pval(D, p);
1876
D = powiu(p, v);
1877
Dp = mulii(D, Dp);
1878
gel(w, 1) = remii(w1, Dp);
1879
}
1880
for (i=2; i<=N; i++) gel(w,i) = FpX_red(gel(w,i), Dp);
1881
S->M = NULL;
1882
S->D = D;
1883
S->w = w;
1884
S->T = T;
1885
}
1886
}
1887
/* f = f(pr/p), q = p^(f+1), a in pr.
1888
* Return 1 if v_pr(a) = 1, and 0 otherwise */
1889
static int
1890
is_uniformizer(GEN a, GEN q, norm_S *S) { return !dvdii(get_norm(S,a), q); }
1891
1892
/* Return x * y, x, y are t_MAT (Fp-basis of in O_K/p), assume (x,y)=1.
1893
* Either x or y may be NULL (= O_K), not both */
1894
static GEN
1895
mul_intersect(GEN x, GEN y, GEN p)
1896
{
1897
if (!x) return y;
1898
if (!y) return x;
1899
return FpM_intersect_i(x, y, p);
1900
}
1901
/* Fp-basis of (ZK/pr): applied to the primes found in primedec_aux()
1902
* true nf */
1903
static GEN
1904
Fp_basis(GEN nf, GEN pr)
1905
{
1906
long i, j, l;
1907
GEN x, y;
1908
/* already in basis form (from Buchman-Lenstra) ? */
1909
if (typ(pr) == t_MAT) return pr;
1910
/* ordinary prid (from Kummer) */
1911
x = pr_hnf(nf, pr);
1912
l = lg(x);
1913
y = cgetg(l, t_MAT);
1914
for (i=j=1; i<l; i++)
1915
if (gequal1(gcoeff(x,i,i))) gel(y,j++) = gel(x,i);
1916
setlg(y, j); return y;
1917
}
1918
/* Let Ip = prod_{ P | p } P be the p-radical. The list L contains the
1919
* P (mod Ip) seen as sub-Fp-vector spaces of ZK/Ip.
1920
* Return the list of (Ip / P) (mod Ip).
1921
* N.B: All ideal multiplications are computed as intersections of Fp-vector
1922
* spaces. true nf */
1923
static GEN
1924
get_LV(GEN nf, GEN L, GEN p, long N)
1925
{
1926
long i, l = lg(L)-1;
1927
GEN LV, LW, A, B;
1928
1929
LV = cgetg(l+1, t_VEC);
1930
if (l == 1) { gel(LV,1) = matid(N); return LV; }
1931
LW = cgetg(l+1, t_VEC);
1932
for (i=1; i<=l; i++) gel(LW,i) = Fp_basis(nf, gel(L,i));
1933
1934
/* A[i] = L[1]...L[i-1], i = 2..l */
1935
A = cgetg(l+1, t_VEC); gel(A,1) = NULL;
1936
for (i=1; i < l; i++) gel(A,i+1) = mul_intersect(gel(A,i), gel(LW,i), p);
1937
/* B[i] = L[i+1]...L[l], i = 1..(l-1) */
1938
B = cgetg(l+1, t_VEC); gel(B,l) = NULL;
1939
for (i=l; i>=2; i--) gel(B,i-1) = mul_intersect(gel(B,i), gel(LW,i), p);
1940
for (i=1; i<=l; i++) gel(LV,i) = mul_intersect(gel(A,i), gel(B,i), p);
1941
return LV;
1942
}
1943
1944
static void
1945
errprime(GEN p) { pari_err_PRIME("idealprimedec",p); }
1946
1947
/* P = Fp-basis (over O_K/p) for pr.
1948
* V = Z-basis for I_p/pr. ramif != 0 iff some pr|p is ramified.
1949
* Return a p-uniformizer for pr. Assume pr not inert, i.e. m > 0 */
1950
static GEN
1951
uniformizer(GEN nf, norm_S *S, GEN P, GEN V, GEN p, int ramif)
1952
{
1953
long i, l, f, m = lg(P)-1, N = nf_get_degree(nf);
1954
GEN u, Mv, x, q;
1955
1956
f = N - m; /* we want v_p(Norm(x)) = p^f */
1957
q = powiu(p,f+1);
1958
1959
u = FpM_FpC_invimage(shallowconcat(P, V), col_ei(N,1), p);
1960
setlg(u, lg(P));
1961
u = centermod(ZM_ZC_mul(P, u), p);
1962
if (is_uniformizer(u, q, S)) return u;
1963
if (signe(gel(u,1)) <= 0) /* make sure u[1] in ]-p,p] */
1964
gel(u,1) = addii(gel(u,1), p); /* try u + p */
1965
else
1966
gel(u,1) = subii(gel(u,1), p); /* try u - p */
1967
if (!ramif || is_uniformizer(u, q, S)) return u;
1968
1969
/* P/p ramified, u in P^2, not in Q for all other Q|p */
1970
Mv = zk_multable(nf, Z_ZC_sub(gen_1,u));
1971
l = lg(P);
1972
for (i=1; i<l; i++)
1973
{
1974
x = centermod(ZC_add(u, ZM_ZC_mul(Mv, gel(P,i))), p);
1975
if (is_uniformizer(x, q, S)) return x;
1976
}
1977
errprime(p);
1978
return NULL; /* LCOV_EXCL_LINE */
1979
}
1980
1981
/*******************************************************************/
1982
/* */
1983
/* BUCHMANN-LENSTRA ALGORITHM */
1984
/* */
1985
/*******************************************************************/
1986
static GEN
1987
mk_pr(GEN p, GEN u, long e, long f, GEN t)
1988
{ return mkvec5(p, u, utoipos(e), utoipos(f), t); }
1989
1990
/* nf a true nf, u in Z[X]/(T); pr = p Z_K + u Z_K of ramification index e */
1991
GEN
1992
idealprimedec_kummer(GEN nf,GEN u,long e,GEN p)
1993
{
1994
GEN t, T = nf_get_pol(nf);
1995
long f = degpol(u), N = degpol(T);
1996
1997
if (f == N)
1998
{ /* inert */
1999
u = scalarcol_shallow(p,N);
2000
t = gen_1;
2001
}
2002
else
2003
{
2004
t = centermod(poltobasis(nf, FpX_div(T, u, p)), p);
2005
u = centermod(poltobasis(nf, u), p);
2006
if (e == 1)
2007
{ /* make sure v_pr(u) = 1 (automatic if e>1) */
2008
GEN cw, w = Q_primitive_part(nf_to_scalar_or_alg(nf, u), &cw);
2009
long v = cw? f - Q_pval(cw, p) * N: f;
2010
if (ZpX_resultant_val(T, w, p, v + 1) > v)
2011
{
2012
GEN c = gel(u,1);
2013
gel(u,1) = signe(c) > 0? subii(c, p): addii(c, p);
2014
}
2015
}
2016
t = zk_multable(nf, t);
2017
}
2018
return mk_pr(p,u,e,f,t);
2019
}
2020
2021
typedef struct {
2022
GEN nf, p;
2023
long I;
2024
} eltmod_muldata;
2025
2026
static GEN
2027
sqr_mod(void *data, GEN x)
2028
{
2029
eltmod_muldata *D = (eltmod_muldata*)data;
2030
return FpC_red(nfsqri(D->nf, x), D->p);
2031
}
2032
static GEN
2033
ei_msqr_mod(void *data, GEN x)
2034
{
2035
GEN x2 = sqr_mod(data, x);
2036
eltmod_muldata *D = (eltmod_muldata*)data;
2037
return FpC_red(zk_ei_mul(D->nf, x2, D->I), D->p);
2038
}
2039
/* nf a true nf; compute lift(nf.zk[I]^p mod p) */
2040
static GEN
2041
pow_ei_mod_p(GEN nf, long I, GEN p)
2042
{
2043
pari_sp av = avma;
2044
eltmod_muldata D;
2045
long N = nf_get_degree(nf);
2046
GEN y = col_ei(N,I);
2047
if (I == 1) return y;
2048
D.nf = nf;
2049
D.p = p;
2050
D.I = I;
2051
y = gen_pow_fold(y, p, (void*)&D, &sqr_mod, &ei_msqr_mod);
2052
return gerepileupto(av,y);
2053
}
2054
2055
/* nf a true nf; return a Z basis of Z_K's p-radical, phi = x--> x^p-x */
2056
static GEN
2057
pradical(GEN nf, GEN p, GEN *phi)
2058
{
2059
long i, N = nf_get_degree(nf);
2060
GEN q,m,frob,rad;
2061
2062
/* matrix of Frob: x->x^p over Z_K/p */
2063
frob = cgetg(N+1,t_MAT);
2064
for (i=1; i<=N; i++) gel(frob,i) = pow_ei_mod_p(nf,i,p);
2065
2066
m = frob; q = p;
2067
while (abscmpiu(q,N) < 0) { q = mulii(q,p); m = FpM_mul(m, frob, p); }
2068
rad = FpM_ker(m, p); /* m = Frob^k, s.t p^k >= N */
2069
for (i=1; i<=N; i++) gcoeff(frob,i,i) = subiu(gcoeff(frob,i,i), 1);
2070
*phi = frob; return rad;
2071
}
2072
2073
/* return powers of a: a^0, ... , a^d, d = dim A */
2074
static GEN
2075
get_powers(GEN mul, GEN p)
2076
{
2077
long i, d = lgcols(mul);
2078
GEN z, pow = cgetg(d+2,t_MAT), P = pow+1;
2079
2080
gel(P,0) = scalarcol_shallow(gen_1, d-1);
2081
z = gel(mul,1);
2082
for (i=1; i<=d; i++)
2083
{
2084
gel(P,i) = z; /* a^i */
2085
if (i!=d) z = FpM_FpC_mul(mul, z, p);
2086
}
2087
return pow;
2088
}
2089
2090
/* minimal polynomial of a in A (dim A = d).
2091
* mul = multiplication table by a in A */
2092
static GEN
2093
pol_min(GEN mul, GEN p)
2094
{
2095
pari_sp av = avma;
2096
GEN z = FpM_deplin(get_powers(mul, p), p);
2097
return gerepilecopy(av, RgV_to_RgX(z,0));
2098
}
2099
2100
static GEN
2101
get_pr(GEN nf, norm_S *S, GEN p, GEN P, GEN V, int ramif, long N, long flim)
2102
{
2103
GEN u, t;
2104
long e, f;
2105
2106
if (typ(P) == t_VEC)
2107
{ /* already done (Kummer) */
2108
f = pr_get_f(P);
2109
if (flim > 0 && f > flim) return NULL;
2110
if (flim == -2) return (GEN)f;
2111
return P;
2112
}
2113
f = N - (lg(P)-1);
2114
if (flim > 0 && f > flim) return NULL;
2115
if (flim == -2) return (GEN)f;
2116
/* P = (p,u) prime. t is an anti-uniformizer: Z_K + t/p Z_K = P^(-1),
2117
* so that v_P(t) = e(P/p)-1 */
2118
if (f == N) {
2119
u = scalarcol_shallow(p,N);
2120
t = gen_1;
2121
e = 1;
2122
} else {
2123
GEN mt;
2124
u = uniformizer(nf, S, P, V, p, ramif);
2125
t = FpM_deplin(zk_multable(nf,u), p);
2126
mt = zk_multable(nf, t);
2127
e = ramif? 1 + ZC_nfval(t,mk_pr(p,u,0,0,mt)): 1;
2128
t = mt;
2129
}
2130
return mk_pr(p,u,e,f,t);
2131
}
2132
2133
/* true nf */
2134
static GEN
2135
primedec_end(GEN nf, GEN L, GEN p, long flim)
2136
{
2137
long i, j, l = lg(L), N = nf_get_degree(nf);
2138
GEN LV = get_LV(nf, L,p,N);
2139
int ramif = dvdii(nf_get_disc(nf), p);
2140
norm_S S; init_norm(&S, nf, p);
2141
for (i = j = 1; i < l; i++)
2142
{
2143
GEN P = get_pr(nf, &S, p, gel(L,i), gel(LV,i), ramif, N, flim);
2144
if (!P) continue;
2145
gel(L,j++) = P;
2146
if (flim == -1) return P;
2147
}
2148
setlg(L, j); return L;
2149
}
2150
2151
/* prime ideal decomposition of p; if flim>0, restrict to f(P,p) <= flim
2152
* if flim = -1 return only the first P
2153
* if flim = -2 return only the f(P/p) in a t_VECSMALL */
2154
static GEN
2155
primedec_aux(GEN nf, GEN p, long flim)
2156
{
2157
const long TYP = (flim == -2)? t_VECSMALL: t_VEC;
2158
GEN E, F, L, Ip, phi, f, g, h, UN, T = nf_get_pol(nf);
2159
long i, k, c, iL, N;
2160
int kummer;
2161
2162
F = FpX_factor(T, p);
2163
E = gel(F,2);
2164
F = gel(F,1);
2165
2166
k = lg(F); if (k == 1) errprime(p);
2167
if ( !dvdii(nf_get_index(nf),p) ) /* p doesn't divide index */
2168
{
2169
L = cgetg(k, TYP);
2170
for (i=1; i<k; i++)
2171
{
2172
GEN t = gel(F,i);
2173
long f = degpol(t);
2174
if (flim > 0 && f > flim) { setlg(L, i); break; }
2175
if (flim == -2)
2176
L[i] = f;
2177
else
2178
gel(L,i) = idealprimedec_kummer(nf, t, E[i],p);
2179
if (flim == -1) return gel(L,1);
2180
}
2181
return L;
2182
}
2183
2184
kummer = 0;
2185
g = FpXV_prod(F, p);
2186
h = FpX_div(T,g,p);
2187
f = FpX_red(ZX_Z_divexact(ZX_sub(ZX_mul(g,h), T), p), p);
2188
2189
N = degpol(T);
2190
L = cgetg(N+1,TYP);
2191
iL = 1;
2192
for (i=1; i<k; i++)
2193
if (E[i] == 1 || signe(FpX_rem(f,gel(F,i),p)))
2194
{
2195
GEN t = gel(F,i);
2196
kummer = 1;
2197
gel(L,iL++) = idealprimedec_kummer(nf, t, E[i],p);
2198
if (flim == -1) return gel(L,1);
2199
}
2200
else /* F[i] | (f,g,h), happens at least once by Dedekind criterion */
2201
E[i] = 0;
2202
2203
/* phi matrix of x -> x^p - x in algebra Z_K/p */
2204
Ip = pradical(nf,p,&phi);
2205
2206
/* split etale algebra Z_K / (p,Ip) */
2207
h = cgetg(N+1,t_VEC);
2208
if (kummer)
2209
{ /* split off Kummer factors */
2210
GEN mb, b = NULL;
2211
for (i=1; i<k; i++)
2212
if (!E[i]) b = b? FpX_mul(b, gel(F,i), p): gel(F,i);
2213
if (!b) errprime(p);
2214
b = FpC_red(poltobasis(nf,b), p);
2215
mb = FpM_red(zk_multable(nf,b), p);
2216
/* Fp-base of ideal (Ip, b) in ZK/p */
2217
gel(h,1) = FpM_image(shallowconcat(mb,Ip), p);
2218
}
2219
else
2220
gel(h,1) = Ip;
2221
2222
UN = col_ei(N, 1);
2223
for (c=1; c; c--)
2224
{ /* Let A:= (Z_K/p) / Ip etale; split A2 := A / Im H ~ Im M2
2225
H * ? + M2 * Mi2 = Id_N ==> M2 * Mi2 projector A --> A2 */
2226
GEN M, Mi, M2, Mi2, phi2, mat1, H = gel(h,c); /* maximal rank */
2227
long dim, r = lg(H)-1;
2228
2229
M = FpM_suppl(shallowconcat(H,UN), p);
2230
Mi = FpM_inv(M, p);
2231
M2 = vecslice(M, r+1,N); /* M = (H|M2) invertible */
2232
Mi2 = rowslice(Mi,r+1,N);
2233
/* FIXME: FpM_mul(,M2) could be done with vecpermute */
2234
phi2 = FpM_mul(Mi2, FpM_mul(phi,M2, p), p);
2235
mat1 = FpM_ker(phi2, p);
2236
dim = lg(mat1)-1; /* A2 product of 'dim' fields */
2237
if (dim > 1)
2238
{ /* phi2 v = 0 => a = M2 v in Ker phi, a not in Fp.1 + H */
2239
GEN R, a, mula, mul2, v = gel(mat1,2);
2240
long n;
2241
2242
a = FpM_FpC_mul(M2,v, p); /* not a scalar */
2243
mula = FpM_red(zk_multable(nf,a), p);
2244
mul2 = FpM_mul(Mi2, FpM_mul(mula,M2, p), p);
2245
R = FpX_roots(pol_min(mul2,p), p); /* totally split mod p */
2246
n = lg(R)-1;
2247
for (i=1; i<=n; i++)
2248
{
2249
GEN I = RgM_Rg_sub_shallow(mula, gel(R,i));
2250
gel(h,c++) = FpM_image(shallowconcat(H, I), p);
2251
}
2252
if (n == dim)
2253
for (i=1; i<=n; i++) gel(L,iL++) = gel(h,--c);
2254
}
2255
else /* A2 field ==> H maximal, f = N-r = dim(A2) */
2256
gel(L,iL++) = H;
2257
}
2258
setlg(L, iL);
2259
return primedec_end(nf, L, p, flim);
2260
}
2261
2262
GEN
2263
idealprimedec_limit_f(GEN nf, GEN p, long f)
2264
{
2265
pari_sp av = avma;
2266
GEN v;
2267
if (typ(p) != t_INT) pari_err_TYPE("idealprimedec",p);
2268
if (f < 0) pari_err_DOMAIN("idealprimedec", "f", "<", gen_0, stoi(f));
2269
v = primedec_aux(checknf(nf), p, f);
2270
v = gen_sort(v, (void*)&cmp_prime_over_p, &cmp_nodata);
2271
return gerepileupto(av,v);
2272
}
2273
GEN
2274
idealprimedec_galois(GEN nf, GEN p)
2275
{
2276
pari_sp av = avma;
2277
GEN v = primedec_aux(checknf(nf), p, -1);
2278
return gerepilecopy(av,v);
2279
}
2280
GEN
2281
idealprimedec_degrees(GEN nf, GEN p)
2282
{
2283
pari_sp av = avma;
2284
GEN v = primedec_aux(checknf(nf), p, -2);
2285
vecsmall_sort(v); return gerepileuptoleaf(av, v);
2286
}
2287
GEN
2288
idealprimedec_limit_norm(GEN nf, GEN p, GEN B)
2289
{ return idealprimedec_limit_f(nf, p, logint(B,p)); }
2290
GEN
2291
idealprimedec(GEN nf, GEN p)
2292
{ return idealprimedec_limit_f(nf, p, 0); }
2293
GEN
2294
nf_pV_to_prV(GEN nf, GEN P)
2295
{
2296
long i, l;
2297
GEN Q = cgetg_copy(P,&l);
2298
if (l == 1) return Q;
2299
for (i = 1; i < l; i++) gel(Q,i) = idealprimedec(nf, gel(P,i));
2300
return shallowconcat1(Q);
2301
}
2302
2303
/* return [Fp[x]: Fp] */
2304
static long
2305
ffdegree(GEN x, GEN frob, GEN p)
2306
{
2307
pari_sp av = avma;
2308
long d, f = lg(frob)-1;
2309
GEN y = x;
2310
2311
for (d=1; d < f; d++)
2312
{
2313
y = FpM_FpC_mul(frob, y, p);
2314
if (ZV_equal(y, x)) break;
2315
}
2316
return gc_long(av,d);
2317
}
2318
2319
static GEN
2320
lift_to_zk(GEN v, GEN c, long N)
2321
{
2322
GEN w = zerocol(N);
2323
long i, l = lg(c);
2324
for (i=1; i<l; i++) gel(w,c[i]) = gel(v,i);
2325
return w;
2326
}
2327
2328
/* return t = 1 mod pr, t = 0 mod p / pr^e(pr/p) */
2329
static GEN
2330
anti_uniformizer(GEN nf, GEN pr)
2331
{
2332
long N = nf_get_degree(nf), e = pr_get_e(pr);
2333
GEN p, b, z;
2334
2335
if (e * pr_get_f(pr) == N) return gen_1;
2336
p = pr_get_p(pr);
2337
b = pr_get_tau(pr); /* ZM */
2338
if (e != 1)
2339
{
2340
GEN q = powiu(pr_get_p(pr), e-1);
2341
b = ZM_Z_divexact(ZM_powu(b,e), q);
2342
}
2343
/* b = tau^e / p^(e-1), v_pr(b) = 0, v_Q(b) >= e(Q/p) for other Q | p */
2344
z = ZM_hnfmodid(FpM_red(b,p), p); /* ideal (p) / pr^e, coprime to pr */
2345
z = idealaddtoone_raw(nf, pr, z);
2346
return Z_ZC_sub(gen_1, FpC_center(FpC_red(z,p), p, shifti(p,-1)));
2347
}
2348
2349
#define mpr_TAU 1
2350
#define mpr_FFP 2
2351
#define mpr_NFP 5
2352
#define SMALLMODPR 4
2353
#define LARGEMODPR 6
2354
static GEN
2355
modpr_TAU(GEN modpr)
2356
{
2357
GEN tau = gel(modpr,mpr_TAU);
2358
return isintzero(tau)? NULL: tau;
2359
}
2360
2361
/* prh = HNF matrix, which is identity but for the first line. Return a
2362
* projector to ZK / prh ~ Z/prh[1,1] */
2363
GEN
2364
dim1proj(GEN prh)
2365
{
2366
long i, N = lg(prh)-1;
2367
GEN ffproj = cgetg(N+1, t_VEC);
2368
GEN x, q = gcoeff(prh,1,1);
2369
gel(ffproj,1) = gen_1;
2370
for (i=2; i<=N; i++)
2371
{
2372
x = gcoeff(prh,1,i);
2373
if (signe(x)) x = subii(q,x);
2374
gel(ffproj,i) = x;
2375
}
2376
return ffproj;
2377
}
2378
2379
/* p not necessarily prime, but coprime to denom(basis) */
2380
GEN
2381
QXQV_to_FpM(GEN basis, GEN T, GEN p)
2382
{
2383
long i, l = lg(basis), f = degpol(T);
2384
GEN z = cgetg(l, t_MAT);
2385
for (i = 1; i < l; i++)
2386
{
2387
GEN w = gel(basis,i);
2388
if (typ(w) == t_INT)
2389
w = scalarcol_shallow(w, f);
2390
else
2391
{
2392
GEN dx;
2393
w = Q_remove_denom(w, &dx);
2394
w = FpXQ_red(w, T, p);
2395
if (dx)
2396
{
2397
dx = Fp_inv(dx, p);
2398
if (!equali1(dx)) w = FpX_Fp_mul(w, dx, p);
2399
}
2400
w = RgX_to_RgC(w, f);
2401
}
2402
gel(z,i) = w; /* w_i mod (T,p) */
2403
}
2404
return z;
2405
}
2406
2407
/* initialize reduction mod pr; if zk = 1, will only init data required to
2408
* reduce *integral* element. Realize (O_K/pr) as Fp[X] / (T), for a
2409
* *monic* T; use variable vT for varn(T) */
2410
static GEN
2411
modprinit(GEN nf, GEN pr, int zk, long vT)
2412
{
2413
pari_sp av = avma;
2414
GEN res, tau, mul, x, p, T, pow, ffproj, nfproj, prh, c;
2415
long N, i, k, f;
2416
2417
nf = checknf(nf); checkprid(pr);
2418
if (vT < 0) vT = nf_get_varn(nf);
2419
f = pr_get_f(pr);
2420
N = nf_get_degree(nf);
2421
prh = pr_hnf(nf, pr);
2422
tau = zk? gen_0: anti_uniformizer(nf, pr);
2423
p = pr_get_p(pr);
2424
2425
if (f == 1)
2426
{
2427
res = cgetg(SMALLMODPR, t_COL);
2428
gel(res,mpr_TAU) = tau;
2429
gel(res,mpr_FFP) = dim1proj(prh);
2430
gel(res,3) = pr; return gerepilecopy(av, res);
2431
}
2432
2433
c = cgetg(f+1, t_VECSMALL);
2434
ffproj = cgetg(N+1, t_MAT);
2435
for (k=i=1; i<=N; i++)
2436
{
2437
x = gcoeff(prh, i,i);
2438
if (!is_pm1(x)) { c[k] = i; gel(ffproj,i) = col_ei(N, i); k++; }
2439
else
2440
gel(ffproj,i) = ZC_neg(gel(prh,i));
2441
}
2442
ffproj = rowpermute(ffproj, c);
2443
if (! dvdii(nf_get_index(nf), p))
2444
{
2445
GEN basis = nf_get_zkprimpart(nf), D = nf_get_zkden(nf);
2446
if (N == f)
2447
{ /* pr inert */
2448
T = nf_get_pol(nf);
2449
T = FpX_red(T,p);
2450
ffproj = RgV_to_RgM(basis, lg(basis)-1);
2451
}
2452
else
2453
{
2454
T = RgV_RgC_mul(basis, pr_get_gen(pr));
2455
T = FpX_normalize(FpX_red(T,p),p);
2456
basis = FqV_red(vecpermute(basis,c), T, p);
2457
basis = RgV_to_RgM(basis, lg(basis)-1);
2458
ffproj = ZM_mul(basis, ffproj);
2459
}
2460
setvarn(T, vT);
2461
ffproj = FpM_red(ffproj, p);
2462
if (!equali1(D))
2463
{
2464
D = modii(D,p);
2465
if (!equali1(D)) ffproj = FpM_Fp_mul(ffproj, Fp_inv(D,p), p);
2466
}
2467
2468
res = cgetg(SMALLMODPR+1, t_COL);
2469
gel(res,mpr_TAU) = tau;
2470
gel(res,mpr_FFP) = ffproj;
2471
gel(res,3) = pr;
2472
gel(res,4) = T; return gerepilecopy(av, res);
2473
}
2474
2475
if (uisprime(f))
2476
{
2477
mul = ei_multable(nf, c[2]);
2478
mul = vecpermute(mul, c);
2479
}
2480
else
2481
{
2482
GEN v, u, u2, frob;
2483
long deg,deg1,deg2;
2484
2485
/* matrix of Frob: x->x^p over Z_K/pr = < w[c1], ..., w[cf] > over Fp */
2486
frob = cgetg(f+1, t_MAT);
2487
for (i=1; i<=f; i++)
2488
{
2489
x = pow_ei_mod_p(nf,c[i],p);
2490
gel(frob,i) = FpM_FpC_mul(ffproj, x, p);
2491
}
2492
u = col_ei(f,2); k = 2;
2493
deg1 = ffdegree(u, frob, p);
2494
while (deg1 < f)
2495
{
2496
k++; u2 = col_ei(f, k);
2497
deg2 = ffdegree(u2, frob, p);
2498
deg = ulcm(deg1,deg2);
2499
if (deg == deg1) continue;
2500
if (deg == deg2) { deg1 = deg2; u = u2; continue; }
2501
u = ZC_add(u, u2);
2502
while (ffdegree(u, frob, p) < deg) u = ZC_add(u, u2);
2503
deg1 = deg;
2504
}
2505
v = lift_to_zk(u,c,N);
2506
2507
mul = cgetg(f+1,t_MAT);
2508
gel(mul,1) = v; /* assume w_1 = 1 */
2509
for (i=2; i<=f; i++) gel(mul,i) = zk_ei_mul(nf,v,c[i]);
2510
}
2511
2512
/* Z_K/pr = Fp(v), mul = mul by v */
2513
mul = FpM_red(mul, p);
2514
mul = FpM_mul(ffproj, mul, p);
2515
2516
pow = get_powers(mul, p);
2517
T = RgV_to_RgX(FpM_deplin(pow, p), vT);
2518
nfproj = cgetg(f+1, t_MAT);
2519
for (i=1; i<=f; i++) gel(nfproj,i) = lift_to_zk(gel(pow,i), c, N);
2520
2521
setlg(pow, f+1);
2522
ffproj = FpM_mul(FpM_inv(pow, p), ffproj, p);
2523
2524
res = cgetg(LARGEMODPR, t_COL);
2525
gel(res,mpr_TAU) = tau;
2526
gel(res,mpr_FFP) = ffproj;
2527
gel(res,3) = pr;
2528
gel(res,4) = T;
2529
gel(res,mpr_NFP) = nfproj; return gerepilecopy(av, res);
2530
}
2531
2532
GEN
2533
nfmodprinit(GEN nf, GEN pr) { return modprinit(nf, pr, 0, -1); }
2534
GEN
2535
zkmodprinit(GEN nf, GEN pr) { return modprinit(nf, pr, 1, -1); }
2536
GEN
2537
nfmodprinit0(GEN nf, GEN pr, long v) { return modprinit(nf, pr, 0, v); }
2538
2539
/* x may be a modpr */
2540
static int
2541
ok_modpr(GEN x)
2542
{ return typ(x) == t_COL && lg(x) >= SMALLMODPR && lg(x) <= LARGEMODPR; }
2543
void
2544
checkmodpr(GEN x)
2545
{
2546
if (!ok_modpr(x)) pari_err_TYPE("checkmodpr [use nfmodprinit]", x);
2547
checkprid(modpr_get_pr(x));
2548
}
2549
GEN
2550
get_modpr(GEN x)
2551
{ return ok_modpr(x)? x: NULL; }
2552
2553
int
2554
checkprid_i(GEN x)
2555
{
2556
return (typ(x) == t_VEC && lg(x) == 6
2557
&& typ(gel(x,2)) == t_COL && typ(gel(x,3)) == t_INT
2558
&& typ(gel(x,5)) != t_COL); /* tau changed to t_MAT/t_INT in 2.6 */
2559
}
2560
void
2561
checkprid(GEN x)
2562
{ if (!checkprid_i(x)) pari_err_TYPE("checkprid",x); }
2563
GEN
2564
get_prid(GEN x)
2565
{
2566
long lx = lg(x);
2567
if (lx == 3 && typ(x) == t_VEC) x = gel(x,1);
2568
if (checkprid_i(x)) return x;
2569
if (ok_modpr(x)) {
2570
x = modpr_get_pr(x);
2571
if (checkprid_i(x)) return x;
2572
}
2573
return NULL;
2574
}
2575
2576
static GEN
2577
to_ff_init(GEN nf, GEN *pr, GEN *T, GEN *p, int zk)
2578
{
2579
GEN modpr = ok_modpr(*pr)? *pr: modprinit(nf, *pr, zk, -1);
2580
*T = modpr_get_T(modpr);
2581
*pr = modpr_get_pr(modpr);
2582
*p = pr_get_p(*pr); return modpr;
2583
}
2584
2585
/* Return an element of O_K which is set to x Mod T */
2586
GEN
2587
modpr_genFq(GEN modpr)
2588
{
2589
switch(lg(modpr))
2590
{
2591
case SMALLMODPR: /* Fp */
2592
return gen_1;
2593
case LARGEMODPR: /* painful case, p \mid index */
2594
return gmael(modpr,mpr_NFP, 2);
2595
default: /* trivial case : p \nmid index */
2596
{
2597
long v = varn( modpr_get_T(modpr) );
2598
return pol_x(v);
2599
}
2600
}
2601
}
2602
2603
GEN
2604
nf_to_Fq_init(GEN nf, GEN *pr, GEN *T, GEN *p) {
2605
GEN modpr = to_ff_init(nf,pr,T,p,0);
2606
GEN tau = modpr_TAU(modpr);
2607
if (!tau) gel(modpr,mpr_TAU) = anti_uniformizer(nf, *pr);
2608
return modpr;
2609
}
2610
GEN
2611
zk_to_Fq_init(GEN nf, GEN *pr, GEN *T, GEN *p) {
2612
return to_ff_init(nf,pr,T,p,1);
2613
}
2614
2615
/* assume x in 'basis' form (t_COL) */
2616
GEN
2617
zk_to_Fq(GEN x, GEN modpr)
2618
{
2619
GEN pr = modpr_get_pr(modpr), p = pr_get_p(pr);
2620
GEN ffproj = gel(modpr,mpr_FFP);
2621
GEN T = modpr_get_T(modpr);
2622
return T? FpM_FpC_mul_FpX(ffproj,x, p, varn(T)): FpV_dotproduct(ffproj,x, p);
2623
}
2624
2625
/* REDUCTION Modulo a prime ideal */
2626
2627
/* nf a true nf */
2628
static GEN
2629
Rg_to_ff(GEN nf, GEN x0, GEN modpr)
2630
{
2631
GEN x = x0, den, pr = modpr_get_pr(modpr), p = pr_get_p(pr);
2632
long tx = typ(x);
2633
2634
if (tx == t_POLMOD) { x = gel(x,2); tx = typ(x); }
2635
switch(tx)
2636
{
2637
case t_INT: return modii(x, p);
2638
case t_FRAC: return Rg_to_Fp(x, p);
2639
case t_POL:
2640
switch(lg(x))
2641
{
2642
case 2: return gen_0;
2643
case 3: return Rg_to_Fp(gel(x,2), p);
2644
}
2645
x = Q_remove_denom(x, &den);
2646
x = poltobasis(nf, x);
2647
/* content(x) and den may not be coprime */
2648
break;
2649
case t_COL:
2650
x = Q_remove_denom(x, &den);
2651
/* content(x) and den are coprime */
2652
if (lg(x)-1 == nf_get_degree(nf)) break;
2653
default: pari_err_TYPE("Rg_to_ff",x);
2654
return NULL;/*LCOV_EXCL_LINE*/
2655
}
2656
if (den)
2657
{
2658
long v = Z_pvalrem(den, p, &den);
2659
if (v)
2660
{
2661
if (tx == t_POL) v -= ZV_pvalrem(x, p, &x);
2662
/* now v = valuation(true denominator of x) */
2663
if (v > 0)
2664
{
2665
GEN tau = modpr_TAU(modpr);
2666
if (!tau) pari_err_TYPE("zk_to_ff", x0);
2667
x = nfmuli(nf,x, nfpow_u(nf, tau, v));
2668
v -= ZV_pvalrem(x, p, &x);
2669
}
2670
if (v > 0) pari_err_INV("Rg_to_ff", mkintmod(gen_0,p));
2671
if (v) return gen_0;
2672
if (is_pm1(den)) den = NULL;
2673
}
2674
x = FpC_red(x, p);
2675
}
2676
x = zk_to_Fq(x, modpr);
2677
if (den)
2678
{
2679
GEN c = Fp_inv(den, p);
2680
x = typ(x) == t_INT? Fp_mul(x,c,p): FpX_Fp_mul(x,c,p);
2681
}
2682
return x;
2683
}
2684
2685
GEN
2686
nfreducemodpr(GEN nf, GEN x, GEN modpr)
2687
{
2688
pari_sp av = avma;
2689
nf = checknf(nf); checkmodpr(modpr);
2690
return gerepileupto(av, algtobasis(nf, Fq_to_nf(Rg_to_ff(nf,x,modpr),modpr)));
2691
}
2692
2693
GEN
2694
nfmodpr(GEN nf, GEN x, GEN pr)
2695
{
2696
pari_sp av = avma;
2697
GEN T, p, modpr;
2698
nf = checknf(nf);
2699
modpr = nf_to_Fq_init(nf, &pr, &T, &p);
2700
if (typ(x) == t_MAT && lg(x) == 3)
2701
{
2702
GEN y, v = famat_nfvalrem(nf, x, pr, &y);
2703
long s = signe(v);
2704
if (s < 0) pari_err_INV("Rg_to_ff", mkintmod(gen_0,p));
2705
if (s > 0) return gc_const(av, gen_0);
2706
x = FqV_factorback(nfV_to_FqV(gel(y,1), nf, modpr), gel(y,2), T, p);
2707
return gerepileupto(av, x);
2708
}
2709
x = Rg_to_ff(nf, x, modpr);
2710
x = Fq_to_FF(x, Tp_to_FF(T,p));
2711
return gerepilecopy(av, x);
2712
}
2713
GEN
2714
nfmodprlift(GEN nf, GEN x, GEN pr)
2715
{
2716
pari_sp av = avma;
2717
GEN y, T, p, modpr;
2718
long i, l, d;
2719
nf = checknf(nf);
2720
switch(typ(x))
2721
{
2722
case t_INT: return icopy(x);
2723
case t_FFELT: break;
2724
case t_VEC: case t_COL: case t_MAT:
2725
y = cgetg_copy(x,&l);
2726
for (i = 1; i < l; i++) gel(y,i) = nfmodprlift(nf,gel(x,i),pr);
2727
return y;
2728
default: pari_err_TYPE("nfmodprlit",x);
2729
}
2730
x = FF_to_FpXQ_i(x);
2731
d = degpol(x);
2732
if (d <= 0) { set_avma(av); return d? gen_0: icopy(gel(x,2)); }
2733
modpr = nf_to_Fq_init(nf, &pr, &T, &p);
2734
return gerepilecopy(av, Fq_to_nf(x, modpr));
2735
}
2736
2737
/* lift A from residue field to nf */
2738
GEN
2739
Fq_to_nf(GEN A, GEN modpr)
2740
{
2741
long dA;
2742
if (typ(A) == t_INT || lg(modpr) < LARGEMODPR) return A;
2743
dA = degpol(A);
2744
if (dA <= 0) return dA ? gen_0: gel(A,2);
2745
return ZM_ZX_mul(gel(modpr,mpr_NFP), A);
2746
}
2747
GEN
2748
FqV_to_nfV(GEN x, GEN modpr)
2749
{ pari_APPLY_same(Fq_to_nf(gel(x,i), modpr)) }
2750
GEN
2751
FqM_to_nfM(GEN A, GEN modpr)
2752
{
2753
long i,j,h,l = lg(A);
2754
GEN B = cgetg(l, t_MAT);
2755
2756
if (l == 1) return B;
2757
h = lgcols(A);
2758
for (j=1; j<l; j++)
2759
{
2760
GEN Aj = gel(A,j), Bj = cgetg(h,t_COL); gel(B,j) = Bj;
2761
for (i=1; i<h; i++) gel(Bj,i) = Fq_to_nf(gel(Aj,i), modpr);
2762
}
2763
return B;
2764
}
2765
GEN
2766
FqX_to_nfX(GEN A, GEN modpr)
2767
{
2768
long i, l;
2769
GEN B;
2770
2771
if (typ(A)!=t_POL) return icopy(A); /* scalar */
2772
B = cgetg_copy(A, &l); B[1] = A[1];
2773
for (i=2; i<l; i++) gel(B,i) = Fq_to_nf(gel(A,i), modpr);
2774
return B;
2775
}
2776
2777
/* reduce A to residue field */
2778
GEN
2779
nf_to_Fq(GEN nf, GEN A, GEN modpr)
2780
{
2781
pari_sp av = avma;
2782
return gerepileupto(av, Rg_to_ff(checknf(nf), A, modpr));
2783
}
2784
/* A t_VEC/t_COL */
2785
GEN
2786
nfV_to_FqV(GEN A, GEN nf,GEN modpr)
2787
{
2788
long i,l = lg(A);
2789
GEN B = cgetg(l,typ(A));
2790
for (i=1; i<l; i++) gel(B,i) = nf_to_Fq(nf,gel(A,i), modpr);
2791
return B;
2792
}
2793
/* A t_MAT */
2794
GEN
2795
nfM_to_FqM(GEN A, GEN nf,GEN modpr)
2796
{
2797
long i,j,h,l = lg(A);
2798
GEN B = cgetg(l,t_MAT);
2799
2800
if (l == 1) return B;
2801
h = lgcols(A);
2802
for (j=1; j<l; j++)
2803
{
2804
GEN Aj = gel(A,j), Bj = cgetg(h,t_COL); gel(B,j) = Bj;
2805
for (i=1; i<h; i++) gel(Bj,i) = nf_to_Fq(nf, gel(Aj,i), modpr);
2806
}
2807
return B;
2808
}
2809
/* A t_POL */
2810
GEN
2811
nfX_to_FqX(GEN A, GEN nf,GEN modpr)
2812
{
2813
long i,l = lg(A);
2814
GEN B = cgetg(l,t_POL); B[1] = A[1];
2815
for (i=2; i<l; i++) gel(B,i) = nf_to_Fq(nf,gel(A,i),modpr);
2816
return normalizepol_lg(B, l);
2817
}
2818
2819
/*******************************************************************/
2820
/* */
2821
/* RELATIVE ROUND 2 */
2822
/* */
2823
/*******************************************************************/
2824
/* Shallow functions */
2825
/* FIXME: use a bb_field and export the nfX_* routines */
2826
static GEN
2827
nfX_sub(GEN nf, GEN x, GEN y)
2828
{
2829
long i, lx = lg(x), ly = lg(y);
2830
GEN z;
2831
if (ly <= lx) {
2832
z = cgetg(lx,t_POL); z[1] = x[1];
2833
for (i=2; i < ly; i++) gel(z,i) = nfsub(nf,gel(x,i),gel(y,i));
2834
for ( ; i < lx; i++) gel(z,i) = gel(x,i);
2835
z = normalizepol_lg(z, lx);
2836
} else {
2837
z = cgetg(ly,t_POL); z[1] = y[1];
2838
for (i=2; i < lx; i++) gel(z,i) = nfsub(nf,gel(x,i),gel(y,i));
2839
for ( ; i < ly; i++) gel(z,i) = gneg(gel(y,i));
2840
z = normalizepol_lg(z, ly);
2841
}
2842
return z;
2843
}
2844
/* FIXME: quadratic multiplication */
2845
static GEN
2846
nfX_mul(GEN nf, GEN a, GEN b)
2847
{
2848
long da = degpol(a), db = degpol(b), dc, lc, k;
2849
GEN c;
2850
if (da < 0 || db < 0) return gen_0;
2851
dc = da + db;
2852
if (dc == 0) return nfmul(nf, gel(a,2),gel(b,2));
2853
lc = dc+3;
2854
c = cgetg(lc, t_POL); c[1] = a[1];
2855
for (k = 0; k <= dc; k++)
2856
{
2857
long i, I = minss(k, da);
2858
GEN d = NULL;
2859
for (i = maxss(k-db, 0); i <= I; i++)
2860
{
2861
GEN e = nfmul(nf, gel(a, i+2), gel(b, k-i+2));
2862
d = d? nfadd(nf, d, e): e;
2863
}
2864
gel(c, k+2) = d;
2865
}
2866
return normalizepol_lg(c, lc);
2867
}
2868
/* assume b monic */
2869
static GEN
2870
nfX_rem(GEN nf, GEN a, GEN b)
2871
{
2872
long da = degpol(a), db = degpol(b);
2873
if (da < 0) return gen_0;
2874
a = leafcopy(a);
2875
while (da >= db)
2876
{
2877
long i, k = da;
2878
GEN A = gel(a, k+2);
2879
for (i = db-1, k--; i >= 0; i--, k--)
2880
gel(a,k+2) = nfsub(nf, gel(a,k+2), nfmul(nf, A, gel(b,i+2)));
2881
a = normalizepol_lg(a, lg(a)-1);
2882
da = degpol(a);
2883
}
2884
return a;
2885
}
2886
static GEN
2887
nfXQ_mul(GEN nf, GEN a, GEN b, GEN T)
2888
{
2889
GEN c = nfX_mul(nf, a, b);
2890
if (typ(c) != t_POL) return c;
2891
return nfX_rem(nf, c, T);
2892
}
2893
2894
static void
2895
fill(long l, GEN H, GEN Hx, GEN I, GEN Ix)
2896
{
2897
long i;
2898
if (typ(Ix) == t_VEC) /* standard */
2899
for (i=1; i<l; i++) { gel(H,i) = gel(Hx,i); gel(I,i) = gel(Ix,i); }
2900
else /* constant ideal */
2901
for (i=1; i<l; i++) { gel(H,i) = gel(Hx,i); gel(I,i) = Ix; }
2902
}
2903
2904
/* given MODULES x and y by their pseudo-bases, returns a pseudo-basis of the
2905
* module generated by x and y. */
2906
static GEN
2907
rnfjoinmodules_i(GEN nf, GEN Hx, GEN Ix, GEN Hy, GEN Iy)
2908
{
2909
long lx = lg(Hx), ly = lg(Hy), l = lx+ly-1;
2910
GEN H = cgetg(l, t_MAT), I = cgetg(l, t_VEC);
2911
fill(lx, H , Hx, I , Ix);
2912
fill(ly, H+lx-1, Hy, I+lx-1, Iy); return nfhnf(nf, mkvec2(H, I));
2913
}
2914
static GEN
2915
rnfjoinmodules(GEN nf, GEN x, GEN y)
2916
{
2917
if (!x) return y;
2918
if (!y) return x;
2919
return rnfjoinmodules_i(nf, gel(x,1), gel(x,2), gel(y,1), gel(y,2));
2920
}
2921
2922
typedef struct {
2923
GEN multab, T,p;
2924
long h;
2925
} rnfeltmod_muldata;
2926
2927
static GEN
2928
_sqr(void *data, GEN x)
2929
{
2930
rnfeltmod_muldata *D = (rnfeltmod_muldata *) data;
2931
GEN z = x? tablesqr(D->multab,x)
2932
: tablemul_ei_ej(D->multab,D->h,D->h);
2933
return FqV_red(z,D->T,D->p);
2934
}
2935
static GEN
2936
_msqr(void *data, GEN x)
2937
{
2938
GEN x2 = _sqr(data, x), z;
2939
rnfeltmod_muldata *D = (rnfeltmod_muldata *) data;
2940
z = tablemul_ei(D->multab, x2, D->h);
2941
return FqV_red(z,D->T,D->p);
2942
}
2943
2944
/* Compute W[h]^n mod (T,p) in the extension, assume n >= 0. T a ZX */
2945
static GEN
2946
rnfeltid_powmod(GEN multab, long h, GEN n, GEN T, GEN p)
2947
{
2948
pari_sp av = avma;
2949
GEN y;
2950
rnfeltmod_muldata D;
2951
2952
if (!signe(n)) return gen_1;
2953
2954
D.multab = multab;
2955
D.h = h;
2956
D.T = T;
2957
D.p = p;
2958
y = gen_pow_fold(NULL, n, (void*)&D, &_sqr, &_msqr);
2959
return gerepilecopy(av, y);
2960
}
2961
2962
/* P != 0 has at most degpol(P) roots. Look for an element in Fq which is not
2963
* a root, cf repres() */
2964
static GEN
2965
FqX_non_root(GEN P, GEN T, GEN p)
2966
{
2967
long dP = degpol(P), f, vT;
2968
long i, j, k, pi, pp;
2969
GEN v;
2970
2971
if (dP == 0) return gen_1;
2972
pp = is_bigint(p) ? dP+1: itos(p);
2973
v = cgetg(dP + 2, t_VEC);
2974
gel(v,1) = gen_0;
2975
if (T)
2976
{ f = degpol(T); vT = varn(T); }
2977
else
2978
{ f = 1; vT = 0; }
2979
for (i=pi=1; i<=f; i++,pi*=pp)
2980
{
2981
GEN gi = i == 1? gen_1: pol_xn(i-1, vT), jgi = gi;
2982
for (j=1; j<pp; j++)
2983
{
2984
for (k=1; k<=pi; k++)
2985
{
2986
GEN z = Fq_add(gel(v,k), jgi, T,p);
2987
if (!gequal0(FqX_eval(P, z, T,p))) return z;
2988
gel(v, j*pi+k) = z;
2989
}
2990
if (j < pp-1) jgi = Fq_add(jgi, gi, T,p); /* j*g[i] */
2991
}
2992
}
2993
return NULL;
2994
}
2995
2996
/* Relative Dedekind criterion over (true) nf, applied to the order defined by a
2997
* root of monic irreducible polynomial P, modulo the prime ideal pr. Assume
2998
* vdisc = v_pr( disc(P) ).
2999
* Return NULL if nf[X]/P is pr-maximal. Otherwise, return [flag, O, v]:
3000
* O = enlarged order, given by a pseudo-basis
3001
* flag = 1 if O is proven pr-maximal (may be 0 and O nevertheless pr-maximal)
3002
* v = v_pr(disc(O)). */
3003
static GEN
3004
rnfdedekind_i(GEN nf, GEN P, GEN pr, long vdisc, long only_maximal)
3005
{
3006
GEN Ppr, A, I, p, tau, g, h, k, base, T, gzk, hzk, prinvp, pal, nfT, modpr;
3007
long m, vt, r, d, i, j, mpr;
3008
3009
if (vdisc < 0) pari_err_TYPE("rnfdedekind [non integral pol]", P);
3010
if (vdisc == 1) return NULL; /* pr-maximal */
3011
if (!only_maximal && !gequal1(leading_coeff(P)))
3012
pari_err_IMPL( "the full Dedekind criterion in the nonmonic case");
3013
/* either monic OR only_maximal = 1 */
3014
m = degpol(P);
3015
nfT = nf_get_pol(nf);
3016
modpr = nf_to_Fq_init(nf,&pr, &T, &p);
3017
Ppr = nfX_to_FqX(P, nf, modpr);
3018
mpr = degpol(Ppr);
3019
if (mpr < m) /* nonmonic => only_maximal = 1 */
3020
{
3021
if (mpr < 0) return NULL;
3022
if (! RgX_valrem(Ppr, &Ppr))
3023
{ /* nonzero constant coefficient */
3024
Ppr = RgX_shift_shallow(RgX_recip_i(Ppr), m - mpr);
3025
P = RgX_recip_i(P);
3026
}
3027
else
3028
{
3029
GEN z = FqX_non_root(Ppr, T, p);
3030
if (!z) pari_err_IMPL( "Dedekind in the difficult case");
3031
z = Fq_to_nf(z, modpr);
3032
if (typ(z) == t_INT)
3033
P = RgX_translate(P, z);
3034
else
3035
P = RgXQX_translate(P, z, T);
3036
P = RgX_recip_i(P);
3037
Ppr = nfX_to_FqX(P, nf, modpr); /* degpol(P) = degpol(Ppr) = m */
3038
}
3039
}
3040
A = gel(FqX_factor(Ppr,T,p),1);
3041
r = lg(A); /* > 1 */
3042
g = gel(A,1);
3043
for (i=2; i<r; i++) g = FqX_mul(g, gel(A,i), T, p);
3044
h = FqX_div(Ppr,g, T, p);
3045
gzk = FqX_to_nfX(g, modpr);
3046
hzk = FqX_to_nfX(h, modpr);
3047
k = nfX_sub(nf, P, nfX_mul(nf, gzk,hzk));
3048
tau = pr_get_tau(pr);
3049
switch(typ(tau))
3050
{
3051
case t_INT: k = gdiv(k, p); break;
3052
case t_MAT: k = RgX_Rg_div(tablemulvec(NULL,tau, k), p); break;
3053
}
3054
k = nfX_to_FqX(k, nf, modpr);
3055
k = FqX_normalize(FqX_gcd(FqX_gcd(g,h, T,p), k, T,p), T,p);
3056
d = degpol(k); /* <= m */
3057
if (!d) return NULL; /* pr-maximal */
3058
if (only_maximal) return gen_0; /* not maximal */
3059
3060
A = cgetg(m+d+1,t_MAT);
3061
I = cgetg(m+d+1,t_VEC); base = mkvec2(A, I);
3062
/* base[2] temporarily multiplied by p, for the final nfhnfmod,
3063
* which requires integral ideals */
3064
prinvp = pr_inv_p(pr); /* again multiplied by p */
3065
for (j=1; j<=m; j++)
3066
{
3067
gel(A,j) = col_ei(m, j);
3068
gel(I,j) = p;
3069
}
3070
pal = FqX_to_nfX(FqX_div(Ppr,k, T,p), modpr);
3071
for ( ; j<=m+d; j++)
3072
{
3073
gel(A,j) = RgX_to_RgC(pal,m);
3074
gel(I,j) = prinvp;
3075
if (j < m+d) pal = RgXQX_rem(RgX_shift_shallow(pal,1),P,nfT);
3076
}
3077
/* the modulus is integral */
3078
base = nfhnfmod(nf,base, idealmulpowprime(nf, powiu(p,m), pr, utoineg(d)));
3079
gel(base,2) = gdiv(gel(base,2), p); /* cancel the factor p */
3080
vt = vdisc - 2*d;
3081
return mkvec3(vt < 2? gen_1: gen_0, base, stoi(vt));
3082
}
3083
3084
/* [L:K] = n */
3085
static GEN
3086
triv_order(long n)
3087
{
3088
GEN z = cgetg(3, t_VEC);
3089
gel(z,1) = matid(n);
3090
gel(z,2) = const_vec(n, gen_1); return z;
3091
}
3092
3093
/* if flag is set, return gen_1 (resp. gen_0) if the order K[X]/(P)
3094
* is pr-maximal (resp. not pr-maximal). */
3095
GEN
3096
rnfdedekind(GEN nf, GEN P, GEN pr, long flag)
3097
{
3098
pari_sp av = avma;
3099
GEN z, dP;
3100
long v;
3101
3102
nf = checknf(nf);
3103
P = RgX_nffix("rnfdedekind", nf_get_pol(nf), P, 1);
3104
dP = nfX_disc(nf, P);
3105
if (!pr)
3106
{
3107
GEN fa = idealfactor(nf, dP);
3108
GEN Q = gel(fa,1), E = gel(fa,2);
3109
pari_sp av2 = avma;
3110
long i, l = lg(Q);
3111
for (i = 1; i < l; i++, set_avma(av2))
3112
{
3113
v = itos(gel(E,i));
3114
if (rnfdedekind_i(nf,P,gel(Q,i),v,1)) { set_avma(av); return gen_0; }
3115
set_avma(av2);
3116
}
3117
set_avma(av); return gen_1;
3118
}
3119
else if (typ(pr) == t_VEC)
3120
{ /* flag = 1 is implicit */
3121
if (lg(pr) == 1) { set_avma(av); return gen_1; }
3122
if (typ(gel(pr,1)) == t_VEC)
3123
{ /* list of primes */
3124
GEN Q = pr;
3125
pari_sp av2 = avma;
3126
long i, l = lg(Q);
3127
for (i = 1; i < l; i++, set_avma(av2))
3128
{
3129
v = nfval(nf, dP, gel(Q,i));
3130
if (rnfdedekind_i(nf,P,gel(Q,i),v,1)) { set_avma(av); return gen_0; }
3131
}
3132
set_avma(av); return gen_1;
3133
}
3134
}
3135
/* single prime */
3136
v = nfval(nf, dP, pr);
3137
z = rnfdedekind_i(nf, P, pr, v, flag);
3138
if (z)
3139
{
3140
if (flag) { set_avma(av); return gen_0; }
3141
z = gerepilecopy(av, z);
3142
}
3143
else
3144
{
3145
set_avma(av); if (flag) return gen_1;
3146
z = cgetg(4, t_VEC);
3147
gel(z,1) = gen_1;
3148
gel(z,2) = triv_order(degpol(P));
3149
gel(z,3) = stoi(v);
3150
}
3151
return z;
3152
}
3153
3154
static int
3155
ideal_is1(GEN x) {
3156
switch(typ(x))
3157
{
3158
case t_INT: return is_pm1(x);
3159
case t_MAT: return RgM_isidentity(x);
3160
}
3161
return 0;
3162
}
3163
3164
/* return a in ideal A such that v_pr(a) = v_pr(A) */
3165
static GEN
3166
minval(GEN nf, GEN A, GEN pr)
3167
{
3168
GEN ab = idealtwoelt(nf,A), a = gel(ab,1), b = gel(ab,2);
3169
if (nfval(nf,a,pr) > nfval(nf,b,pr)) a = b;
3170
return a;
3171
}
3172
3173
/* nf a true nf. Return NULL if power order if pr-maximal */
3174
static GEN
3175
rnfmaxord(GEN nf, GEN pol, GEN pr, long vdisc)
3176
{
3177
pari_sp av = avma, av1;
3178
long i, j, k, n, nn, vpol, cnt, sep;
3179
GEN q, q1, p, T, modpr, W, I, p1;
3180
GEN prhinv, mpi, Id;
3181
3182
if (DEBUGLEVEL>1) err_printf(" treating %Ps^%ld\n", pr, vdisc);
3183
modpr = nf_to_Fq_init(nf,&pr,&T,&p);
3184
av1 = avma;
3185
p1 = rnfdedekind_i(nf, pol, modpr, vdisc, 0);
3186
if (!p1) return gc_NULL(av);
3187
if (is_pm1(gel(p1,1))) return gerepilecopy(av,gel(p1,2));
3188
sep = itos(gel(p1,3));
3189
W = gmael(p1,2,1);
3190
I = gmael(p1,2,2);
3191
gerepileall(av1, 2, &W, &I);
3192
3193
mpi = zk_multable(nf, pr_get_gen(pr));
3194
n = degpol(pol); nn = n*n;
3195
vpol = varn(pol);
3196
q1 = q = pr_norm(pr);
3197
while (abscmpiu(q1,n) < 0) q1 = mulii(q1,q);
3198
Id = matid(n);
3199
prhinv = pr_inv(pr);
3200
av1 = avma;
3201
for(cnt=1;; cnt++)
3202
{
3203
GEN I0 = leafcopy(I), W0 = leafcopy(W);
3204
GEN Wa, Winv, Ip, A, MW, MWmod, F, pseudo, C, G;
3205
GEN Tauinv = cgetg(n+1, t_VEC), Tau = cgetg(n+1, t_VEC);
3206
3207
if (DEBUGLEVEL>1) err_printf(" pass no %ld\n",cnt);
3208
for (j=1; j<=n; j++)
3209
{
3210
GEN tau, tauinv;
3211
if (ideal_is1(gel(I,j)))
3212
{
3213
gel(I,j) = gel(Tau,j) = gel(Tauinv,j) = gen_1;
3214
continue;
3215
}
3216
gel(Tau,j) = tau = minval(nf, gel(I,j), pr);
3217
gel(Tauinv,j) = tauinv = nfinv(nf, tau);
3218
gel(W,j) = nfC_nf_mul(nf, gel(W,j), tau);
3219
gel(I,j) = idealmul(nf, tauinv, gel(I,j)); /* v_pr(I[j]) = 0 */
3220
}
3221
/* W = (Z_K/pr)-basis of O/pr. O = (W0,I0) ~ (W, I) */
3222
3223
/* compute MW: W_i*W_j = sum MW_k,(i,j) W_k */
3224
Wa = RgM_to_RgXV(W,vpol);
3225
Winv = nfM_inv(nf, W);
3226
MW = cgetg(nn+1, t_MAT);
3227
/* W_1 = 1 */
3228
for (j=1; j<=n; j++) gel(MW, j) = gel(MW, (j-1)*n+1) = gel(Id,j);
3229
for (i=2; i<=n; i++)
3230
for (j=i; j<=n; j++)
3231
{
3232
GEN z = nfXQ_mul(nf, gel(Wa,i), gel(Wa,j), pol);
3233
if (typ(z) != t_POL)
3234
z = nfC_nf_mul(nf, gel(Winv,1), z);
3235
else
3236
{
3237
z = RgX_to_RgC(z, lg(Winv)-1);
3238
z = nfM_nfC_mul(nf, Winv, z);
3239
}
3240
gel(MW, (i-1)*n+j) = gel(MW, (j-1)*n+i) = z;
3241
}
3242
3243
/* compute Ip = pr-radical [ could use Ker(trace) if q large ] */
3244
MWmod = nfM_to_FqM(MW,nf,modpr);
3245
F = cgetg(n+1, t_MAT); gel(F,1) = gel(Id,1);
3246
for (j=2; j<=n; j++) gel(F,j) = rnfeltid_powmod(MWmod, j, q1, T,p);
3247
Ip = FqM_ker(F,T,p);
3248
if (lg(Ip) == 1) { W = W0; I = I0; break; }
3249
3250
/* Fill C: W_k A_j = sum_i C_(i,j),k A_i */
3251
A = FqM_to_nfM(FqM_suppl(Ip,T,p), modpr);
3252
for (j = lg(Ip); j<=n; j++) gel(A,j) = nfC_multable_mul(gel(A,j), mpi);
3253
MW = nfM_mul(nf, nfM_inv(nf,A), MW);
3254
C = cgetg(n+1, t_MAT);
3255
for (k=1; k<=n; k++)
3256
{
3257
GEN mek = vecslice(MW, (k-1)*n+1, k*n), Ck;
3258
gel(C,k) = Ck = cgetg(nn+1, t_COL);
3259
for (j=1; j<=n; j++)
3260
{
3261
GEN z = nfM_nfC_mul(nf, mek, gel(A,j));
3262
for (i=1; i<=n; i++) gel(Ck, (j-1)*n+i) = nf_to_Fq(nf,gel(z,i),modpr);
3263
}
3264
}
3265
G = FqM_to_nfM(FqM_ker(C,T,p), modpr);
3266
3267
pseudo = rnfjoinmodules_i(nf, G,prhinv, Id,I);
3268
/* express W in terms of the power basis */
3269
W = nfM_mul(nf, W, gel(pseudo,1));
3270
I = gel(pseudo,2);
3271
/* restore the HNF property W[i,i] = 1. NB: W upper triangular, with
3272
* W[i,i] = Tau[i] */
3273
for (j=1; j<=n; j++)
3274
if (gel(Tau,j) != gen_1)
3275
{
3276
gel(W,j) = nfC_nf_mul(nf, gel(W,j), gel(Tauinv,j));
3277
gel(I,j) = idealmul(nf, gel(Tau,j), gel(I,j));
3278
}
3279
if (DEBUGLEVEL>3) err_printf(" new order:\n%Ps\n%Ps\n", W, I);
3280
if (sep <= 3 || gequal(I,I0)) break;
3281
3282
if (gc_needed(av1,2))
3283
{
3284
if(DEBUGMEM>1) pari_warn(warnmem,"rnfmaxord");
3285
gerepileall(av1,2, &W,&I);
3286
}
3287
}
3288
return gerepilecopy(av, mkvec2(W, I));
3289
}
3290
3291
GEN
3292
Rg_nffix(const char *f, GEN T, GEN c, int lift)
3293
{
3294
switch(typ(c))
3295
{
3296
case t_INT: case t_FRAC: return c;
3297
case t_POL:
3298
if (lg(c) >= lg(T)) c = RgX_rem(c,T);
3299
break;
3300
case t_POLMOD:
3301
if (!RgX_equal_var(gel(c,1), T)) pari_err_MODULUS(f, gel(c,1),T);
3302
c = gel(c,2);
3303
switch(typ(c))
3304
{
3305
case t_POL: break;
3306
case t_INT: case t_FRAC: return c;
3307
default: pari_err_TYPE(f, c);
3308
}
3309
break;
3310
default: pari_err_TYPE(f,c);
3311
}
3312
/* typ(c) = t_POL */
3313
if (varn(c) != varn(T)) pari_err_VAR(f, c,T);
3314
switch(lg(c))
3315
{
3316
case 2: return gen_0;
3317
case 3:
3318
c = gel(c,2); if (is_rational_t(typ(c))) return c;
3319
pari_err_TYPE(f,c);
3320
}
3321
RgX_check_QX(c, f);
3322
return lift? c: mkpolmod(c, T);
3323
}
3324
/* check whether P is a polynomials with coeffs in number field Q[y]/(T) */
3325
GEN
3326
RgX_nffix(const char *f, GEN T, GEN P, int lift)
3327
{
3328
long i, l, vT = varn(T);
3329
GEN Q = cgetg_copy(P, &l);
3330
if (typ(P) != t_POL) pari_err_TYPE(stack_strcat(f," [t_POL expected]"), P);
3331
if (varncmp(varn(P), vT) >= 0) pari_err_PRIORITY(f, P, ">=", vT);
3332
Q[1] = P[1];
3333
for (i=2; i<l; i++) gel(Q,i) = Rg_nffix(f, T, gel(P,i), lift);
3334
return normalizepol_lg(Q, l);
3335
}
3336
GEN
3337
RgV_nffix(const char *f, GEN T, GEN P, int lift)
3338
{
3339
long i, l;
3340
GEN Q = cgetg_copy(P, &l);
3341
for (i=1; i<l; i++) gel(Q,i) = Rg_nffix(f, T, gel(P,i), lift);
3342
return Q;
3343
}
3344
3345
static GEN
3346
get_d(GEN nf, GEN d)
3347
{
3348
GEN b = idealredmodpower(nf, d, 2, 100000);
3349
return nfmul(nf, d, nfsqr(nf,b));
3350
}
3351
3352
static GEN
3353
pr_factorback(GEN nf, GEN fa)
3354
{
3355
GEN P = gel(fa,1), E = gel(fa,2), z = gen_1;
3356
long i, l = lg(P);
3357
for (i = 1; i < l; i++) z = idealmulpowprime(nf, z, gel(P,i), gel(E,i));
3358
return z;
3359
}
3360
static GEN
3361
pr_factorback_scal(GEN nf, GEN fa)
3362
{
3363
GEN D = pr_factorback(nf,fa);
3364
if (typ(D) == t_MAT && RgM_isscalar(D,NULL)) D = gcoeff(D,1,1);
3365
return D;
3366
}
3367
3368
/* nf = base field K
3369
* pol= monic polynomial in Z_K[X] defining a relative extension L = K[X]/(pol).
3370
* Returns a pseudo-basis [A,I] of Z_L, set *pD to [D,d] and *pf to the
3371
* index-ideal; rnf is used when lim != 0 and may be NULL */
3372
GEN
3373
rnfallbase(GEN nf, GEN pol, GEN lim, GEN rnf, GEN *pD, GEN *pf, GEN *pDKP)
3374
{
3375
long i, j, jf, l;
3376
GEN fa, E, P, Ef, Pf, z, disc;
3377
3378
nf = checknf(nf); pol = liftpol_shallow(pol);
3379
if (!gequal1(leading_coeff(pol)))
3380
pari_err_IMPL("nonmonic relative polynomials in rnfallbase");
3381
disc = nf_to_scalar_or_basis(nf, nfX_disc(nf, pol));
3382
if (lim)
3383
{
3384
GEN rnfeq, zknf, dzknf, U, vU, dA, A, MB, dB, BdB, vj, B, Tabs;
3385
GEN D = idealhnf(nf, disc);
3386
long rU, m = nf_get_degree(nf), n = degpol(pol), N = n*m;
3387
nfmaxord_t S;
3388
3389
if (typ(lim) == t_INT)
3390
P = ZV_union_shallow(nf_get_ramified_primes(nf),
3391
gel(Z_factor_limit(gcoeff(D,1,1), itou(lim)), 1));
3392
else
3393
{
3394
P = cgetg_copy(lim, &l);
3395
for (i = 1; i < l; i++)
3396
{
3397
GEN p = gel(lim,i);
3398
if (typ(p) != t_INT) p = pr_get_p(p);
3399
gel(P,i) = p;
3400
}
3401
P = ZV_sort_uniq(P);
3402
}
3403
if (rnf)
3404
{
3405
rnfeq = rnf_get_map(rnf);
3406
zknf = rnf_get_nfzk(rnf);
3407
}
3408
else
3409
{
3410
rnfeq = nf_rnfeq(nf, pol);
3411
zknf = nf_nfzk(nf, rnfeq);
3412
}
3413
dzknf = gel(zknf,1);
3414
if (gequal1(dzknf)) dzknf = NULL;
3415
Tabs = gel(rnfeq,1);
3416
nfmaxord(&S, mkvec2(Tabs,P), 0);
3417
B = RgXV_unscale(S.basis, S.unscale);
3418
BdB = Q_remove_denom(B, &dB);
3419
MB = RgXV_to_RgM(BdB, N); /* HNF */
3420
3421
vU = cgetg(N+1, t_VEC);
3422
vj = cgetg(N+1, t_VECSMALL);
3423
gel(vU,1) = U = cgetg(m+1, t_MAT);
3424
gel(U,1) = col_ei(N, 1);
3425
A = dB? (dzknf? gdiv(dB,dzknf): dB): NULL;
3426
if (A && gequal1(A)) A = NULL;
3427
for (j = 2; j <= m; j++)
3428
{
3429
GEN t = gel(zknf,j);
3430
if (A) t = ZX_Z_mul(t, A);
3431
gel(U,j) = hnf_solve(MB, RgX_to_RgC(t, N));
3432
}
3433
for (i = 2; i <= N; i++)
3434
{
3435
GEN b = gel(BdB,i);
3436
gel(vU,i) = U = cgetg(m+1, t_MAT);
3437
gel(U,1) = hnf_solve(MB, RgX_to_RgC(b, N));
3438
for (j = 2; j <= m; j++)
3439
{
3440
GEN t = ZX_rem(ZX_mul(b, gel(zknf,j)), Tabs);
3441
if (dzknf) t = gdiv(t, dzknf);
3442
gel(U,j) = hnf_solve(MB, RgX_to_RgC(t, N));
3443
}
3444
}
3445
vj[1] = 1; U = gel(vU,1); rU = m;
3446
for (i = j = 2; i <= N; i++)
3447
{
3448
GEN V = shallowconcat(U, gel(vU,i));
3449
if (ZM_rank(V) != rU)
3450
{
3451
U = V; rU += m; vj[j++] = i;
3452
if (rU == N) break;
3453
}
3454
}
3455
if (dB) for(;;)
3456
{
3457
GEN c = gen_1, H = ZM_hnfmodid(U, dB);
3458
long ic = 0;
3459
for (i = 1; i <= N; i++)
3460
if (cmpii(gcoeff(H,i,i), c) > 0) { c = gcoeff(H,i,i); ic = i; }
3461
if (!ic) break;
3462
vj[j++] = ic;
3463
U = shallowconcat(H, gel(vU, ic));
3464
}
3465
setlg(vj, j);
3466
B = vecpermute(B, vj);
3467
3468
l = lg(B);
3469
A = cgetg(l,t_MAT);
3470
for (j = 1; j < l; j++)
3471
{
3472
GEN t = eltabstorel_lift(rnfeq, gel(B,j));
3473
gel(A,j) = Rg_to_RgC(t, n);
3474
}
3475
A = RgM_to_nfM(nf, A);
3476
A = Q_remove_denom(A, &dA);
3477
if (!dA)
3478
{ /* order is maximal */
3479
z = triv_order(n);
3480
if (pf) *pf = gen_1;
3481
}
3482
else
3483
{
3484
GEN fi;
3485
/* the first n columns of A are probably in HNF already */
3486
A = shallowconcat(vecslice(A,n+1,lg(A)-1), vecslice(A,1,n));
3487
A = mkvec2(A, const_vec(l-1,gen_1));
3488
if (DEBUGLEVEL > 2) err_printf("rnfallbase: nfhnf in dim %ld\n", l-1);
3489
z = nfhnfmod(nf, A, nfdetint(nf,A));
3490
gel(z,2) = gdiv(gel(z,2), dA);
3491
fi = idealprod(nf,gel(z,2));
3492
D = idealmul(nf, D, idealsqr(nf, fi));
3493
if (pf) *pf = idealinv(nf, fi);
3494
}
3495
if (RgM_isscalar(D,NULL)) D = gcoeff(D,1,1);
3496
if (pDKP) { settyp(S.dKP, t_VEC); *pDKP = S.dKP; }
3497
*pD = mkvec2(D, get_d(nf, disc)); return z;
3498
}
3499
fa = idealfactor(nf, disc);
3500
P = gel(fa,1); l = lg(P); z = NULL;
3501
E = gel(fa,2);
3502
Pf = cgetg(l, t_COL);
3503
Ef = cgetg(l, t_COL);
3504
for (i = j = jf = 1; i < l; i++)
3505
{
3506
GEN pr = gel(P,i);
3507
long e = itos(gel(E,i));
3508
if (e > 1)
3509
{
3510
GEN vD = rnfmaxord(nf, pol, pr, e);
3511
if (vD)
3512
{
3513
long ef = idealprodval(nf, gel(vD,2), pr);
3514
z = rnfjoinmodules(nf, z, vD);
3515
if (ef) { gel(Pf, jf) = pr; gel(Ef, jf++) = stoi(-ef); }
3516
e += 2 * ef;
3517
}
3518
}
3519
if (e) { gel(P, j) = pr; gel(E, j++) = stoi(e); }
3520
}
3521
setlg(P,j);
3522
setlg(E,j);
3523
if (pDKP)
3524
{
3525
GEN v = cgetg(j, t_VEC);
3526
for (i = 1; i < j; i++) gel(v,i) = pr_get_p(gel(P,i));
3527
*pDKP = ZV_sort_uniq(v);
3528
}
3529
if (pf)
3530
{
3531
setlg(Pf, jf);
3532
setlg(Ef, jf); *pf = pr_factorback_scal(nf, mkmat2(Pf,Ef));
3533
}
3534
*pD = mkvec2(pr_factorback_scal(nf,fa), get_d(nf, disc));
3535
return z? z: triv_order(degpol(pol));
3536
}
3537
3538
static GEN
3539
RgX_to_algX(GEN nf, GEN x)
3540
{
3541
long i, l;
3542
GEN y = cgetg_copy(x, &l); y[1] = x[1];
3543
for (i=2; i<l; i++) gel(y,i) = nf_to_scalar_or_alg(nf, gel(x,i));
3544
return y;
3545
}
3546
3547
GEN
3548
nfX_to_monic(GEN nf, GEN T, GEN *pL)
3549
{
3550
GEN lT, g, a;
3551
long i, l = lg(T);
3552
if (l == 2) return pol_0(varn(T));
3553
if (l == 3) return pol_1(varn(T));
3554
nf = checknf(nf);
3555
T = Q_primpart(RgX_to_nfX(nf, T));
3556
lT = leading_coeff(T); if (pL) *pL = lT;
3557
if (isint1(T)) return T;
3558
g = cgetg_copy(T, &l); g[1] = T[1]; a = lT;
3559
gel(g, l-1) = gen_1;
3560
gel(g, l-2) = gel(T,l-2);
3561
if (l == 4) return g;
3562
if (typ(lT) == t_INT)
3563
{
3564
gel(g, l-3) = gmul(a, gel(T,l-3));
3565
for (i = l-4; i > 1; i--) { a = mulii(a,lT); gel(g,i) = gmul(a, gel(T,i)); }
3566
}
3567
else
3568
{
3569
gel(g, l-3) = nfmul(nf, a, gel(T,l-3));
3570
for (i = l-3; i > 1; i--)
3571
{
3572
a = nfmul(nf,a,lT);
3573
gel(g,i) = nfmul(nf, a, gel(T,i));
3574
}
3575
}
3576
return RgX_to_algX(nf, g);
3577
}
3578
3579
GEN
3580
rnfdisc_factored(GEN nf, GEN pol, GEN *pd)
3581
{
3582
long i, j, l;
3583
GEN fa, E, P, disc, lim;
3584
3585
nf = checknf(nf);
3586
pol = rnfdisc_get_T(nf, pol, &lim);
3587
disc = nf_to_scalar_or_basis(nf, nfX_disc(nf, pol));
3588
pol = nfX_to_monic(nf, pol, NULL);
3589
fa = idealfactor_partial(nf, disc, lim);
3590
P = gel(fa,1); l = lg(P);
3591
E = gel(fa,2);
3592
for (i = j = 1; i < l; i++)
3593
{
3594
long e = itos(gel(E,i));
3595
GEN pr = gel(P,i);
3596
if (e > 1)
3597
{
3598
GEN vD = rnfmaxord(nf, pol, pr, e);
3599
if (vD) e += 2*idealprodval(nf, gel(vD,2), pr);
3600
}
3601
if (e) { gel(P, j) = pr; gel(E, j++) = stoi(e); }
3602
}
3603
if (pd) *pd = get_d(nf, disc);
3604
setlg(P, j);
3605
setlg(E, j); return fa;
3606
}
3607
GEN
3608
rnfdiscf(GEN nf, GEN pol)
3609
{
3610
pari_sp av = avma;
3611
GEN d, fa = rnfdisc_factored(nf, pol, &d);
3612
return gerepilecopy(av, mkvec2(pr_factorback_scal(nf,fa), d));
3613
}
3614
3615
GEN
3616
gen_if_principal(GEN bnf, GEN x)
3617
{
3618
pari_sp av = avma;
3619
GEN z = bnfisprincipal0(bnf,x, nf_GEN_IF_PRINCIPAL | nf_FORCE);
3620
return isintzero(z)? gc_NULL(av): z;
3621
}
3622
3623
static int
3624
is_pseudo_matrix(GEN O)
3625
{
3626
return (typ(O) ==t_VEC && lg(O) >= 3
3627
&& typ(gel(O,1)) == t_MAT
3628
&& typ(gel(O,2)) == t_VEC
3629
&& lgcols(O) == lg(gel(O,2)));
3630
}
3631
3632
/* given bnf and a pseudo-basis of an order in HNF [A,I], tries to simplify
3633
* the HNF as much as possible. The resulting matrix will be upper triangular
3634
* but the diagonal coefficients will not be equal to 1. The ideals are
3635
* guaranteed to be integral and primitive. */
3636
GEN
3637
rnfsimplifybasis(GEN bnf, GEN x)
3638
{
3639
pari_sp av = avma;
3640
long i, l;
3641
GEN y, Az, Iz, nf, A, I;
3642
3643
bnf = checkbnf(bnf); nf = bnf_get_nf(bnf);
3644
if (!is_pseudo_matrix(x)) pari_err_TYPE("rnfsimplifybasis",x);
3645
A = gel(x,1);
3646
I = gel(x,2); l = lg(I);
3647
y = cgetg(3, t_VEC);
3648
Az = cgetg(l, t_MAT); gel(y,1) = Az;
3649
Iz = cgetg(l, t_VEC); gel(y,2) = Iz;
3650
for (i = 1; i < l; i++)
3651
{
3652
GEN c, d;
3653
if (ideal_is1(gel(I,i))) {
3654
gel(Iz,i) = gen_1;
3655
gel(Az,i) = gel(A,i);
3656
continue;
3657
}
3658
3659
gel(Iz,i) = Q_primitive_part(gel(I,i), &c);
3660
gel(Az,i) = c? RgC_Rg_mul(gel(A,i),c): gel(A,i);
3661
if (c && ideal_is1(gel(Iz,i))) continue;
3662
3663
d = gen_if_principal(bnf, gel(Iz,i));
3664
if (d)
3665
{
3666
gel(Iz,i) = gen_1;
3667
gel(Az,i) = nfC_nf_mul(nf, gel(Az,i), d);
3668
}
3669
}
3670
return gerepilecopy(av, y);
3671
}
3672
3673
static GEN
3674
get_order(GEN nf, GEN O, const char *s)
3675
{
3676
if (typ(O) == t_POL)
3677
return rnfpseudobasis(nf, O);
3678
if (!is_pseudo_matrix(O)) pari_err_TYPE(s, O);
3679
return O;
3680
}
3681
3682
GEN
3683
rnfdet(GEN nf, GEN order)
3684
{
3685
pari_sp av = avma;
3686
GEN A, I, D;
3687
nf = checknf(nf);
3688
order = get_order(nf, order, "rnfdet");
3689
A = gel(order,1);
3690
I = gel(order,2);
3691
D = idealmul(nf, nfM_det(nf,A), idealprod(nf,I));
3692
return gerepileupto(av, D);
3693
}
3694
3695
/* Given two fractional ideals a and b, gives x in a, y in b, z in b^-1,
3696
t in a^-1 such that xt-yz=1. In the present version, z is in Z. */
3697
static void
3698
nfidealdet1(GEN nf, GEN a, GEN b, GEN *px, GEN *py, GEN *pz, GEN *pt)
3699
{
3700
GEN x, uv, y, da, db;
3701
3702
a = idealinv(nf,a);
3703
a = Q_remove_denom(a, &da);
3704
b = Q_remove_denom(b, &db);
3705
x = idealcoprime(nf,a,b);
3706
uv = idealaddtoone(nf, idealmul(nf,x,a), b);
3707
y = gel(uv,2);
3708
if (da) x = gmul(x,da);
3709
if (db) y = gdiv(y,db);
3710
*px = x;
3711
*py = y;
3712
*pz = db ? negi(db): gen_m1;
3713
*pt = nfdiv(nf, gel(uv,1), x);
3714
}
3715
3716
/* given a pseudo-basis of an order in HNF [A,I] (or [A,I,D,d]), gives an
3717
* n x n matrix (not in HNF) of a pseudo-basis and an ideal vector
3718
* [1,1,...,1,I] such that order = Z_K^(n-1) x I.
3719
* Uses the approximation theorem ==> slow. */
3720
GEN
3721
rnfsteinitz(GEN nf, GEN order)
3722
{
3723
pari_sp av = avma;
3724
long i, n, l;
3725
GEN A, I, p1;
3726
3727
nf = checknf(nf);
3728
order = get_order(nf, order, "rnfsteinitz");
3729
A = RgM_to_nfM(nf, gel(order,1));
3730
I = leafcopy(gel(order,2)); n=lg(A)-1;
3731
for (i=1; i<n; i++)
3732
{
3733
GEN c1, c2, b, a = gel(I,i);
3734
gel(I,i) = gen_1;
3735
if (ideal_is1(a)) continue;
3736
3737
c1 = gel(A,i);
3738
c2 = gel(A,i+1);
3739
b = gel(I,i+1);
3740
if (ideal_is1(b))
3741
{
3742
gel(A,i) = c2;
3743
gel(A,i+1) = gneg(c1);
3744
gel(I,i+1) = a;
3745
}
3746
else
3747
{
3748
pari_sp av2 = avma;
3749
GEN x, y, z, t;
3750
nfidealdet1(nf,a,b, &x,&y,&z,&t);
3751
x = RgC_add(nfC_nf_mul(nf, c1, x), nfC_nf_mul(nf, c2, y));
3752
y = RgC_add(nfC_nf_mul(nf, c1, z), nfC_nf_mul(nf, c2, t));
3753
gerepileall(av2, 2, &x,&y);
3754
gel(A,i) = x;
3755
gel(A,i+1) = y;
3756
gel(I,i+1) = Q_primitive_part(idealmul(nf,a,b), &p1);
3757
if (p1) gel(A,i+1) = nfC_nf_mul(nf, gel(A,i+1), p1);
3758
}
3759
}
3760
l = lg(order);
3761
p1 = cgetg(l,t_VEC);
3762
gel(p1,1) = A;
3763
gel(p1,2) = I; for (i=3; i<l; i++) gel(p1,i) = gel(order,i);
3764
return gerepilecopy(av, p1);
3765
}
3766
3767
/* Given bnf and either an order as output by rnfpseudobasis or a polynomial,
3768
* and outputs a basis if it is free, an n+1-generating set if it is not */
3769
GEN
3770
rnfbasis(GEN bnf, GEN order)
3771
{
3772
pari_sp av = avma;
3773
long j, n;
3774
GEN nf, A, I, cl, col, a;
3775
3776
bnf = checkbnf(bnf); nf = bnf_get_nf(bnf);
3777
order = get_order(nf, order, "rnfbasis");
3778
I = gel(order,2); n = lg(I)-1;
3779
j=1; while (j<n && ideal_is1(gel(I,j))) j++;
3780
if (j<n)
3781
{
3782
order = rnfsteinitz(nf,order);
3783
I = gel(order,2);
3784
}
3785
A = gel(order,1);
3786
col= gel(A,n); A = vecslice(A, 1, n-1);
3787
cl = gel(I,n);
3788
a = gen_if_principal(bnf, cl);
3789
if (!a)
3790
{
3791
GEN v = idealtwoelt(nf, cl);
3792
A = shallowconcat(A, gmul(gel(v,1), col));
3793
a = gel(v,2);
3794
}
3795
A = shallowconcat(A, nfC_nf_mul(nf, col, a));
3796
return gerepilecopy(av, A);
3797
}
3798
3799
/* Given bnf and either an order as output by rnfpseudobasis or a polynomial,
3800
* and outputs a basis (not pseudo) in Hermite Normal Form if it exists, zero
3801
* if not
3802
*/
3803
GEN
3804
rnfhnfbasis(GEN bnf, GEN order)
3805
{
3806
pari_sp av = avma;
3807
long j, n;
3808
GEN nf, A, I, a;
3809
3810
bnf = checkbnf(bnf); nf = bnf_get_nf(bnf);
3811
order = get_order(nf, order, "rnfbasis");
3812
A = gel(order,1); A = RgM_shallowcopy(A);
3813
I = gel(order,2); n = lg(A)-1;
3814
for (j=1; j<=n; j++)
3815
{
3816
if (ideal_is1(gel(I,j))) continue;
3817
a = gen_if_principal(bnf, gel(I,j));
3818
if (!a) { set_avma(av); return gen_0; }
3819
gel(A,j) = nfC_nf_mul(nf, gel(A,j), a);
3820
}
3821
return gerepilecopy(av,A);
3822
}
3823
3824
static long
3825
rnfisfree_aux(GEN bnf, GEN order)
3826
{
3827
long n, j;
3828
GEN nf, P, I;
3829
3830
bnf = checkbnf(bnf);
3831
if (is_pm1( bnf_get_no(bnf) )) return 1;
3832
nf = bnf_get_nf(bnf);
3833
order = get_order(nf, order, "rnfisfree");
3834
I = gel(order,2); n = lg(I)-1;
3835
j=1; while (j<=n && ideal_is1(gel(I,j))) j++;
3836
if (j>n) return 1;
3837
3838
P = gel(I,j);
3839
for (j++; j<=n; j++)
3840
if (!ideal_is1(gel(I,j))) P = idealmul(nf,P,gel(I,j));
3841
return gequal0( isprincipal(bnf,P) );
3842
}
3843
3844
long
3845
rnfisfree(GEN bnf, GEN order)
3846
{ pari_sp av = avma; return gc_long(av, rnfisfree_aux(bnf,order)); }
3847
3848
/**********************************************************************/
3849
/** **/
3850
/** COMPOSITUM OF TWO NUMBER FIELDS **/
3851
/** **/
3852
/**********************************************************************/
3853
static GEN
3854
compositum_fix(GEN nf, GEN A)
3855
{
3856
int ok;
3857
if (nf)
3858
{
3859
A = Q_primpart(liftpol_shallow(A)); RgX_check_ZXX(A,"polcompositum");
3860
ok = nfissquarefree(nf,A);
3861
}
3862
else
3863
{
3864
A = Q_primpart(A); RgX_check_ZX(A,"polcompositum");
3865
ok = ZX_is_squarefree(A);
3866
}
3867
if (!ok) pari_err_DOMAIN("polcompositum","issquarefree(arg)","=",gen_0,A);
3868
return A;
3869
}
3870
#define next_lambda(a) (a>0 ? -a : 1-a)
3871
3872
static long
3873
nfcompositum_lambda(GEN nf, GEN A, GEN B, long lambda)
3874
{
3875
pari_sp av = avma;
3876
forprime_t S;
3877
GEN T = nf_get_pol(nf);
3878
long vT = varn(T);
3879
ulong p;
3880
init_modular_big(&S);
3881
p = u_forprime_next(&S);
3882
while (1)
3883
{
3884
GEN Hp, Tp, a;
3885
if (DEBUGLEVEL>4) err_printf("Trying lambda = %ld\n", lambda);
3886
a = ZXX_to_FlxX(RgX_rescale(A, stoi(-lambda)), p, vT);
3887
Tp = ZX_to_Flx(T, p);
3888
Hp = FlxqX_direct_compositum(a, ZXX_to_FlxX(B, p, vT), Tp, p);
3889
if (!FlxqX_is_squarefree(Hp, Tp, p))
3890
{ lambda = next_lambda(lambda); continue; }
3891
if (DEBUGLEVEL>4) err_printf("Final lambda = %ld\n", lambda);
3892
return gc_long(av, lambda);
3893
}
3894
}
3895
3896
/* modular version */
3897
GEN
3898
nfcompositum(GEN nf, GEN A, GEN B, long flag)
3899
{
3900
pari_sp av = avma;
3901
int same;
3902
long v, k;
3903
GEN C, D, LPRS;
3904
3905
if (typ(A)!=t_POL) pari_err_TYPE("polcompositum",A);
3906
if (typ(B)!=t_POL) pari_err_TYPE("polcompositum",B);
3907
if (degpol(A)<=0 || degpol(B)<=0) pari_err_CONSTPOL("polcompositum");
3908
v = varn(A);
3909
if (varn(B) != v) pari_err_VAR("polcompositum", A,B);
3910
if (nf)
3911
{
3912
nf = checknf(nf);
3913
if (varncmp(v,nf_get_varn(nf))>=0) pari_err_PRIORITY("polcompositum", nf, ">=", v);
3914
}
3915
same = (A == B || RgX_equal(A,B));
3916
A = compositum_fix(nf,A);
3917
B = same ? A: compositum_fix(nf,B);
3918
3919
D = LPRS = NULL; /* -Wall */
3920
k = same? -1: 1;
3921
if (nf)
3922
{
3923
long v0 = fetch_var();
3924
GEN q, T = nf_get_pol(nf);
3925
A = liftpol_shallow(A);
3926
B = liftpol_shallow(B);
3927
k = nfcompositum_lambda(nf, A, B, k);
3928
if (flag&1)
3929
{
3930
GEN H0, H1;
3931
GEN chgvar = deg1pol_shallow(stoi(k),pol_x(v0),v);
3932
GEN B1 = poleval(QXQX_to_mod_shallow(B, T), chgvar);
3933
C = RgX_resultant_all(QXQX_to_mod_shallow(A, T), B1, &q);
3934
C = gsubst(C,v0,pol_x(v));
3935
C = lift_if_rational(C);
3936
H0 = gsubst(gel(q,2),v0,pol_x(v));
3937
H1 = gsubst(gel(q,3),v0,pol_x(v));
3938
if (typ(H0) != t_POL) H0 = scalarpol_shallow(H0,v);
3939
if (typ(H1) != t_POL) H1 = scalarpol_shallow(H1,v);
3940
H0 = lift_if_rational(H0);
3941
H1 = lift_if_rational(H1);
3942
LPRS = mkvec2(H0,H1);
3943
}
3944
else
3945
{
3946
C = nf_direct_compositum(nf, RgX_rescale(A,stoi(-k)), B);
3947
setvarn(C, v); C = QXQX_to_mod_shallow(C, T);
3948
}
3949
}
3950
else
3951
{
3952
B = leafcopy(B); setvarn(B,fetch_var_higher());
3953
C = (flag&1)? ZX_ZXY_resultant_all(A, B, &k, &LPRS)
3954
: ZX_compositum(A, B, &k);
3955
setvarn(C, v);
3956
}
3957
/* C = Res_Y (A(Y), B(X + kY)) guaranteed squarefree */
3958
if (flag & 2)
3959
C = mkvec(C);
3960
else
3961
{
3962
if (same)
3963
{
3964
D = RgX_rescale(A, stoi(1 - k));
3965
if (nf) D = QXQX_to_mod_shallow(D, nf_get_pol(nf));
3966
C = RgX_div(C, D);
3967
if (degpol(C) <= 0)
3968
C = mkvec(D);
3969
else
3970
C = shallowconcat(nf? gel(nffactor(nf,C),1): ZX_DDF(C), D);
3971
}
3972
else
3973
C = nf? gel(nffactor(nf,C),1): ZX_DDF(C);
3974
}
3975
gen_sort_inplace(C, (void*)(nf?&cmp_RgX: &cmpii), &gen_cmp_RgX, NULL);
3976
if (flag&1)
3977
{ /* a,b,c root of A,B,C = compositum, c = b - k a */
3978
long i, l = lg(C);
3979
GEN a, b, mH0 = RgX_neg(gel(LPRS,1)), H1 = gel(LPRS,2);
3980
setvarn(mH0,v);
3981
setvarn(H1,v);
3982
for (i=1; i<l; i++)
3983
{
3984
GEN D = gel(C,i);
3985
a = RgXQ_mul(mH0, nf? RgXQ_inv(H1,D): QXQ_inv(H1,D), D);
3986
b = gadd(pol_x(v), gmulsg(k,a));
3987
if (degpol(D) == 1) b = RgX_rem(b,D);
3988
gel(C,i) = mkvec4(D, mkpolmod(a,D), mkpolmod(b,D), stoi(-k));
3989
}
3990
}
3991
(void)delete_var();
3992
settyp(C, t_VEC);
3993
if (flag&2) C = gel(C,1);
3994
return gerepilecopy(av, C);
3995
}
3996
GEN
3997
polcompositum0(GEN A, GEN B, long flag)
3998
{ return nfcompositum(NULL,A,B,flag); }
3999
4000
GEN
4001
compositum(GEN pol1,GEN pol2) { return polcompositum0(pol1,pol2,0); }
4002
GEN
4003
compositum2(GEN pol1,GEN pol2) { return polcompositum0(pol1,pol2,1); }
4004
4005