Book a Demo!
CoCalc Logo Icon
StoreFeaturesDocsShareSupportNewsAboutPoliciesSign UpSign In
Download

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

28486 views
License: GPL3
ubuntu2004
1
/* Copyright (C) 2000 The PARI group.
2
3
This file is part of the PARI/GP package.
4
5
PARI/GP is free software; you can redistribute it and/or modify it under the
6
terms of the GNU General Public License as published by the Free Software
7
Foundation; either version 2 of the License, or (at your option) any later
8
version. It is distributed in the hope that it will be useful, but WITHOUT
9
ANY WARRANTY WHATSOEVER.
10
11
Check the License for details. You should have received a copy of it, along
12
with the package; see the file 'COPYING'. If not, write to the Free Software
13
Foundation, Inc., 51 Franklin Street, Fifth Floor, Boston, MA 02110-1301 USA. */
14
#include "pari.h"
15
#include "paripriv.h"
16
17
#define DEBUGLEVEL DEBUGLEVEL_hensel
18
19
/***********************************************************************/
20
/** **/
21
/** QUADRATIC HENSEL LIFT (adapted from V. Shoup's NTL) **/
22
/** **/
23
/***********************************************************************/
24
/* Setup for divide/conquer quadratic Hensel lift
25
* a = set of k t_POL in Z[X] = factors over Fp (T=NULL) or Fp[Y]/(T)
26
* V = set of products of factors built as follows
27
* 1) V[1..k] = initial a
28
* 2) iterate:
29
* append to V the two smallest factors (minimal degree) in a, remove them
30
* from a and replace them by their product [net loss for a = 1 factor]
31
*
32
* W = bezout coeffs W[i]V[i] + W[i+1]V[i+1] = 1
33
*
34
* link[i] = -j if V[i] = a[j]
35
* j if V[i] = V[j] * V[j+1]
36
* Arrays (link, V, W) pre-allocated for 2k - 2 elements */
37
static void
38
BuildTree(GEN link, GEN V, GEN W, GEN a, GEN T, GEN p)
39
{
40
long k = lg(a)-1;
41
long i, j, s, minp, mind;
42
43
for (i=1; i<=k; i++) { gel(V,i) = gel(a,i); link[i] = -i; }
44
45
for (j=1; j <= 2*k-5; j+=2,i++)
46
{
47
minp = j;
48
mind = degpol(gel(V,j));
49
for (s=j+1; s<i; s++)
50
if (degpol(gel(V,s)) < mind) { minp = s; mind = degpol(gel(V,s)); }
51
52
swap(gel(V,j), gel(V,minp));
53
lswap(link[j], link[minp]);
54
55
minp = j+1;
56
mind = degpol(gel(V,j+1));
57
for (s=j+2; s<i; s++)
58
if (degpol(gel(V,s)) < mind) { minp = s; mind = degpol(gel(V,s)); }
59
60
swap(gel(V,j+1), gel(V,minp));
61
lswap(link[j+1], link[minp]);
62
63
gel(V,i) = FqX_mul(gel(V,j), gel(V,j+1), T, p);
64
link[i] = j;
65
}
66
67
for (j=1; j <= 2*k-3; j+=2)
68
{
69
GEN d, u, v;
70
d = FqX_extgcd(gel(V,j), gel(V,j+1), T, p, &u, &v);
71
if (degpol(d) > 0) pari_err_COPRIME("BuildTree", gel(V,j), gel(V,j+1));
72
d = gel(d,2);
73
if (!gequal1(d))
74
{
75
if (typ(d)==t_POL)
76
{
77
d = FpXQ_inv(d, T, p);
78
u = FqX_Fq_mul(u, d, T, p);
79
v = FqX_Fq_mul(v, d, T, p);
80
}
81
else
82
{
83
d = Fp_inv(d, p);
84
u = FqX_Fp_mul(u, d, T,p);
85
v = FqX_Fp_mul(v, d, T,p);
86
}
87
}
88
gel(W,j) = u;
89
gel(W,j+1) = v;
90
}
91
}
92
93
/* au + bv = 1 (p0), ab = f (p0). Lift mod p1 = p0 pd (<= p0^2).
94
* If noinv is set, don't lift the inverses u and v */
95
static void
96
ZpX_HenselLift(GEN V, GEN W, long j, GEN f, GEN pd, GEN p0, GEN p1, int noinv)
97
{
98
pari_sp av = avma;
99
long space = lg(f) * lgefint(p1);
100
GEN a2, b2, g, z, s, t;
101
GEN a = gel(V,j), b = gel(V,j+1);
102
GEN u = gel(W,j), v = gel(W,j+1);
103
104
(void)new_chunk(space); /* HACK */
105
g = ZX_sub(f, ZX_mul(a,b));
106
g = ZX_Z_divexact(g, p0);
107
g = FpX_red(g, pd);
108
z = FpX_mul(v,g, pd);
109
t = FpX_divrem(z,a, pd, &s);
110
t = ZX_add(ZX_mul(u,g), ZX_mul(t,b));
111
t = FpX_red(t, pd);
112
t = ZX_Z_mul(t,p0);
113
s = ZX_Z_mul(s,p0);
114
set_avma(av);
115
a2 = ZX_add(a,s);
116
b2 = ZX_add(b,t);
117
118
/* already reduced mod p1 = pd p0 */
119
gel(V,j) = a2;
120
gel(V,j+1) = b2;
121
if (noinv) return;
122
123
av = avma;
124
(void)new_chunk(space); /* HACK */
125
g = ZX_add(ZX_mul(u,a2), ZX_mul(v,b2));
126
g = Z_ZX_sub(gen_1, g);
127
g = ZX_Z_divexact(g, p0);
128
g = FpX_red(g, pd);
129
z = FpX_mul(v,g, pd);
130
t = FpX_divrem(z,a, pd, &s);
131
t = ZX_add(ZX_mul(u,g), ZX_mul(t,b));
132
t = FpX_red(t, pd);
133
t = ZX_Z_mul(t,p0);
134
s = ZX_Z_mul(s,p0);
135
set_avma(av);
136
gel(W,j) = ZX_add(u,t);
137
gel(W,j+1) = ZX_add(v,s);
138
}
139
140
static void
141
ZpXQ_HenselLift(GEN V, GEN W, long j, GEN f, GEN Td, GEN T1, GEN pd, GEN p0, GEN p1, int noinv)
142
{
143
pari_sp av = avma;
144
const long n = degpol(T1), vT = varn(T1);
145
long space = lg(f) * lgefint(p1) * lg(T1);
146
GEN a2, b2, g, z, s, t;
147
GEN a = gel(V,j), b = gel(V,j+1);
148
GEN u = gel(W,j), v = gel(W,j+1);
149
150
(void)new_chunk(space); /* HACK */
151
g = RgX_sub(f, Kronecker_to_ZXX(ZXX_mul_Kronecker(a,b,n), n, vT));
152
g = FpXQX_red(g, T1, p1);
153
g = RgX_Rg_divexact(g, p0);
154
z = FpXQX_mul(v,g, Td,pd);
155
t = FpXQX_divrem(z,a, Td,pd, &s);
156
t = ZX_add(ZXX_mul_Kronecker(u,g,n), ZXX_mul_Kronecker(t,b,n));
157
t = Kronecker_to_ZXX(t, n, vT);
158
t = FpXQX_red(t, Td, pd);
159
t = RgX_Rg_mul(t,p0);
160
s = RgX_Rg_mul(s,p0);
161
set_avma(av);
162
163
a2 = RgX_add(a,s);
164
b2 = RgX_add(b,t);
165
/* already reduced mod p1 = pd p0 */
166
gel(V,j) = a2;
167
gel(V,j+1) = b2;
168
if (noinv) return;
169
170
av = avma;
171
(void)new_chunk(space); /* HACK */
172
g = ZX_add(ZXX_mul_Kronecker(u,a2,n), ZXX_mul_Kronecker(v,b2,n));
173
g = Kronecker_to_ZXX(g, n, vT);
174
g = Rg_RgX_sub(gen_1, g);
175
g = FpXQX_red(g, T1, p1);
176
g = RgX_Rg_divexact(g, p0);
177
z = FpXQX_mul(v,g, Td,pd);
178
t = FpXQX_divrem(z,a, Td,pd, &s);
179
t = ZX_add(ZXX_mul_Kronecker(u,g,n), ZXX_mul_Kronecker(t,b,n));
180
t = Kronecker_to_ZXX(t, n, vT);
181
t = FpXQX_red(t, Td, pd);
182
t = RgX_Rg_mul(t,p0);
183
s = RgX_Rg_mul(s,p0);
184
set_avma(av);
185
gel(W,j) = RgX_add(u,t);
186
gel(W,j+1) = RgX_add(v,s);
187
}
188
189
/* v list of factors, w list of inverses. f = v[j] v[j+1]
190
* Lift v[j] and v[j+1] mod p0 pd (possibly mod T), then all their divisors */
191
static void
192
ZpX_RecTreeLift(GEN link, GEN v, GEN w, GEN pd, GEN p0, GEN p1,
193
GEN f, long j, int noinv)
194
{
195
if (j < 0) return;
196
ZpX_HenselLift(v, w, j, f, pd, p0,p1, noinv);
197
ZpX_RecTreeLift(link, v, w, pd, p0,p1, gel(v,j) , link[j ], noinv);
198
ZpX_RecTreeLift(link, v, w, pd, p0,p1, gel(v,j+1), link[j+1], noinv);
199
}
200
static void
201
ZpXQ_RecTreeLift(GEN link, GEN v, GEN w, GEN Td, GEN T1, GEN pd, GEN p0, GEN p1,
202
GEN f, long j, int noinv)
203
{
204
if (j < 0) return;
205
ZpXQ_HenselLift(v, w, j, f, Td,T1, pd, p0,p1, noinv);
206
ZpXQ_RecTreeLift(link, v, w, Td,T1, pd, p0,p1, gel(v,j) , link[j ], noinv);
207
ZpXQ_RecTreeLift(link, v, w, Td,T1, pd, p0,p1, gel(v,j+1), link[j+1], noinv);
208
}
209
210
/* Assume n > 0. We want to go to accuracy n, starting from accuracy 1, using
211
* a quadratically convergent algorithm. Goal: 9 -> 1,2,3,5,9 instead of
212
* 1,2,4,8,9 (sequence of accuracies).
213
*
214
* Let a0 = 1, a1 = 2, a2, ... ak = n, the sequence of accuracies. To obtain
215
* it, work backwards:
216
* a(k) = n, a(i-1) = (a(i) + 1) \ 2,
217
* but we do not want to store a(i) explicitly, even as a t_VECSMALL, since
218
* this would leave an object on the stack. We store a(i) implicitly in a
219
* MASK: let a(0) = 1, if the i-bit of MASK is set, set a(i+1) = 2 a(i) - 1,
220
* and 2a(i) otherwise.
221
*
222
* In fact, we do something a little more complicated to simplify the
223
* function interface and avoid returning k and MASK separately: we return
224
* MASK + 2^(k+1), so the highest bit of the mask indicates the length of the
225
* sequence, and the following ones are as above. */
226
ulong
227
quadratic_prec_mask(long n)
228
{
229
long a = n, i;
230
ulong mask = 0;
231
for(i = 1;; i++, mask <<= 1)
232
{
233
mask |= (a&1); a = (a+1)>>1;
234
if (a==1) return mask | (1UL << i);
235
}
236
}
237
238
/* Lift to precision p^e0.
239
* a = modular factors of f mod (p,T) [possibly T=NULL]
240
* OR a TreeLift structure [e, link, v, w]: go on lifting
241
* flag = 0: standard.
242
* flag = 1: return TreeLift structure */
243
static GEN
244
MultiLift(GEN f, GEN a, GEN T, GEN p, long e0, long flag)
245
{
246
long i, eold, e, k = lg(a) - 1;
247
GEN E, v, w, link, penew, Tnew;
248
ulong mask;
249
pari_timer Ti;
250
251
if (k < 2) pari_err_DOMAIN("MultiLift", "#(modular factors)", "<", gen_2,a);
252
if (e0 < 1) pari_err_DOMAIN("MultiLift", "precision", "<", gen_1,stoi(e0));
253
if (e0 == 1) return a;
254
255
if (DEBUGLEVEL > 3) timer_start(&Ti);
256
if (typ(gel(a,1)) == t_INT)
257
{ /* a = TreeLift structure */
258
e = itos(gel(a,1));
259
link = gel(a,2);
260
v = gel(a,3);
261
w = gel(a,4);
262
}
263
else
264
{
265
e = 1;
266
v = cgetg(2*k-2 + 1, t_VEC);
267
w = cgetg(2*k-2 + 1, t_VEC);
268
link=cgetg(2*k-2 + 1, t_VECSMALL);
269
BuildTree(link, v, w, a, T? FpX_red(T,p): NULL, p);
270
if (DEBUGLEVEL > 3) timer_printf(&Ti, "building tree");
271
}
272
mask = quadratic_prec_mask(e0);
273
eold = 1;
274
penew = NULL;
275
Tnew = NULL;
276
if (DEBUGLEVEL > 3) err_printf("lifting to prec %ld\n", e0);
277
while (mask > 1)
278
{
279
long enew = eold << 1;
280
if (mask & 1) enew--;
281
mask >>= 1;
282
if (enew >= e) { /* mask == 1: last iteration */
283
GEN peold = penew? penew: powiu(p, eold);
284
GEN Td = NULL, pd;
285
long d = enew - eold; /* = eold or eold-1 */
286
/* lift from p^eold to p^enew */
287
pd = (d == eold)? peold: diviiexact(peold, p); /* p^d */
288
penew = mulii(peold,pd);
289
if (T) {
290
if (Tnew)
291
Td = (d == eold)? Tnew: FpX_red(Tnew,pd);
292
else
293
Td = FpX_red(T, peold);
294
Tnew = FpX_red(T, penew);
295
ZpXQ_RecTreeLift(link, v, w, Td, Tnew, pd, peold, penew, f, lg(v)-2,
296
(flag == 0 && mask == 1));
297
}
298
else
299
ZpX_RecTreeLift(link, v, w, pd, peold, penew, f, lg(v)-2,
300
(flag == 0 && mask == 1));
301
if (DEBUGLEVEL > 3) timer_printf(&Ti, "reaching prec %ld", enew);
302
}
303
eold = enew;
304
}
305
306
if (flag)
307
E = mkvec4(utoipos(e0), link, v, w);
308
else
309
{
310
E = cgetg(k+1, t_VEC);
311
for (i = 1; i <= 2*k-2; i++)
312
{
313
long t = link[i];
314
if (t < 0) gel(E,-t) = gel(v,i);
315
}
316
}
317
return E;
318
}
319
320
/* Q list of (coprime, monic) factors of pol mod (T,p). Lift mod p^e = pe.
321
* T may be NULL */
322
GEN
323
ZpX_liftfact(GEN pol, GEN Q, GEN pe, GEN p, long e)
324
{
325
pari_sp av = avma;
326
pol = FpX_normalize(pol, pe);
327
if (lg(Q) == 2) return mkvec(pol);
328
return gerepilecopy(av, MultiLift(pol, Q, NULL, p, e, 0));
329
}
330
331
GEN
332
ZpXQX_liftfact(GEN pol, GEN Q, GEN T, GEN pe, GEN p, long e)
333
{
334
pari_sp av = avma;
335
pol = FpXQX_normalize(pol, T, pe);
336
if (lg(Q) == 2) return mkvec(pol);
337
return gerepilecopy(av, MultiLift(pol, Q, T, p, e, 0));
338
}
339
340
GEN
341
ZqX_liftfact(GEN f, GEN a, GEN T, GEN pe, GEN p, long e)
342
{ return T ? ZpXQX_liftfact(f, a, T, pe, p, e): ZpX_liftfact(f, a, pe, p, e); }
343
GEN
344
ZqX_liftroot(GEN f, GEN a, GEN T, GEN p, long e)
345
{ return T ? ZpXQX_liftroot(f, a,T , p, e): ZpX_liftroot(f, a, p, e); }
346
347
/* U = NULL treated as 1 */
348
static void
349
BezoutPropagate(GEN link, GEN v, GEN w, GEN pe, GEN U, GEN f, long j)
350
{
351
GEN Q, R;
352
if (j < 0) return;
353
354
Q = FpX_mul(gel(v,j), gel(w,j), pe);
355
if (U)
356
{
357
Q = FpXQ_mul(Q, U, f, pe);
358
R = FpX_sub(U, Q, pe);
359
}
360
else
361
R = Fp_FpX_sub(gen_1, Q, pe);
362
gel(w,j+1) = Q; /* 0 mod U v[j], 1 mod (1-U) v[j+1] */
363
gel(w,j) = R; /* 1 mod U v[j], 0 mod (1-U) v[j+1] */
364
BezoutPropagate(link, v, w, pe, R, f, link[j ]);
365
BezoutPropagate(link, v, w, pe, Q, f, link[j+1]);
366
}
367
368
/* as above, but return the Bezout coefficients for the lifted modular factors
369
* U[i] = 1 mod Qlift[i]
370
* 0 mod Qlift[j], j != i */
371
GEN
372
bezout_lift_fact(GEN pol, GEN Q, GEN p, long e)
373
{
374
pari_sp av = avma;
375
GEN E, link, v, w, pe;
376
long i, k = lg(Q)-1;
377
if (k == 1) retmkvec(pol_1(varn(pol)));
378
pe = powiu(p, e);
379
pol = FpX_normalize(pol, pe);
380
E = MultiLift(pol, Q, NULL, p, e, 1);
381
link = gel(E,2);
382
v = gel(E,3);
383
w = gel(E,4);
384
BezoutPropagate(link, v, w, pe, NULL, pol, lg(v)-2);
385
E = cgetg(k+1, t_VEC);
386
for (i = 1; i <= 2*k-2; i++)
387
{
388
long t = link[i];
389
if (t < 0) E[-t] = w[i];
390
}
391
return gerepilecopy(av, E);
392
}
393
394
/* Front-end for ZpX_liftfact:
395
lift the factorization of pol mod p given by L to p^N (if possible) */
396
GEN
397
polhensellift(GEN pol, GEN L, GEN Tp, long N)
398
{
399
GEN T, p;
400
long i, l;
401
pari_sp av = avma;
402
void (*chk)(GEN, const char*);
403
404
if (typ(pol) != t_POL) pari_err_TYPE("polhensellift",pol);
405
RgX_check_ZXX(pol, "polhensellift");
406
if (!is_vec_t(typ(L)) || lg(L) < 3) pari_err_TYPE("polhensellift",L);
407
if (N < 1) pari_err_DOMAIN("polhensellift", "precision", "<", gen_1,stoi(N));
408
if (!ff_parse_Tp(Tp, &T, &p, 0)) pari_err_TYPE("polhensellift",Tp);
409
chk = T? RgX_check_ZXX: RgX_check_ZX;
410
l = lg(L); L = leafcopy(L);
411
for (i = 1; i < l; i++)
412
{
413
GEN q = gel(L,i);
414
if (typ(q) != t_POL) gel(L,i) = scalar_ZX_shallow(q, varn(pol));
415
else chk(q, "polhensellift");
416
}
417
return gerepilecopy(av, ZqX_liftfact(pol, L, T, powiu(p,N), p, N));
418
}
419
420
static GEN
421
FqV_roots_from_deg1(GEN x, GEN T, GEN p)
422
{
423
long i,l = lg(x);
424
GEN r = cgetg(l,t_COL);
425
for (i=1; i<l; i++) { GEN P = gel(x,i); gel(r,i) = Fq_neg(gel(P,2), T, p); }
426
return r;
427
}
428
429
static GEN
430
ZpX_liftroots_full(GEN f, GEN S, GEN q, GEN p, long e)
431
{
432
pari_sp av = avma;
433
GEN y = ZpX_liftfact(f, deg1_from_roots(S, varn(f)), q, p, e);
434
return gerepileupto(av, FqV_roots_from_deg1(y, NULL, q));
435
}
436
437
GEN
438
ZpX_roots(GEN F, GEN p, long e)
439
{
440
pari_sp av = avma;
441
GEN q = powiu(p, e);
442
GEN f = FpX_normalize(F, p);
443
GEN g = FpX_normalize(FpX_split_part(f, p), p);
444
GEN S;
445
if (degpol(g) < degpol(f))
446
{
447
GEN h = FpX_div(f, g, p);
448
F = gel(ZpX_liftfact(F, mkvec2(g, h), q, p, e), 1);
449
}
450
S = FpX_roots(g, p);
451
return gerepileupto(av, ZpX_liftroots_full(F, S, q, p, e));
452
}
453
454
static GEN
455
ZpXQX_liftroots_full(GEN f, GEN S, GEN T, GEN q, GEN p, long e)
456
{
457
pari_sp av = avma;
458
GEN y = ZpXQX_liftfact(f, deg1_from_roots(S, varn(f)), T, q, p, e);
459
return gerepileupto(av, FqV_roots_from_deg1(y, T, q));
460
}
461
462
GEN
463
ZpXQX_roots(GEN F, GEN T, GEN p, long e)
464
{
465
pari_sp av = avma;
466
GEN q = powiu(p, e);
467
GEN f = FpXQX_normalize(F, T, p);
468
GEN g = FpXQX_normalize(FpXQX_split_part(f, T, p), T, p);
469
GEN S;
470
if (degpol(g) < degpol(f))
471
{
472
GEN h = FpXQX_div(f, g, T, p);
473
F = gel(ZpXQX_liftfact(F, mkvec2(g, h), T, q, p, e), 1);
474
}
475
S = FpXQX_roots(g, T, p);
476
return gerepileupto(av, ZpXQX_liftroots_full(F, S, T, q, p, e));
477
}
478
479
GEN
480
ZqX_roots(GEN F, GEN T, GEN p, long e)
481
{
482
return T ? ZpXQX_roots(F, T, p, e): ZpX_roots(F, p, e);
483
}
484
485
/* SPEC:
486
* p is a t_INT > 1, e >= 1
487
* f is a ZX with leading term prime to p.
488
* a is a simple root mod l for all l|p.
489
* Return roots of f mod p^e, as integers (implicitly mod p^e)
490
* STANDARD USE: p is a prime power */
491
GEN
492
ZpX_liftroot(GEN f, GEN a, GEN p, long e)
493
{
494
pari_sp av = avma;
495
GEN q = p, fr, W;
496
ulong mask;
497
498
a = modii(a,q);
499
if (e == 1) return a;
500
mask = quadratic_prec_mask(e);
501
fr = FpX_red(f,q);
502
W = Fp_inv(FpX_eval(ZX_deriv(fr), a, q), q); /* 1/f'(a) mod p */
503
for(;;)
504
{
505
q = sqri(q);
506
if (mask & 1) q = diviiexact(q, p);
507
mask >>= 1;
508
fr = FpX_red(f,q);
509
a = Fp_sub(a, Fp_mul(W, FpX_eval(fr, a,q), q), q);
510
if (mask == 1) return gerepileuptoint(av, a);
511
W = Fp_sub(shifti(W,1), Fp_mul(Fp_sqr(W,q), FpX_eval(ZX_deriv(fr),a,q), q), q);
512
}
513
}
514
515
GEN
516
ZpX_liftroots(GEN f, GEN S, GEN p, long e)
517
{
518
long i, n = lg(S)-1, d = degpol(f);
519
GEN r;
520
if (n == d) return ZpX_liftroots_full(f, S, powiu(p, e), p, e);
521
r = cgetg(n+1, typ(S));
522
for (i=1; i <= n; i++)
523
gel(r,i) = ZpX_liftroot(f, gel(S,i), p, e);
524
return r;
525
}
526
527
GEN
528
ZpXQX_liftroot_vald(GEN f, GEN a, long v, GEN T, GEN p, long e)
529
{
530
pari_sp av = avma, av2;
531
GEN pv = p, q, W, df, Tq, fr, dfr;
532
ulong mask;
533
a = Fq_red(a, T, p);
534
if (e <= v+1) return a;
535
df = RgX_deriv(f);
536
if (v) { pv = powiu(p,v); df = ZXX_Z_divexact(df, pv); }
537
mask = quadratic_prec_mask(e-v);
538
Tq = FpXT_red(T, p); dfr = FpXQX_red(df, Tq, p);
539
W = Fq_inv(FqX_eval(dfr, a, Tq, p), Tq, p); /* 1/f'(a) mod (T,p) */
540
q = p;
541
av2 = avma;
542
for (;;)
543
{
544
GEN u, fa, qv, q2v, q2, Tq2;
545
q2 = q; q = sqri(q);
546
if (mask & 1) q = diviiexact(q,p);
547
mask >>= 1;
548
if (v) { qv = mulii(q, pv); q2v = mulii(q2, pv); }
549
else { qv = q; q2v = q2; }
550
Tq2 = FpXT_red(T, q2v); Tq = FpXT_red(T, qv);
551
fr = FpXQX_red(f, Tq, qv);
552
fa = FqX_eval(fr, a, Tq, qv);
553
fa = typ(fa)==t_INT? diviiexact(fa,q2v): ZX_Z_divexact(fa, q2v);
554
a = Fq_sub(a, Fq_mul(Fq_mul(W,fa,Tq2,q2v),q2, Tq,qv), Tq, qv);
555
if (mask == 1) return gerepileupto(av, a);
556
dfr = FpXQX_red(df, Tq, q);
557
u = Fq_sub(Fq_mul(W,FqX_eval(dfr,a,Tq,q),Tq,q),gen_1,Tq,q);
558
u = typ(u)==t_INT? diviiexact(u,q2): ZX_Z_divexact(u,q2);
559
W = Fq_sub(W, Fq_mul(Fq_mul(u,W,Tq2,q2),q2, Tq,q),Tq,q);
560
if (gc_needed(av2,2))
561
{
562
if(DEBUGMEM>1) pari_warn(warnmem,"ZpXQX_liftroot, e = %ld", e);
563
gerepileall(av2, 3, &a, &W, &q);
564
}
565
}
566
}
567
568
GEN
569
ZpXQX_liftroot(GEN f, GEN a, GEN T, GEN p, long e) { return ZpXQX_liftroot_vald(f,a,0,T,p,e); }
570
571
GEN
572
ZpXQX_liftroots(GEN f, GEN S, GEN T, GEN p, long e)
573
{
574
long i, n = lg(S)-1, d = degpol(f);
575
GEN r;
576
if (n == d) return ZpXQX_liftroots_full(f, S, T, powiu(p, e), p, e);
577
r = cgetg(n+1, typ(S));
578
for (i=1; i <= n; i++)
579
gel(r,i) = ZpXQX_liftroot(f, gel(S,i), T, p, e);
580
return r;
581
}
582
583
/* Same as ZpX_liftroot for the polynomial X^n-b*/
584
GEN
585
Zp_sqrtnlift(GEN b, GEN n, GEN a, GEN p, long e)
586
{
587
pari_sp ltop=avma;
588
GEN q, w, n_1;
589
ulong mask;
590
long pis2 = equalii(n, gen_2)? 1: 0;
591
if (e == 1) return icopy(a);
592
n_1 = subiu(n,1);
593
mask = quadratic_prec_mask(e);
594
w = Fp_inv(pis2 ? shifti(a,1): Fp_mul(n,Fp_pow(a,n_1,p), p), p);
595
q = p;
596
for(;;)
597
{
598
q = sqri(q);
599
if (mask & 1) q = diviiexact(q, p);
600
mask >>= 1;
601
if (lgefint(q) == 3 && lgefint(n) == 3)
602
{
603
ulong Q = uel(q,2), N = uel(n,2);
604
ulong A = umodiu(a, Q);
605
ulong B = umodiu(b, Q);
606
ulong W = umodiu(w, Q);
607
608
A = Fl_sub(A, Fl_mul(W, Fl_sub(Fl_powu(A,N,Q), B, Q), Q), Q);
609
a = utoi(A);
610
if (mask == 1) break;
611
W = Fl_sub(Fl_add(W,W,Q),
612
Fl_mul(Fl_sqr(W,Q), Fl_mul(N,Fl_powu(A, N-1, Q), Q), Q), Q);
613
w = utoi(W);
614
}
615
else
616
{
617
/* a -= w (a^n - b) */
618
a = modii(subii(a, mulii(w, subii(Fp_pow(a,n,q),b))), q);
619
if (mask == 1) break;
620
/* w += w - w^2 n a^(n-1)*/
621
w = subii(shifti(w,1), Fp_mul(Fp_sqr(w,q),
622
pis2? shifti(a,1): mulii(n,Fp_pow(a,n_1,q)), q));
623
}
624
}
625
return gerepileuptoint(ltop,a);
626
}
627
628
/* Same as ZpX_liftroot for the polynomial X^2-b */
629
GEN
630
Zp_sqrtlift(GEN b, GEN a, GEN p, long e)
631
{
632
return Zp_sqrtnlift(b, gen_2, a, p, e);
633
}
634
635
GEN
636
Zp_sqrt(GEN x, GEN p, long e)
637
{
638
pari_sp av;
639
GEN z;
640
if (absequaliu(p,2)) return Z2_sqrt(x,e);
641
av = avma;
642
z = Fp_sqrt(Fp_red(x, p), p);
643
if (!z) return NULL;
644
if (e > 1) z = Zp_sqrtlift(x, z, p, e);
645
return gerepileuptoint(av, z);
646
}
647
648
/* Compute (x-1)/(x+1)/p^k */
649
static GEN
650
ZpXQ_log_to_ath(GEN x, long k, GEN T, GEN p, long e, GEN pe)
651
{
652
pari_sp av = avma;
653
long vT = get_FpX_var(T);
654
GEN bn, bdi;
655
GEN bd = ZX_Z_add(x, gen_1);
656
if (absequaliu(p,2)) /*For p=2, we need to simplify by 2*/
657
{
658
bn = ZX_shifti(x,-(k+1));
659
bdi= ZpXQ_invlift(ZX_shifti(bd ,-1), pol_1(vT), T, p, e);
660
}
661
else
662
{
663
bn = ZX_Z_divexact(ZX_Z_sub(x, gen_1),powiu(p,k));
664
bdi= ZpXQ_invlift(bd, scalarpol(Fp_inv(gen_2,p),vT), T, p, e);
665
}
666
return gerepileupto(av, FpXQ_mul(bn, bdi, T, pe));
667
}
668
669
/* Assume p odd, a = 1 [p], return log(a) */
670
GEN
671
ZpXQ_log(GEN a, GEN T, GEN p, long N)
672
{
673
pari_sp av = avma;
674
pari_timer ti;
675
long is2 = absequaliu(p,2);
676
ulong pp = is2 ? 0: itou_or_0(p);
677
double lp = is2 ? 1: pp ? log2(pp): expi(p);
678
long k = maxss(1 , (long) .5+pow((double)(N>>1)/(lp*lp), 1./3));
679
GEN ak, s, b, pol;
680
long e = is2 ? N-1: N;
681
long i, l = (e-2)/(2*(k+is2));
682
GEN pe = powiu(p,e);
683
GEN TNk, pNk = powiu(p,N+k);
684
if( DEBUGLEVEL>=3) timer_start(&ti);
685
TNk = FpX_get_red(get_FpX_mod(T), pNk);
686
ak = FpXQ_pow(a, powiu(p,k), TNk, pNk);
687
if( DEBUGLEVEL>=3) timer_printf(&ti,"FpXQ_pow(%ld)",k);
688
b = ZpXQ_log_to_ath(ak, k, T, p, e, pe);
689
if( DEBUGLEVEL>=3) timer_printf(&ti,"ZpXQ_log_to_ath");
690
pol= cgetg(l+3,t_POL);
691
pol[1] = evalsigne(1)|evalvarn(0);
692
for(i=0; i<=l; i++)
693
{
694
GEN g;
695
ulong z = 2*i+1;
696
if (pp)
697
{
698
long w = u_lvalrem(z, pp, &z);
699
g = powuu(pp,2*i*k-w);
700
}
701
else g = powiu(p,2*i*k);
702
gel(pol,i+2) = Fp_divu(g, z,pe);
703
}
704
if( DEBUGLEVEL>=3) timer_printf(&ti,"pol(%ld)",l);
705
s = FpX_FpXQ_eval(pol, FpXQ_sqr(b, T, pe), T, pe);
706
if( DEBUGLEVEL>=3) timer_printf(&ti,"FpX_FpXQ_eval");
707
s = ZX_shifti(FpXQ_mul(b, s, T, pe), 1);
708
return gerepileupto(av, is2? s: FpX_red(s, pe));
709
}
710
711
/***********************************************************************/
712
/** **/
713
/** Generic quadratic hensel lift over Zp[X] **/
714
/** **/
715
/***********************************************************************/
716
/* q = p^N */
717
718
GEN
719
gen_ZpM_Dixon(GEN F, GEN V, GEN q, GEN p, long N, void *E,
720
GEN lin(void *E, GEN F, GEN d, GEN q),
721
GEN invl(void *E, GEN d))
722
{
723
pari_sp av = avma;
724
long N2, M;
725
GEN VN2, V2, VM, bil;
726
GEN q2, qM;
727
V = FpM_red(V, q);
728
if (N == 1) return invl(E, V);
729
N2 = (N + 1)>>1; M = N - N2;
730
F = FpM_red(F, q);
731
qM = powiu(p, M);
732
q2 = M == N2? qM: mulii(qM, p);
733
/* q2 = p^N2, qM = p^M, q = q2 * qM */
734
VN2 = gen_ZpM_Dixon(F, V, q2, p, N2, E, lin, invl);
735
bil = lin(E, F, VN2, q);
736
V2 = ZM_Z_divexact(ZM_sub(V, bil), q2);
737
VM = gen_ZpM_Dixon(F, V2, qM, p, M, E, lin, invl);
738
return gerepileupto(av, FpM_red(ZM_add(VN2, ZM_Z_mul(VM, q2)), q));
739
}
740
741
GEN
742
gen_ZpX_Dixon(GEN F, GEN V, GEN q, GEN p, long N, void *E,
743
GEN lin(void *E, GEN F, GEN d, GEN q),
744
GEN invl(void *E, GEN d))
745
{
746
pari_sp av = avma;
747
long N2, M;
748
GEN VN2, V2, VM, bil;
749
GEN q2, qM;
750
V = FpX_red(V, q);
751
if (N == 1) return invl(E, V);
752
N2 = (N + 1)>>1; M = N - N2;
753
F = FpXT_red(F, q);
754
qM = powiu(p, M);
755
q2 = M == N2? qM: mulii(qM, p);
756
/* q2 = p^N2, qM = p^M, q = q2 * qM */
757
VN2 = gen_ZpX_Dixon(F, V, q2, p, N2, E, lin, invl);
758
bil = lin(E, F, VN2, q);
759
V2 = ZX_Z_divexact(ZX_sub(V, bil), q2);
760
VM = gen_ZpX_Dixon(F, V2, qM, p, M, E, lin, invl);
761
return gerepileupto(av, FpX_red(ZX_add(VN2, ZX_Z_mul(VM, q2)), q));
762
}
763
764
GEN
765
gen_ZpM_Newton(GEN x, GEN p, long n, void *E,
766
GEN eval(void *E, GEN f, GEN q),
767
GEN invd(void *E, GEN V, GEN v, GEN q, long M))
768
{
769
pari_sp ltop = avma, av;
770
long N = 1, N2, M;
771
long mask;
772
GEN q = p;
773
if (n == 1) return gcopy(x);
774
mask = quadratic_prec_mask(n);
775
av = avma;
776
while (mask > 1)
777
{
778
GEN qM, q2, v, V;
779
N2 = N; N <<= 1;
780
q2 = q;
781
if (mask&1UL) { /* can never happen when q2 = p */
782
N--; M = N2-1;
783
qM = diviiexact(q2,p); /* > 1 */
784
q = mulii(qM,q2);
785
} else {
786
M = N2;
787
qM = q2;
788
q = sqri(q2);
789
}
790
/* q2 = p^N2, qM = p^M, q = p^N = q2 * qM */
791
mask >>= 1;
792
v = eval(E, x, q);
793
V = ZM_Z_divexact(gel(v,1), q2);
794
x = FpM_sub(x, ZM_Z_mul(invd(E, V, v, qM, M), q2), q);
795
if (gc_needed(av, 1))
796
{
797
if(DEBUGMEM>1) pari_warn(warnmem,"gen_ZpX_Newton");
798
gerepileall(av, 2, &x, &q);
799
}
800
}
801
return gerepileupto(ltop, x);
802
}
803
804
static GEN
805
_ZpM_invd(void *E, GEN V, GEN v, GEN q, long M/*unused*/)
806
{
807
(void)E; (void)M;
808
return FpM_mul(V, gel(v,2), q);
809
}
810
811
static GEN
812
_ZpM_eval(void *E, GEN x, GEN q)
813
{
814
GEN f = RgM_Rg_sub_shallow(FpM_mul(x, FpM_red((GEN) E, q), q), gen_1);
815
return mkvec2(f, x);
816
}
817
818
GEN
819
ZpM_invlift(GEN M, GEN C, GEN p, long n)
820
{
821
return gen_ZpM_Newton(C, p, n, (void *)M, _ZpM_eval, _ZpM_invd);
822
}
823
824
GEN
825
gen_ZpX_Newton(GEN x, GEN p, long n, void *E,
826
GEN eval(void *E, GEN f, GEN q),
827
GEN invd(void *E, GEN V, GEN v, GEN q, long M))
828
{
829
pari_sp ltop = avma, av;
830
long N = 1, N2, M;
831
long mask;
832
GEN q = p;
833
if (n == 1) return gcopy(x);
834
mask = quadratic_prec_mask(n);
835
av = avma;
836
while (mask > 1)
837
{
838
GEN qM, q2, v, V;
839
N2 = N; N <<= 1;
840
q2 = q;
841
if (mask&1UL) { /* can never happen when q2 = p */
842
N--; M = N2-1;
843
qM = diviiexact(q2,p); /* > 1 */
844
q = mulii(qM,q2);
845
} else {
846
M = N2;
847
qM = q2;
848
q = sqri(q2);
849
}
850
/* q2 = p^N2, qM = p^M, q = p^N = q2 * qM */
851
mask >>= 1;
852
v = eval(E, x, q);
853
V = ZX_Z_divexact(gel(v,1), q2);
854
x = FpX_sub(x, ZX_Z_mul(invd(E, V, v, qM, M), q2), q);
855
if (gc_needed(av, 1))
856
{
857
if(DEBUGMEM>1) pari_warn(warnmem,"gen_ZpX_Newton");
858
gerepileall(av, 2, &x, &q);
859
}
860
}
861
return gerepileupto(ltop, x);
862
}
863
864
struct _ZpXQ_inv
865
{
866
GEN T, a, p ,n;
867
};
868
869
static GEN
870
_inv_invd(void *E, GEN V, GEN v, GEN q, long M/*unused*/)
871
{
872
struct _ZpXQ_inv *d = (struct _ZpXQ_inv *) E;
873
GEN Tq = FpXT_red(d->T, q);
874
(void)M;
875
return FpXQ_mul(V, gel(v,2), Tq, q);
876
}
877
878
static GEN
879
_inv_eval(void *E, GEN x, GEN q)
880
{
881
struct _ZpXQ_inv *d = (struct _ZpXQ_inv *) E;
882
GEN Tq = FpXT_red(d->T, q);
883
GEN f = FpX_Fp_sub(FpXQ_mul(x, FpX_red(d->a, q), Tq, q), gen_1, q);
884
return mkvec2(f, x);
885
}
886
887
GEN
888
ZpXQ_invlift(GEN a, GEN x, GEN T, GEN p, long e)
889
{
890
struct _ZpXQ_inv d;
891
d.a = a; d.T = T; d.p = p;
892
return gen_ZpX_Newton(x, p, e, &d, _inv_eval, _inv_invd);
893
}
894
895
GEN
896
ZpXQ_inv(GEN a, GEN T, GEN p, long e)
897
{
898
pari_sp av=avma;
899
GEN ai;
900
if (lgefint(p)==3)
901
{
902
ulong pp = p[2];
903
ai = Flx_to_ZX(Flxq_inv(ZX_to_Flx(a,pp), ZXT_to_FlxT(T, pp), pp));
904
} else
905
ai = FpXQ_inv(FpX_red(a,p), FpXT_red(T,p),p);
906
return gerepileupto(av, ZpXQ_invlift(a, ai, T, p, e));
907
}
908
909
GEN
910
ZpXQ_div(GEN a, GEN b, GEN T, GEN q, GEN p, long e)
911
{
912
return FpXQ_mul(a, ZpXQ_inv(b, T, p, e), T, q);
913
}
914
915
GEN
916
ZpXQX_divrem(GEN x, GEN Sp, GEN T, GEN q, GEN p, long e, GEN *pr)
917
{
918
pari_sp av = avma;
919
GEN S = get_FpXQX_mod(Sp);
920
GEN b = leading_coeff(S), bi;
921
GEN S2, Q;
922
if (typ(b)==t_INT) return FpXQX_divrem(x, Sp, T, q, pr);
923
bi = ZpXQ_inv(b, T, p, e);
924
S2 = FqX_Fq_mul_to_monic(S, bi, T, q);
925
Q = FpXQX_divrem(x, S2, T, q, pr);
926
if (pr==ONLY_DIVIDES && !Q) { set_avma(av); return NULL; }
927
if (pr==ONLY_REM || pr==ONLY_DIVIDES) return gerepileupto(av, Q);
928
Q = FpXQX_FpXQ_mul(Q, bi, T, q);
929
gerepileall(av, 2, &Q, pr);
930
return Q;
931
}
932
933
GEN
934
ZpXQX_digits(GEN x, GEN B, GEN T, GEN q, GEN p, long e)
935
{
936
pari_sp av = avma;
937
GEN b = leading_coeff(B), bi;
938
GEN B2, P, V, W;
939
long i, lV;
940
if (typ(b)==t_INT) return FpXQX_digits(x, B, T, q);
941
bi = ZpXQ_inv(b, T, p, e);
942
B2 = FqX_Fq_mul_to_monic(B, bi, T, q);
943
V = FpXQX_digits(x, B2, T, q);
944
lV = lg(V)-1;
945
P = FpXQ_powers(bi, lV-1, T, q);
946
W = cgetg(lV+1, t_VEC);
947
for(i=1; i<=lV; i++)
948
gel(W, i) = FpXQX_FpXQ_mul(gel(V,i), gel(P, i), T, q);
949
return gerepileupto(av, W);
950
}
951
952
struct _ZpXQ_sqrtn
953
{
954
GEN T, a, n, ai;
955
};
956
957
static GEN
958
_sqrtn_invd(void *E, GEN V, GEN v, GEN q, long M)
959
{
960
struct _ZpXQ_sqrtn *d = (struct _ZpXQ_sqrtn *) E;
961
GEN Tq = FpX_red(d->T, q), aiq = FpX_red(d->ai, q);
962
(void)M;
963
return FpXQ_mul(FpXQ_mul(V, gel(v,2), Tq, q), aiq, Tq, q);
964
}
965
966
static GEN
967
_sqrtn_eval(void *E, GEN x, GEN q)
968
{
969
struct _ZpXQ_sqrtn *d = (struct _ZpXQ_sqrtn *) E;
970
GEN Tq = FpX_red(d->T, q);
971
GEN f = FpX_sub(FpXQ_pow(x, d->n, Tq, q), d->a, q);
972
return mkvec2(f, x);
973
}
974
975
GEN
976
ZpXQ_sqrtnlift(GEN a, GEN n, GEN x, GEN T, GEN p, long e)
977
{
978
struct _ZpXQ_sqrtn d;
979
d.a = a; d.T = T; d.n = n;
980
d.ai = ZpXQ_inv(ZX_Z_mul(a, n),T,p,(e+1)>>1);
981
return gen_ZpX_Newton(x, p, e, &d, _sqrtn_eval, _sqrtn_invd);
982
}
983
984
static GEN
985
to_ZX(GEN a, long v) { return typ(a)==t_INT? scalarpol_shallow(a,v): a; }
986
987
GEN
988
Zq_sqrtnlift(GEN a, GEN n, GEN x, GEN T, GEN p, long e)
989
{
990
return T? ZpXQ_sqrtnlift(to_ZX(a,varn(T)), n, to_ZX(x,varn(T)), T, p, e)
991
: Zp_sqrtnlift(a, n, x, p, e);
992
}
993
994
GEN
995
ZpXQ_sqrt(GEN a, GEN T, GEN p, long e)
996
{
997
pari_sp av = avma;
998
GEN z = FpXQ_sqrt(FpX_red(a, p), T, p);
999
if (!z) return NULL;
1000
if (e <= 1) return gerepileupto(av, z);
1001
return gerepileupto(av, ZpXQ_sqrtnlift(a, gen_2, z, T, p, e));
1002
}
1003
1004
GEN
1005
ZpX_ZpXQ_liftroot_ea(GEN P, GEN S, GEN T, GEN p, long n, void *E,
1006
GEN early(void *E, GEN x, GEN q))
1007
{
1008
pari_sp ltop = avma, av;
1009
long N, r;
1010
long mask;
1011
GEN q2, q, W, Q, Tq2, Tq, Pq;
1012
pari_timer ti;
1013
T = FpX_get_red(T, powiu(p, n));
1014
if (n == 1) return gcopy(S);
1015
mask = quadratic_prec_mask(n);
1016
av = avma;
1017
q2 = p; q = sqri(p); mask >>= 1; N = 2;
1018
if (DEBUGLEVEL > 3) timer_start(&ti);
1019
Tq = FpXT_red(T,q);
1020
Tq2 = FpXT_red(Tq,q2);
1021
Pq = FpX_red(P,q);
1022
W = FpXQ_inv(FpX_FpXQ_eval(FpX_deriv(P,q2), S, Tq2, q2), Tq2, q2);
1023
Q = ZX_Z_divexact(FpX_FpXQ_eval(Pq, S, Tq, q), q2);
1024
r = brent_kung_optpow(degpol(P), 4, 3);
1025
if (DEBUGLEVEL > 3)
1026
err_printf("ZpX_ZpXQ_liftroot: lifting to prec %ld\n",n);
1027
for (;;)
1028
{
1029
GEN H, Sq, Wq, Spow, dP, qq, Pqq, Tqq;
1030
H = FpXQ_mul(W, Q, Tq2, q2);
1031
Sq = FpX_sub(S, ZX_Z_mul(H, q2), q);
1032
if (DEBUGLEVEL > 3)
1033
timer_printf(&ti,"ZpX_ZpXQ_liftroot: reaching prec %ld",N);
1034
if (mask==1)
1035
return gerepileupto(ltop, Sq);
1036
if (early)
1037
{
1038
GEN Se = early(E, Sq, q);
1039
if (Se) return gerepileupto(ltop, Se);
1040
}
1041
qq = sqri(q); N <<= 1;
1042
if (mask&1UL) { qq = diviiexact(qq, p); N--; }
1043
mask >>= 1;
1044
Pqq = FpX_red(P, qq);
1045
Tqq = FpXT_red(T, qq);
1046
Spow = FpXQ_powers(Sq, r, Tqq, qq);
1047
Q = ZX_Z_divexact(FpX_FpXQV_eval(Pqq, Spow, Tqq, qq), q);
1048
dP = FpX_FpXQV_eval(FpX_deriv(Pq, q), FpXV_red(Spow, q), Tq, q);
1049
Wq = ZX_Z_divexact(FpX_Fp_sub(FpXQ_mul(W, dP, Tq, q), gen_1, q), q2);
1050
Wq = ZX_Z_mul(FpXQ_mul(W, Wq, Tq2, q2), q2);
1051
Wq = FpX_sub(W, Wq, q);
1052
S = Sq; W = Wq; q2 = q; q = qq; Tq2 = Tq; Tq = Tqq; Pq = Pqq;
1053
if (gc_needed(av, 1))
1054
{
1055
if(DEBUGMEM>1) pari_warn(warnmem,"ZpX_ZpXQ_Newton");
1056
gerepileall(av, 8, &S, &W, &Q, &Tq2, &Tq, &Pq, &q, &q2);
1057
}
1058
}
1059
}
1060
1061
GEN
1062
ZpX_ZpXQ_liftroot(GEN P, GEN S, GEN T, GEN p, long n)
1063
{
1064
return ZpX_ZpXQ_liftroot_ea(P, S, T, p, n, NULL, NULL);
1065
}
1066
1067
GEN
1068
ZpX_Frobenius(GEN T, GEN p, long e)
1069
{
1070
return ZpX_ZpXQ_liftroot(get_FpX_mod(T), FpX_Frobenius(T, p), T, p, e);
1071
}
1072
1073
GEN
1074
ZpXQM_prodFrobenius(GEN M, GEN T, GEN p, long e)
1075
{
1076
pari_sp av = avma;
1077
GEN xp = ZpX_Frobenius(T, p, e);
1078
GEN z = FpXQM_autsum(mkvec2(xp, M), get_FpX_degree(T), T, powiu(p,e));
1079
return gerepilecopy(av, gel(z,2));
1080
}
1081
1082
/* Canonical lift of polynomial */
1083
1084
static GEN _can_invl(void *E, GEN V) {(void) E; return V; }
1085
1086
static GEN _can_lin(void *E, GEN F, GEN V, GEN q)
1087
{
1088
GEN v = RgX_splitting(V, 3);
1089
(void) E;
1090
return FpX_sub(V,ZXV_dotproduct(v, F), q);
1091
}
1092
1093
static GEN
1094
_can_iter(void *E, GEN f, GEN q)
1095
{
1096
GEN h = RgX_splitting(f,3);
1097
GEN h1s = ZX_sqr(gel(h,1)), h2s = ZX_sqr(gel(h,2)), h3s = ZX_sqr(gel(h,3));
1098
GEN h12 = ZX_mul(gel(h,1), gel(h,2));
1099
GEN h13 = ZX_mul(gel(h,1), gel(h,3));
1100
GEN h23 = ZX_mul(gel(h,2), gel(h,3));
1101
GEN h1c = ZX_mul(gel(h,1), h1s);
1102
GEN h3c = ZX_mul(gel(h,3), h3s);
1103
GEN th = ZX_mul(ZX_sub(h2s,ZX_mulu(h13,3)),gel(h,2));
1104
GEN y = FpX_sub(f,ZX_add(RgX_shift_shallow(h3c,2),ZX_add(RgX_shift_shallow(th,1),h1c)),q);
1105
(void) E;
1106
return mkvecn(7,y,h1s,h2s,h3s,h12,h13,h23);
1107
}
1108
1109
static GEN
1110
_can_invd(void *E, GEN V, GEN v, GEN qM, long M)
1111
{
1112
GEN h1s=gel(v,2), h2s=gel(v,3), h3s=gel(v,4);
1113
GEN h12=gel(v,5), h13=gel(v,6), h23=gel(v,7);
1114
GEN F = mkvec3(ZX_sub(h1s,RgX_shift_shallow(h23,1)),RgX_shift_shallow(ZX_sub(h2s,h13),1),
1115
ZX_sub(RgX_shift_shallow(h3s,2),RgX_shift_shallow(h12,1)));
1116
(void)E;
1117
return gen_ZpX_Dixon(ZXV_Z_mul(F, utoi(3)), V, qM, utoi(3), M, NULL,
1118
_can_lin, _can_invl);
1119
}
1120
1121
static GEN
1122
F3x_frobeniuslift(GEN P, long n)
1123
{ return gen_ZpX_Newton(Flx_to_ZX(P),utoi(3), n, NULL, _can_iter, _can_invd); }
1124
1125
static GEN _can5_invl(void *E, GEN V) {(void) E; return V; }
1126
1127
static GEN _can5_lin(void *E, GEN F, GEN V, GEN q)
1128
{
1129
ulong p = *(ulong*)E;
1130
GEN v = RgX_splitting(V, p);
1131
return FpX_sub(V,ZXV_dotproduct(v, F), q);
1132
}
1133
1134
/* P(X,t) -> P(X*t^n,t) mod (t^p-1) */
1135
static GEN
1136
_shift(GEN P, long n, ulong p, long v)
1137
{
1138
long i, l=lg(P);
1139
GEN r = cgetg(l,t_POL); r[1] = P[1];
1140
for(i=2;i<l;i++)
1141
{
1142
long s = n*(i-2)%p;
1143
GEN ci = gel(P,i);
1144
if (typ(ci)==t_INT)
1145
gel(r,i) = monomial(ci, s, v);
1146
else
1147
gel(r,i) = RgX_rotate_shallow(ci, s, p);
1148
}
1149
return FpXX_renormalize(r, l);
1150
}
1151
1152
struct _can_mul
1153
{
1154
GEN T, q;
1155
ulong p;
1156
};
1157
1158
static GEN
1159
_can5_mul(void *E, GEN A, GEN B)
1160
{
1161
struct _can_mul *d = (struct _can_mul *)E;
1162
GEN a = gel(A,1), b = gel(B,1);
1163
long n = itos(gel(A,2));
1164
GEN bn = _shift(b, n, d->p, get_FpX_var(d->T));
1165
GEN c = FpXQX_mul(a, bn, d->T, d->q);
1166
return mkvec2(c, addii(gel(A,2), gel(B,2)));
1167
}
1168
1169
static GEN
1170
_can5_sqr(void *E, GEN A)
1171
{
1172
return _can5_mul(E,A,A);
1173
}
1174
1175
static GEN
1176
_can5_iter(void *E, GEN f, GEN q)
1177
{
1178
pari_sp av = avma;
1179
struct _can_mul D;
1180
ulong p = *(ulong*)E;
1181
long i, vT = fetch_var();
1182
GEN N, P, d, V, fs;
1183
D.q = q; D.T = ZX_Z_sub(pol_xn(p,vT),gen_1);
1184
D.p = p;
1185
fs = mkvec2(_shift(f, 1, p, vT), gen_1);
1186
N = gel(gen_powu_i(fs,p-1,(void*)&D,_can5_sqr,_can5_mul),1);
1187
N = ZXX_evalx0(FpXQX_red(N,polcyclo(p,vT),q));
1188
P = FpX_mul(N,f,q);
1189
P = RgX_deflate(P, p);
1190
d = RgX_splitting(N, p);
1191
V = cgetg(p+1,t_VEC);
1192
gel(V,1) = ZX_mulu(gel(d,1), p);
1193
for(i=2; i<= (long)p; i++)
1194
gel(V,i) = ZX_mulu(RgX_shift_shallow(gel(d,p+2-i), 1), p);
1195
(void)delete_var(); return gerepilecopy(av, mkvec2(ZX_sub(f,P),V));
1196
}
1197
1198
static GEN
1199
_can5_invd(void *E, GEN H, GEN v, GEN qM, long M)
1200
{
1201
ulong p = *(long*)E;
1202
return gen_ZpX_Dixon(gel(v,2), H, qM, utoipos(p), M, E, _can5_lin, _can5_invl);
1203
}
1204
1205
GEN
1206
Flx_Teichmuller(GEN P, ulong p, long n)
1207
{
1208
return p==3 ? F3x_frobeniuslift(P,n):
1209
gen_ZpX_Newton(Flx_to_ZX(P),utoipos(p), n, &p, _can5_iter, _can5_invd);
1210
}
1211
1212
GEN
1213
polteichmuller(GEN P, ulong p, long n)
1214
{
1215
pari_sp av = avma;
1216
GEN q = NULL;
1217
if (typ(P)!=t_POL || !RgX_is_FpX(P,&q)) pari_err_TYPE("polteichmuller",P);
1218
if (q && !equaliu(q,p)) pari_err_MODULUS("polteichmuller",q,utoi(p));
1219
if (n <= 0)
1220
pari_err_DOMAIN("polteichmuller", "precision", "<=",gen_0,stoi(n));
1221
return gerepileupto(av, p==2 ? F2x_Teichmuller(RgX_to_F2x(P), n)
1222
: Flx_Teichmuller(RgX_to_Flx(P, p), p, n));
1223
}
1224
1225