Line data Source code
1 : /* Copyright (C) 2000-2005 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 :
15 : /***********************************************************************/
16 : /** **/
17 : /** ARITHMETIC OPERATIONS ON POLYNOMIALS **/
18 : /** (third part) **/
19 : /** **/
20 : /***********************************************************************/
21 : #include "pari.h"
22 : #include "paripriv.h"
23 :
24 : #define DEBUGLEVEL DEBUGLEVEL_pol
25 :
26 : /************************************************************************
27 : ** **
28 : ** Ring membership **
29 : ** **
30 : ************************************************************************/
31 : struct charact {
32 : GEN q;
33 : int isprime;
34 : };
35 : static void
36 1225 : char_update_prime(struct charact *S, GEN p)
37 : {
38 1225 : if (!S->isprime) { S->isprime = 1; S->q = p; }
39 1225 : if (!equalii(p, S->q)) pari_err_MODULUS("characteristic", S->q, p);
40 1218 : }
41 : static void
42 6594 : char_update_int(struct charact *S, GEN n)
43 : {
44 6594 : if (S->isprime)
45 : {
46 7 : if (dvdii(n, S->q)) return;
47 7 : pari_err_MODULUS("characteristic", S->q, n);
48 : }
49 6587 : S->q = gcdii(S->q, n);
50 : }
51 : static void
52 179046 : charact(struct charact *S, GEN x)
53 : {
54 179046 : const long tx = typ(x);
55 : long i, l;
56 179046 : switch(tx)
57 : {
58 5145 : case t_INTMOD:char_update_int(S, gel(x,1)); break;
59 1134 : case t_FFELT: char_update_prime(S, gel(x,4)); break;
60 26761 : case t_COMPLEX: case t_QUAD:
61 : case t_POLMOD: case t_POL: case t_SER: case t_RFRAC:
62 : case t_VEC: case t_COL: case t_MAT:
63 26761 : l = lg(x);
64 178171 : for (i=lontyp[tx]; i < l; i++) charact(S,gel(x,i));
65 26747 : break;
66 7 : case t_LIST:
67 7 : x = list_data(x);
68 7 : if (x) charact(S, x);
69 7 : break;
70 : }
71 179018 : }
72 : static void
73 4739 : charact_res(struct charact *S, GEN x)
74 : {
75 4739 : const long tx = typ(x);
76 : long i, l;
77 4739 : switch(tx)
78 : {
79 1449 : case t_INTMOD:char_update_int(S, gel(x,1)); break;
80 0 : case t_FFELT: char_update_prime(S, gel(x,4)); break;
81 91 : case t_PADIC: char_update_prime(S, padic_p(x)); break;
82 1722 : case t_COMPLEX: case t_QUAD:
83 : case t_POLMOD: case t_POL: case t_SER: case t_RFRAC:
84 : case t_VEC: case t_COL: case t_MAT:
85 1722 : l = lg(x);
86 6132 : for (i=lontyp[tx]; i < l; i++) charact_res(S,gel(x,i));
87 1722 : break;
88 0 : case t_LIST:
89 0 : x = list_data(x);
90 0 : if (x) charact_res(S, x);
91 0 : break;
92 : }
93 4739 : }
94 : GEN
95 27622 : characteristic(GEN x)
96 : {
97 : struct charact S;
98 27622 : S.q = gen_0; S.isprime = 0;
99 27622 : charact(&S, x); return S.q;
100 : }
101 : GEN
102 329 : residual_characteristic(GEN x)
103 : {
104 : struct charact S;
105 329 : S.q = gen_0; S.isprime = 0;
106 329 : charact_res(&S, x); return S.q;
107 : }
108 :
109 : int
110 71020556 : Rg_is_Fp(GEN x, GEN *pp)
111 : {
112 : GEN mod;
113 71020556 : switch(typ(x))
114 : {
115 2482676 : case t_INTMOD:
116 2482676 : mod = gel(x,1);
117 2482676 : if (!*pp) *pp = mod;
118 2341976 : else if (mod != *pp && !equalii(mod, *pp))
119 : {
120 0 : if (DEBUGLEVEL) pari_warn(warner,"different moduli in Rg_is_Fp");
121 0 : return 0;
122 : }
123 2482676 : return 1;
124 57129405 : case t_INT:
125 57129405 : return 1;
126 11408475 : default: return 0;
127 : }
128 : }
129 :
130 : int
131 28180284 : RgX_is_FpX(GEN x, GEN *pp)
132 : {
133 28180284 : long i, lx = lg(x);
134 87766228 : for (i=2; i<lx; i++)
135 70994422 : if (!Rg_is_Fp(gel(x, i), pp))
136 11408470 : return 0;
137 16771806 : return 1;
138 : }
139 :
140 : int
141 0 : RgV_is_FpV(GEN x, GEN *pp)
142 : {
143 0 : long i, lx = lg(x);
144 0 : for (i=1; i<lx; i++)
145 0 : if (!Rg_is_Fp(gel(x,i), pp)) return 0;
146 0 : return 1;
147 : }
148 :
149 : int
150 0 : RgM_is_FpM(GEN x, GEN *pp)
151 : {
152 0 : long i, lx = lg(x);
153 0 : for (i=1; i<lx; i++)
154 0 : if (!RgV_is_FpV(gel(x, i), pp)) return 0;
155 0 : return 1;
156 : }
157 :
158 : int
159 60802 : Rg_is_FpXQ(GEN x, GEN *pT, GEN *pp)
160 : {
161 : GEN pol, mod, p;
162 60802 : switch(typ(x))
163 : {
164 26131 : case t_INTMOD:
165 26131 : return Rg_is_Fp(x, pp);
166 8561 : case t_INT:
167 8561 : return 1;
168 21 : case t_POL:
169 21 : return RgX_is_FpX(x, pp);
170 21350 : case t_FFELT:
171 21350 : mod = x; p = FF_p_i(x);
172 21350 : if (!*pp) *pp = p;
173 21350 : if (!*pT) *pT = mod;
174 19824 : else if (typ(*pT)!=t_FFELT || !FF_samefield(*pT,mod))
175 : {
176 42 : if (DEBUGLEVEL) pari_warn(warner,"different moduli in Rg_is_FpXQ");
177 42 : return 0;
178 : }
179 21308 : return 1;
180 4585 : case t_POLMOD:
181 4585 : mod = gel(x,1); pol = gel(x, 2);
182 4585 : if (!RgX_is_FpX(mod, pp)) return 0;
183 4585 : if (typ(pol)==t_POL)
184 : {
185 4578 : if (!RgX_is_FpX(pol, pp)) return 0;
186 : }
187 7 : else if (!Rg_is_Fp(pol, pp)) return 0;
188 4585 : if (!*pT) *pT = mod;
189 4431 : else if (mod != *pT && !gequal(mod, *pT))
190 : {
191 0 : if (DEBUGLEVEL) pari_warn(warner,"different moduli in Rg_is_FpXQ");
192 0 : return 0;
193 : }
194 4585 : return 1;
195 :
196 154 : default: return 0;
197 : }
198 : }
199 :
200 : int
201 3381 : RgX_is_FpXQX(GEN x, GEN *pT, GEN *pp)
202 : {
203 3381 : long i, lx = lg(x);
204 63427 : for (i = 2; i < lx; i++)
205 60144 : if (!Rg_is_FpXQ(gel(x,i), pT, pp)) return 0;
206 3283 : return 1;
207 : }
208 :
209 : /************************************************************************
210 : ** **
211 : ** Ring conversion **
212 : ** **
213 : ************************************************************************/
214 :
215 : /* p > 0 a t_INT, return lift(x * Mod(1,p)).
216 : * If x is an INTMOD, assume modulus is a multiple of p. */
217 : GEN
218 52285018 : Rg_to_Fp(GEN x, GEN p)
219 : {
220 52285018 : if (lgefint(p) == 3) return utoi(Rg_to_Fl(x, uel(p,2)));
221 4555004 : switch(typ(x))
222 : {
223 288479 : case t_INT: return modii(x, p);
224 18790 : case t_FRAC: {
225 18790 : pari_sp av = avma;
226 18790 : GEN z = modii(gel(x,1), p);
227 18790 : if (z == gen_0) return gen_0;
228 18785 : return gc_INT(av, remii(mulii(z, Fp_inv(gel(x,2), p)), p));
229 : }
230 70 : case t_PADIC: return padic_to_Fp(x, p);
231 4247668 : case t_INTMOD: {
232 4247668 : GEN q = gel(x,1), a = gel(x,2);
233 4247668 : if (equalii(q, p)) return icopy(a);
234 14 : if (!dvdii(q,p)) pari_err_MODULUS("Rg_to_Fp", q, p);
235 0 : return remii(a, p);
236 : }
237 0 : default: pari_err_TYPE("Rg_to_Fp",x);
238 : return NULL; /* LCOV_EXCL_LINE */
239 : }
240 : }
241 : /* If x is a POLMOD, assume modulus is a multiple of T. */
242 : GEN
243 1291958 : Rg_to_FpXQ(GEN x, GEN T, GEN p)
244 : {
245 1291958 : long ta, tx = typ(x), v = get_FpX_var(T);
246 : GEN a, b;
247 1291958 : if (is_const_t(tx))
248 : {
249 59175 : if (tx == t_FFELT)
250 : {
251 17355 : GEN z = FF_to_FpXQ(x);
252 17355 : setvarn(z, v);
253 17355 : return z;
254 : }
255 41820 : return scalar_ZX(degpol(T)? Rg_to_Fp(x, p): gen_0, v);
256 : }
257 1232783 : switch(tx)
258 : {
259 1230676 : case t_POLMOD:
260 1230676 : b = gel(x,1);
261 1230676 : a = gel(x,2); ta = typ(a);
262 1230676 : if (is_const_t(ta))
263 3885 : return scalar_ZX(degpol(T)? Rg_to_Fp(a, p): gen_0, v);
264 1226791 : b = RgX_to_FpX(b, p); if (varn(b) != v) break;
265 1226791 : a = RgX_to_FpX(a, p);
266 1226791 : if (ZX_equal(b,get_FpX_mod(T)) || signe(FpX_rem(b,T,p))==0)
267 1226791 : return FpX_rem(a, T, p);
268 0 : break;
269 2107 : case t_POL:
270 2107 : if (varn(x) != v) break;
271 2100 : return FpX_rem(RgX_to_FpX(x,p), T, p);
272 0 : case t_RFRAC:
273 0 : a = Rg_to_FpXQ(gel(x,1), T,p);
274 0 : b = Rg_to_FpXQ(gel(x,2), T,p);
275 0 : return FpXQ_div(a,b, T,p);
276 : }
277 7 : pari_err_TYPE("Rg_to_FpXQ",x);
278 : return NULL; /* LCOV_EXCL_LINE */
279 : }
280 : GEN
281 3335006 : RgX_to_FpX(GEN x, GEN p)
282 : {
283 : long i, l;
284 3335006 : GEN z = cgetg_copy(x, &l); z[1] = x[1];
285 14762943 : for (i = 2; i < l; i++) gel(z,i) = Rg_to_Fp(gel(x,i), p);
286 3335006 : return FpX_renormalize(z, l);
287 : }
288 :
289 : GEN
290 140 : RgV_to_FpV(GEN x, GEN p)
291 483 : { pari_APPLY_type(t_VEC, Rg_to_Fp(gel(x,i), p)) }
292 :
293 : GEN
294 1751090 : RgC_to_FpC(GEN x, GEN p)
295 28499485 : { pari_APPLY_type(t_COL, Rg_to_Fp(gel(x,i), p)) }
296 :
297 : GEN
298 222349 : RgM_to_FpM(GEN x, GEN p)
299 1973397 : { pari_APPLY_same(RgC_to_FpC(gel(x,i), p)) }
300 :
301 : GEN
302 369001 : RgV_to_Flv(GEN x, ulong p)
303 1676894 : { pari_APPLY_ulong(Rg_to_Fl(gel(x,i), p)) }
304 :
305 : GEN
306 118296 : RgM_to_Flm(GEN x, ulong p)
307 422998 : { pari_APPLY_same(RgV_to_Flv(gel(x,i), p)) }
308 :
309 : GEN
310 5098 : RgX_to_FpXQX(GEN x, GEN T, GEN p)
311 : {
312 5098 : long i, l = lg(x);
313 5098 : GEN z = cgetg(l, t_POL); z[1] = x[1];
314 43366 : for (i = 2; i < l; i++) gel(z,i) = Rg_to_FpXQ(gel(x,i), T,p);
315 5098 : return FpXQX_renormalize(z, l);
316 : }
317 : GEN
318 49186 : RgX_to_FqX(GEN x, GEN T, GEN p)
319 : {
320 49186 : long i, l = lg(x);
321 49186 : GEN z = cgetg(l, t_POL); z[1] = x[1];
322 49186 : if (T)
323 10990 : for (i = 2; i < l; i++) gel(z,i) = Rg_to_FpXQ(gel(x,i), T, p);
324 : else
325 791394 : for (i = 2; i < l; i++) gel(z,i) = Rg_to_Fp(gel(x,i), p);
326 49186 : return FpXQX_renormalize(z, l);
327 : }
328 :
329 : GEN
330 219128 : RgC_to_FqC(GEN x, GEN T, GEN p)
331 : {
332 219128 : long i, l = lg(x);
333 219128 : GEN z = cgetg(l, t_COL);
334 219128 : if (T)
335 1430310 : for (i = 1; i < l; i++) gel(z,i) = Rg_to_FpXQ(gel(x,i), T, p);
336 : else
337 0 : for (i = 1; i < l; i++) gel(z,i) = Rg_to_Fp(gel(x,i), p);
338 219128 : return z;
339 : }
340 :
341 : GEN
342 52325 : RgM_to_FqM(GEN x, GEN T, GEN p)
343 271418 : { pari_APPLY_same(RgC_to_FqC(gel(x, i), T, p)) }
344 :
345 : /* lg(V) > 1 */
346 : GEN
347 851487 : FpXV_FpC_mul(GEN V, GEN W, GEN p)
348 : {
349 851487 : pari_sp av = avma;
350 851487 : long i, l = lg(V);
351 851487 : GEN z = ZX_Z_mul(gel(V,1),gel(W,1));
352 4201029 : for(i=2; i<l; i++)
353 : {
354 3349542 : z = ZX_add(z, ZX_Z_mul(gel(V,i),gel(W,i)));
355 3349542 : if ((i & 7) == 0) z = gc_upto(av, z);
356 : }
357 851487 : return gc_upto(av, FpX_red(z,p));
358 : }
359 :
360 : GEN
361 55832 : FqX_Fq_add(GEN y, GEN x, GEN T, GEN p)
362 : {
363 55832 : long i, lz = lg(y);
364 : GEN z;
365 55832 : if (!T) return FpX_Fp_add(y, x, p);
366 8666 : if (lz == 2) return scalarpol(x, varn(y));
367 7952 : z = cgetg(lz,t_POL); z[1] = y[1];
368 7952 : gel(z,2) = Fq_add(gel(y,2),x, T, p);
369 7952 : if (lz == 3) z = FpXX_renormalize(z,lz);
370 : else
371 33145 : for(i=3;i<lz;i++) gel(z,i) = gcopy(gel(y,i));
372 7952 : return z;
373 : }
374 :
375 : GEN
376 1059 : FqX_Fq_sub(GEN y, GEN x, GEN T, GEN p)
377 : {
378 1059 : long i, lz = lg(y);
379 : GEN z;
380 1059 : if (!T) return FpX_Fp_sub(y, x, p);
381 1059 : if (lz == 2) return scalarpol(x, varn(y));
382 1059 : z = cgetg(lz,t_POL); z[1] = y[1];
383 1059 : gel(z,2) = Fq_sub(gel(y,2), x, T, p);
384 1059 : if (lz == 3) z = FpXX_renormalize(z,lz);
385 : else
386 10278 : for(i=3;i<lz;i++) gel(z,i) = gcopy(gel(y,i));
387 1059 : return z;
388 : }
389 :
390 : GEN
391 149052 : FqX_Fq_mul_to_monic(GEN P, GEN U, GEN T, GEN p)
392 : {
393 : long i, lP;
394 149052 : GEN res = cgetg_copy(P, &lP); res[1] = P[1];
395 918799 : for(i=2; i<lP-1; i++) gel(res,i) = Fq_mul(U,gel(P,i), T,p);
396 149052 : gel(res,lP-1) = gen_1; return res;
397 : }
398 :
399 : GEN
400 38189 : FpXQX_normalize(GEN z, GEN T, GEN p)
401 : {
402 : GEN lc;
403 38189 : if (lg(z) == 2) return z;
404 38175 : lc = leading_coeff(z);
405 38175 : if (typ(lc) == t_POL)
406 : {
407 2167 : if (lg(lc) > 3) /* nonconstant */
408 1902 : return FqX_Fq_mul_to_monic(z, Fq_inv(lc,T,p), T,p);
409 : /* constant */
410 265 : lc = gel(lc,2);
411 265 : z = shallowcopy(z);
412 265 : gel(z, lg(z)-1) = lc;
413 : }
414 : /* lc a t_INT */
415 36273 : if (equali1(lc)) return z;
416 87 : return FqX_Fq_mul_to_monic(z, Fp_inv(lc,p), T,p);
417 : }
418 :
419 : GEN
420 390935 : FqX_eval(GEN x, GEN y, GEN T, GEN p)
421 : {
422 : pari_sp av;
423 : GEN p1, r;
424 390935 : long j, i=lg(x)-1;
425 390935 : if (i<=2)
426 45107 : return (i==2)? Fq_red(gel(x,2), T, p): gen_0;
427 345828 : av=avma; p1=gel(x,i);
428 : /* specific attention to sparse polynomials (see poleval)*/
429 : /*You've guessed it! It's a copy-paste(tm)*/
430 1150444 : for (i--; i>=2; i=j-1)
431 : {
432 805304 : for (j=i; !signe(gel(x,j)); j--)
433 686 : if (j==2)
434 : {
435 490 : if (i!=j) y = Fq_pow(y,utoipos(i-j+1), T, p);
436 490 : return gc_upto(av, Fq_mul(p1,y, T, p));
437 : }
438 804618 : r = (i==j)? y: Fq_pow(y, utoipos(i-j+1), T, p);
439 804618 : p1 = Fq_add(Fq_mul(p1,r,T,p), gel(x,j), T, p);
440 : }
441 345336 : return gc_upto(av, p1);
442 : }
443 :
444 : GEN
445 97321 : FqXY_evalx(GEN Q, GEN x, GEN T, GEN p)
446 : {
447 97321 : long i, lb = lg(Q);
448 : GEN z;
449 97321 : if (!T) return FpXY_evalx(Q, x, p);
450 86961 : z = cgetg(lb, t_POL); z[1] = Q[1];
451 462945 : for (i=2; i<lb; i++)
452 : {
453 375984 : GEN q = gel(Q,i);
454 375984 : gel(z,i) = typ(q) == t_INT? modii(q,p): FqX_eval(q, x, T, p);
455 : }
456 86961 : return FpXQX_renormalize(z, lb);
457 : }
458 :
459 : /* Q an FpXY, evaluate at (X,Y) = (x,y) */
460 : GEN
461 14623 : FqXY_eval(GEN Q, GEN y, GEN x, GEN T, GEN p)
462 : {
463 14623 : pari_sp av = avma;
464 14623 : if (!T) return FpXY_eval(Q, y, x, p);
465 966 : return gc_upto(av, FqX_eval(FqXY_evalx(Q, x, T, p), y, T, p));
466 : }
467 :
468 : /* a X^d */
469 : GEN
470 13104914 : monomial(GEN a, long d, long v)
471 : {
472 : long i, n;
473 : GEN P;
474 13104914 : if (d < 0) {
475 14 : if (isrationalzero(a)) return pol_0(v);
476 14 : retmkrfrac(a, pol_xn(-d, v));
477 : }
478 13104900 : if (gequal0(a))
479 : {
480 93989 : if (isexactzero(a)) return scalarpol_shallow(a,v);
481 0 : n = d+2; P = cgetg(n+1, t_POL);
482 0 : P[1] = evalsigne(0) | evalvarn(v);
483 : }
484 : else
485 : {
486 13010909 : n = d+2; P = cgetg(n+1, t_POL);
487 13010911 : P[1] = evalsigne(1) | evalvarn(v);
488 : }
489 32179365 : for (i = 2; i < n; i++) gel(P,i) = gen_0;
490 13010911 : gel(P,i) = a; return P;
491 : }
492 : GEN
493 157969 : monomialcopy(GEN a, long d, long v)
494 : {
495 : long i, n;
496 : GEN P;
497 157969 : if (d < 0) {
498 14 : if (isrationalzero(a)) return pol_0(v);
499 14 : retmkrfrac(gcopy(a), pol_xn(-d, v));
500 : }
501 157955 : if (gequal0(a))
502 : {
503 14 : if (isexactzero(a)) return scalarpol(a,v);
504 7 : n = d+2; P = cgetg(n+1, t_POL);
505 7 : P[1] = evalsigne(0) | evalvarn(v);
506 : }
507 : else
508 : {
509 157941 : n = d+2; P = cgetg(n+1, t_POL);
510 157941 : P[1] = evalsigne(1) | evalvarn(v);
511 : }
512 314678 : for (i = 2; i < n; i++) gel(P,i) = gen_0;
513 157948 : gel(P,i) = gcopy(a); return P;
514 : }
515 : GEN
516 325824 : pol_x_powers(long N, long v)
517 : {
518 325824 : GEN L = cgetg(N+1,t_VEC);
519 : long i;
520 788717 : for (i=1; i<=N; i++) gel(L,i) = pol_xn(i-1, v);
521 325827 : return L;
522 : }
523 :
524 : GEN
525 0 : FqXQ_powers(GEN x, long l, GEN S, GEN T, GEN p)
526 : {
527 0 : return T ? FpXQXQ_powers(x, l, S, T, p): FpXQ_powers(x, l, S, p);
528 : }
529 :
530 : GEN
531 0 : FqXQ_matrix_pow(GEN y, long n, long m, GEN S, GEN T, GEN p)
532 : {
533 0 : return T ? FpXQXQ_matrix_pow(y, n, m, S, T, p): FpXQ_matrix_pow(y, n, m, S, p);
534 : }
535 :
536 : /*******************************************************************/
537 : /* */
538 : /* Fq */
539 : /* */
540 : /*******************************************************************/
541 :
542 : GEN
543 11592874 : Fq_add(GEN x, GEN y, GEN T/*unused*/, GEN p)
544 : {
545 : (void)T;
546 11592874 : switch((typ(x)==t_POL)|((typ(y)==t_POL)<<1))
547 : {
548 1143687 : case 0: return Fp_add(x,y,p);
549 748136 : case 1: return FpX_Fp_add(x,y,p);
550 91991 : case 2: return FpX_Fp_add(y,x,p);
551 9609060 : case 3: return FpX_add(x,y,p);
552 : }
553 : return NULL;/*LCOV_EXCL_LINE*/
554 : }
555 :
556 : GEN
557 8562688 : Fq_sub(GEN x, GEN y, GEN T/*unused*/, GEN p)
558 : {
559 : (void)T;
560 8562688 : switch((typ(x)==t_POL)|((typ(y)==t_POL)<<1))
561 : {
562 255995 : case 0: return Fp_sub(x,y,p);
563 24540 : case 1: return FpX_Fp_sub(x,y,p);
564 23896 : case 2: return Fp_FpX_sub(x,y,p);
565 8258257 : case 3: return FpX_sub(x,y,p);
566 : }
567 : return NULL;/*LCOV_EXCL_LINE*/
568 : }
569 :
570 : GEN
571 1079556 : Fq_neg(GEN x, GEN T/*unused*/, GEN p)
572 : {
573 : (void)T;
574 1079556 : return (typ(x)==t_POL)? FpX_neg(x,p): Fp_neg(x,p);
575 : }
576 :
577 : GEN
578 81354 : Fq_halve(GEN x, GEN T/*unused*/, GEN p)
579 : {
580 : (void)T;
581 81354 : return (typ(x)==t_POL)? FpX_halve(x,p): Fp_halve(x,p);
582 : }
583 :
584 : /* If T==NULL do not reduce*/
585 : GEN
586 8608648 : Fq_mul(GEN x, GEN y, GEN T, GEN p)
587 : {
588 8608648 : switch((typ(x)==t_POL)|((typ(y)==t_POL)<<1))
589 : {
590 1037917 : case 0: return Fp_mul(x,y,p);
591 128565 : case 1: return FpX_Fp_mul(x,y,p);
592 394686 : case 2: return FpX_Fp_mul(y,x,p);
593 7047481 : case 3: if (T) return FpXQ_mul(x,y,T,p);
594 4478770 : else return FpX_mul(x,y,p);
595 : }
596 : return NULL;/*LCOV_EXCL_LINE*/
597 : }
598 :
599 : /* If T==NULL do not reduce*/
600 : GEN
601 490555 : Fq_mulu(GEN x, ulong y, /*unused*/GEN T, GEN p)
602 : {
603 : (void) T;
604 490555 : return typ(x)==t_POL ? FpX_Fp_mul(x,utoi(y),p): Fp_mulu(x, y, p);
605 : }
606 :
607 : /* y t_INT */
608 : GEN
609 43929 : Fq_Fp_mul(GEN x, GEN y, GEN T/*unused*/, GEN p)
610 : {
611 : (void)T;
612 6822 : return (typ(x) == t_POL)? FpX_Fp_mul(x,y,p)
613 50751 : : Fp_mul(x,y,p);
614 : }
615 : /* If T==NULL do not reduce*/
616 : GEN
617 611173 : Fq_sqr(GEN x, GEN T, GEN p)
618 : {
619 611173 : if (typ(x) == t_POL)
620 : {
621 70585 : if (T) return FpXQ_sqr(x,T,p);
622 0 : else return FpX_sqr(x,p);
623 : }
624 : else
625 540588 : return Fp_sqr(x,p);
626 : }
627 :
628 : GEN
629 0 : Fq_neg_inv(GEN x, GEN T, GEN p)
630 : {
631 0 : if (typ(x) == t_INT) return Fp_inv(Fp_neg(x,p),p);
632 0 : return FpXQ_inv(FpX_neg(x,p),T,p);
633 : }
634 :
635 : GEN
636 0 : Fq_invsafe(GEN x, GEN pol, GEN p)
637 : {
638 0 : if (typ(x) == t_INT) return Fp_invsafe(x,p);
639 0 : return FpXQ_invsafe(x,pol,p);
640 : }
641 :
642 : GEN
643 89311 : Fq_inv(GEN x, GEN pol, GEN p)
644 : {
645 89311 : if (typ(x) == t_INT) return Fp_inv(x,p);
646 81545 : return FpXQ_inv(x,pol,p);
647 : }
648 :
649 : GEN
650 343791 : Fq_div(GEN x, GEN y, GEN pol, GEN p)
651 : {
652 343791 : switch((typ(x)==t_POL)|((typ(y)==t_POL)<<1))
653 : {
654 318402 : case 0: return Fp_div(x,y,p);
655 16702 : case 1: return FpX_Fp_div(x,y,p);
656 406 : case 2: return FpX_Fp_mul(FpXQ_inv(y,pol,p),x,p);
657 8281 : case 3: return FpXQ_div(x,y,pol,p);
658 : }
659 : return NULL;/*LCOV_EXCL_LINE*/
660 : }
661 :
662 : GEN
663 795544 : Fq_pow(GEN x, GEN n, GEN pol, GEN p)
664 : {
665 795544 : if (typ(x) == t_INT) return Fp_pow(x,n,p);
666 136711 : return FpXQ_pow(x,n,pol,p);
667 : }
668 :
669 : GEN
670 15050 : Fq_powu(GEN x, ulong n, GEN pol, GEN p)
671 : {
672 15050 : if (typ(x) == t_INT) return Fp_powu(x,n,p);
673 1267 : return FpXQ_powu(x,n,pol,p);
674 : }
675 :
676 : GEN
677 1892926 : Fq_sqrt(GEN x, GEN T, GEN p)
678 : {
679 1892926 : if (typ(x) == t_INT)
680 : {
681 1825064 : if (!T || odd(get_FpX_degree(T))) return Fp_sqrt(x,p);
682 9621 : x = scalarpol_shallow(x, get_FpX_var(T));
683 : }
684 77483 : return FpXQ_sqrt(x,T,p);
685 : }
686 : GEN
687 170786 : Fq_sqrtn(GEN x, GEN n, GEN T, GEN p, GEN *zeta)
688 : {
689 170786 : if (typ(x) == t_INT)
690 : {
691 : long d;
692 170415 : if (!T) return Fp_sqrtn(x,n,p,zeta);
693 126 : d = get_FpX_degree(T);
694 126 : if (ugcdiu(n,d) == 1)
695 : {
696 56 : if (!zeta) return Fp_sqrtn(x,n,p,NULL);
697 : /* gcd(n,p-1)=gcd(n,q-1): same number of solutions in Fp and F_q */
698 21 : if (equalii(gcdii(subiu(p,1),n), gcdii(subiu(Fp_powu(p,d,n), 1), n)))
699 14 : return Fp_sqrtn(x,n,p,zeta);
700 : }
701 77 : x = scalarpol(x, get_FpX_var(T)); /* left on stack */
702 : }
703 448 : return FpXQ_sqrtn(x,n,T,p,zeta);
704 : }
705 :
706 : struct _Fq_field
707 : {
708 : GEN T, p;
709 : };
710 :
711 : static GEN
712 303247 : _Fq_red(void *E, GEN x)
713 303247 : { struct _Fq_field *s = (struct _Fq_field *)E;
714 303247 : return Fq_red(x, s->T, s->p);
715 : }
716 :
717 : static GEN
718 362523 : _Fq_add(void *E, GEN x, GEN y)
719 : {
720 : (void) E;
721 362523 : switch((typ(x)==t_POL)|((typ(y)==t_POL)<<1))
722 : {
723 14 : case 0: return addii(x,y);
724 0 : case 1: return ZX_Z_add(x,y);
725 15918 : case 2: return ZX_Z_add(y,x);
726 346591 : default: return ZX_add(x,y);
727 : }
728 : }
729 :
730 : static GEN
731 315028 : _Fq_neg(void *E, GEN x) { (void) E; return typ(x)==t_POL?ZX_neg(x):negi(x); }
732 :
733 : static GEN
734 243341 : _Fq_mul(void *E, GEN x, GEN y)
735 : {
736 : (void) E;
737 243341 : switch((typ(x)==t_POL)|((typ(y)==t_POL)<<1))
738 : {
739 679 : case 0: return mulii(x,y);
740 1085 : case 1: return ZX_Z_mul(x,y);
741 1043 : case 2: return ZX_Z_mul(y,x);
742 240534 : default: return ZX_mul(x,y);
743 : }
744 : }
745 :
746 : static GEN
747 9331 : _Fq_inv(void *E, GEN x)
748 9331 : { struct _Fq_field *s = (struct _Fq_field *)E;
749 9331 : return Fq_inv(x,s->T,s->p);
750 : }
751 :
752 : static int
753 38388 : _Fq_equal0(GEN x) { return signe(x)==0; }
754 :
755 : static GEN
756 4151 : _Fq_s(void *E, long x) { (void) E; return stoi(x); }
757 :
758 : static const struct bb_field Fq_field={_Fq_red,_Fq_add,_Fq_mul,_Fq_neg,
759 : _Fq_inv,_Fq_equal0,_Fq_s};
760 :
761 4725 : const struct bb_field *get_Fq_field(void **E, GEN T, GEN p)
762 : {
763 4725 : if (!T)
764 0 : return get_Fp_field(E, p);
765 : else
766 : {
767 4725 : GEN z = new_chunk(sizeof(struct _Fq_field));
768 4725 : struct _Fq_field *e = (struct _Fq_field *) z;
769 4725 : e->T = T; e->p = p; *E = (void*)e;
770 4725 : return &Fq_field;
771 : }
772 : }
773 :
774 : /*******************************************************************/
775 : /* */
776 : /* Fq[X] */
777 : /* */
778 : /*******************************************************************/
779 : /* P(X + c) */
780 : static GEN
781 434 : Fp_XpN_powu(GEN u, long n, GEN p, long v)
782 : {
783 : pari_sp av;
784 : long k;
785 : GEN B, C, V;
786 434 : if (!n) return pol_1(v);
787 434 : if (is_pm1(u))
788 434 : return Xpm1_powu(n, signe(u), v);
789 0 : av = avma;
790 0 : V = Fp_powers(u, n, p);
791 0 : B = FpV_red(vecbinomial(n), p);
792 0 : C = cgetg(n+3, t_POL);
793 0 : C[1] = evalsigne(1)| evalvarn(v);
794 0 : for (k=1; k <= n+1; k++)
795 0 : gel(C,k+1) = Fp_mul(gel(V,n+2-k), gel(B,k), p);
796 0 : return gc_upto(av, C);
797 : }
798 :
799 : static GEN
800 700 : FpX_translate_basecase(GEN P, GEN c, GEN p)
801 : {
802 700 : pari_sp av = avma;
803 : GEN Q, *R;
804 : long i, k, n;
805 :
806 700 : if (!signe(P) || !signe(c)) return ZX_copy(P);
807 560 : Q = leafcopy(P);
808 560 : R = (GEN*)(Q+2); n = degpol(P);
809 1316 : for (i=1; i<=n; i++)
810 : {
811 2016 : for (k=n-i; k<n; k++)
812 1260 : R[k] = Fp_add(R[k], Fp_mul(c, R[k+1], p), p);
813 :
814 756 : if (gc_needed(av,2))
815 : {
816 0 : if(DEBUGMEM>1) pari_warn(warnmem,"FpX_translate, i = %ld/%ld", i,n);
817 0 : Q = gc_GEN(av, Q); R = (GEN*)Q+2;
818 : }
819 : }
820 560 : return gc_GEN(av, FpX_renormalize(Q, lg(Q)));
821 : }
822 :
823 : GEN
824 1134 : FpX_translate(GEN P, GEN c, GEN p)
825 : {
826 1134 : pari_sp av = avma;
827 1134 : long n = degpol(P);
828 1134 : if (n<=3 || 25*(n-3) < expi(p))
829 700 : return FpX_translate_basecase(P, c, p);
830 : else
831 : {
832 434 : long d = n >> 1;
833 434 : GEN Q = FpX_translate(RgX_shift_shallow(P, -d), c, p);
834 434 : GEN R = FpX_translate(RgXn_red_shallow(P, d), c, p);
835 434 : GEN S = Fp_XpN_powu(c, d, p, varn(P));
836 434 : return gc_upto(av, FpX_add(FpX_mul(Q, S, p), R, p));
837 : }
838 : }
839 :
840 : /* P(X + c), c an Fq */
841 : GEN
842 33880 : FqX_translate(GEN P, GEN c, GEN T, GEN p)
843 : {
844 33880 : pari_sp av = avma;
845 : GEN Q, *R;
846 : long i, k, n;
847 :
848 : /* signe works for t_(INT|POL) */
849 33880 : if (!signe(P) || !signe(c)) return RgX_copy(P);
850 33880 : Q = leafcopy(P);
851 33880 : R = (GEN*)(Q+2); n = degpol(P);
852 150059 : for (i=1; i<=n; i++)
853 : {
854 433559 : for (k=n-i; k<n; k++)
855 317380 : R[k] = Fq_add(R[k], Fq_mul(c, R[k+1], T, p), T, p);
856 :
857 116179 : if (gc_needed(av,2))
858 : {
859 0 : if(DEBUGMEM>1) pari_warn(warnmem,"FqX_translate, i = %ld/%ld", i,n);
860 0 : Q = gc_GEN(av, Q); R = (GEN*)Q+2;
861 : }
862 : }
863 33880 : return gc_GEN(av, FpXQX_renormalize(Q, lg(Q)));
864 : }
865 :
866 : GEN
867 40452 : FqV_roots_to_pol(GEN V, GEN T, GEN p, long v)
868 : {
869 40452 : pari_sp ltop = avma;
870 : long k;
871 : GEN W;
872 40452 : if (lgefint(p) == 3)
873 : {
874 31741 : ulong pp = p[2];
875 31741 : GEN Tl = ZX_to_Flx(T, pp);
876 31744 : GEN Vl = FqC_to_FlxqC(V, Tl, pp);
877 31745 : Tl = FlxqV_roots_to_pol(Vl, Tl, pp, v);
878 31746 : return gc_upto(ltop, FlxX_to_ZXX(Tl));
879 : }
880 8711 : W = cgetg(lg(V),t_VEC);
881 78117 : for(k=1; k < lg(V); k++)
882 69406 : gel(W,k) = deg1pol_shallow(gen_1,Fq_neg(gel(V,k),T,p),v);
883 8711 : return gc_upto(ltop, FpXQXV_prod(W, T, p));
884 : }
885 :
886 : GEN
887 109459 : FqV_red(GEN x, GEN T, GEN p)
888 778238 : { pari_APPLY_type(t_VEC, Fq_red(gel(x,i), T, p)) }
889 :
890 : GEN
891 24058 : FqC_red(GEN x, GEN T, GEN p)
892 163384 : { pari_APPLY_type(t_COL, Fq_red(gel(x,i), T, p)) }
893 :
894 : GEN
895 1701 : FqM_red(GEN x, GEN T, GEN p)
896 5411 : { pari_APPLY_same(FqC_red(gel(x,i), T, p)) }
897 :
898 : GEN
899 0 : FqC_add(GEN x, GEN y, GEN T, GEN p)
900 : {
901 0 : if (!T) return FpC_add(x, y, p);
902 0 : pari_APPLY_type(t_COL, Fq_add(gel(x,i), gel(y,i), T, p))
903 : }
904 :
905 : GEN
906 0 : FqC_sub(GEN x, GEN y, GEN T, GEN p)
907 : {
908 0 : if (!T) return FpC_sub(x, y, p);
909 0 : pari_APPLY_type(t_COL, Fq_sub(gel(x,i), gel(y,i), T, p))
910 : }
911 :
912 : GEN
913 0 : FqC_Fq_mul(GEN x, GEN y, GEN T, GEN p)
914 : {
915 0 : if (!T) return FpC_Fp_mul(x, y, p);
916 0 : pari_APPLY_type(t_COL, Fq_mul(gel(x,i),y,T,p))
917 : }
918 :
919 : GEN
920 105 : FqC_FqV_mul(GEN x, GEN y, GEN T, GEN p)
921 : {
922 105 : long i,j, lx=lg(x), ly=lg(y);
923 : GEN z;
924 105 : if (ly==1) return cgetg(1,t_MAT);
925 105 : z = cgetg(ly,t_MAT);
926 819 : for (j=1; j < ly; j++)
927 : {
928 714 : GEN zj = cgetg(lx,t_COL);
929 4200 : for (i=1; i<lx; i++) gel(zj,i) = Fq_mul(gel(x,i),gel(y,j), T, p);
930 714 : gel(z, j) = zj;
931 : }
932 105 : return z;
933 : }
934 :
935 : GEN
936 5467 : FpXC_center(GEN x, GEN p, GEN pov2)
937 41476 : { pari_APPLY_type(t_COL, FpX_center(gel(x,i), p, pov2)) }
938 :
939 : GEN
940 109021 : FqC_to_FlxqC(GEN x, GEN T, ulong p)
941 109021 : { long sv = get_Flx_var(T);
942 4834755 : pari_APPLY_type(t_COL,typ(gel(x,i))==t_INT ?
943 : Z_to_Flx(gel(x,i), p, sv): ZX_to_Flx(gel(x,i), p)) }
944 :
945 : GEN
946 8708 : FqM_to_FlxqM(GEN x, GEN T, ulong p)
947 85985 : { pari_APPLY_same(FqC_to_FlxqC(gel(x,i), T, p)) }
948 :
949 : GEN
950 1800 : FpXM_center(GEN x, GEN p, GEN pov2)
951 7267 : { pari_APPLY_same(FpXC_center(gel(x,i), p, pov2)) }
952 :
953 : /*******************************************************************/
954 : /* */
955 : /* GENERIC CRT */
956 : /* */
957 : /*******************************************************************/
958 : static GEN
959 8309555 : primelist(forprime_t *S, long n, GEN dB)
960 : {
961 8309555 : GEN P = cgetg(n+1, t_VECSMALL);
962 8309542 : long i = 1;
963 : ulong p;
964 20071949 : while (i <= n && (p = u_forprime_next(S)))
965 11762408 : if (!dB || umodiu(dB, p)) P[i++] = p;
966 8309540 : return P;
967 : }
968 :
969 : void
970 7726980 : gen_inccrt_i(const char *str, GEN worker, GEN dB, long n, long mmin,
971 : forprime_t *S, GEN *pH, GEN *pmod, GEN crt(GEN, GEN, GEN*),
972 : GEN center(GEN, GEN, GEN))
973 : {
974 7726980 : long m = mmin? minss(mmin, n): usqrt(n);
975 : GEN H, P, mod;
976 : pari_timer ti;
977 7726977 : if (DEBUGLEVEL > 4)
978 : {
979 0 : timer_start(&ti);
980 0 : err_printf("%s: nb primes: %ld\n",str, n);
981 : }
982 7726976 : if (m == 1)
983 : {
984 7415914 : GEN P = primelist(S, n, dB);
985 7415894 : GEN done = closure_callgen1(worker, P);
986 7415890 : H = gel(done,1);
987 7415890 : mod = gel(done,2);
988 7415890 : if (!*pH && center) H = center(H, mod, shifti(mod,-1));
989 7415828 : if (DEBUGLEVEL>4) timer_printf(&ti,"%s: modular", str);
990 : }
991 : else
992 : {
993 311062 : long i, s = (n+m-1)/m, r = m - (m*s-n), di = 0;
994 : struct pari_mt pt;
995 311062 : long pending = 0;
996 311062 : H = cgetg(m+1, t_VEC); P = cgetg(m+1, t_VEC);
997 311062 : mt_queue_start_lim(&pt, worker, m);
998 1270306 : for (i=1; i<=m || pending; i++)
999 : {
1000 : GEN done;
1001 959242 : GEN pr = i <= m ? mkvec(primelist(S, i<=r ? s: s-1, dB)): NULL;
1002 959244 : mt_queue_submit(&pt, i, pr);
1003 959245 : done = mt_queue_get(&pt, NULL, &pending);
1004 959245 : if (done)
1005 : {
1006 893645 : di++;
1007 893645 : gel(H, di) = gel(done,1);
1008 893645 : gel(P, di) = gel(done,2);
1009 893645 : if (DEBUGLEVEL>5) err_printf("%ld%% ",100*di/m);
1010 : }
1011 : }
1012 311064 : mt_queue_end(&pt);
1013 311064 : if (DEBUGLEVEL>5) err_printf("\n");
1014 311064 : if (DEBUGLEVEL>4) timer_printf(&ti,"%s: modular", str);
1015 311064 : H = crt(H, P, &mod);
1016 311064 : if (DEBUGLEVEL>4) timer_printf(&ti,"%s: chinese", str);
1017 : }
1018 7726892 : if (*pH) H = crt(mkvec2(*pH, H), mkvec2(*pmod, mod), &mod);
1019 7726892 : *pH = H; *pmod = mod;
1020 7726892 : }
1021 : void
1022 2060066 : gen_inccrt(const char *str, GEN worker, GEN dB, long n, long mmin,
1023 : forprime_t *S, GEN *pH, GEN *pmod, GEN crt(GEN, GEN, GEN*),
1024 : GEN center(GEN, GEN, GEN))
1025 : {
1026 2060066 : pari_sp av = avma;
1027 2060066 : gen_inccrt_i(str, worker, dB, n, mmin, S, pH, pmod, crt, center);
1028 2060015 : (void)gc_all(av, 2, pH, pmod);
1029 2060168 : }
1030 :
1031 : GEN
1032 1273428 : gen_crt(const char *str, GEN worker, forprime_t *S, GEN dB, ulong bound, long mmin, GEN *pmod,
1033 : GEN crt(GEN, GEN, GEN*), GEN center(GEN, GEN, GEN))
1034 : {
1035 1273428 : GEN mod = gen_1, H = NULL;
1036 : ulong e;
1037 :
1038 1273428 : bound++;
1039 2546946 : while (bound > (e = expi(mod)))
1040 : {
1041 1273390 : long n = (bound - e) / expu(S->p) + 1;
1042 1273415 : gen_inccrt(str, worker, dB, n, mmin, S, &H, &mod, crt, center);
1043 : }
1044 1273501 : if (pmod) *pmod = mod;
1045 1273501 : return H;
1046 : }
1047 :
1048 : /*******************************************************************/
1049 : /* */
1050 : /* MODULAR GCD */
1051 : /* */
1052 : /*******************************************************************/
1053 : /* return z = a mod q, b mod p (p,q) = 1; qinv = 1/q mod p; a in ]-q,q] */
1054 : static GEN
1055 5155154 : Fl_chinese_coprime(GEN a, ulong b, GEN q, ulong p, ulong qinv, GEN pq, GEN pq2)
1056 : {
1057 5155154 : ulong d, amod = umodiu(a, p);
1058 5155174 : pari_sp av = avma;
1059 : GEN ax;
1060 :
1061 5155174 : if (b == amod) return NULL;
1062 2126315 : d = Fl_mul(Fl_sub(b, amod, p), qinv, p); /* != 0 */
1063 2126937 : if (d >= 1 + (p>>1))
1064 1037825 : ax = subii(a, mului(p-d, q));
1065 : else
1066 : {
1067 1089112 : ax = addii(a, mului(d, q)); /* in ]0, pq[ assuming a in ]-q,q[ */
1068 1088663 : if (cmpii(ax,pq2) > 0) ax = subii(ax,pq);
1069 : }
1070 2126024 : return gc_INT(av, ax);
1071 : }
1072 : GEN
1073 364 : Z_init_CRT(ulong Hp, ulong p) { return stoi(Fl_center(Hp, p, p>>1)); }
1074 : GEN
1075 31794 : ZX_init_CRT(GEN Hp, ulong p, long v)
1076 : {
1077 31794 : long i, l = lg(Hp), lim = (long)(p>>1);
1078 31794 : GEN H = cgetg(l, t_POL);
1079 31794 : H[1] = evalsigne(1) | evalvarn(v);
1080 796087 : for (i=2; i<l; i++)
1081 764294 : gel(H,i) = stoi(Fl_center(Hp[i], p, lim));
1082 31793 : return ZX_renormalize(H,l);
1083 : }
1084 :
1085 : GEN
1086 5789 : ZM_init_CRT(GEN Hp, ulong p)
1087 : {
1088 5789 : long i,j, m, l = lg(Hp), lim = (long)(p>>1);
1089 5789 : GEN c, cp, H = cgetg(l, t_MAT);
1090 5789 : if (l==1) return H;
1091 5789 : m = lgcols(Hp);
1092 18984 : for (j=1; j<l; j++)
1093 : {
1094 13195 : cp = gel(Hp,j);
1095 13195 : c = cgetg(m, t_COL);
1096 13195 : gel(H,j) = c;
1097 166166 : for (i=1; i<m; i++) gel(c,i) = stoi(Fl_center(cp[i],p, lim));
1098 : }
1099 5789 : return H;
1100 : }
1101 :
1102 : int
1103 7560 : Z_incremental_CRT(GEN *H, ulong Hp, GEN *ptq, ulong p)
1104 : {
1105 7560 : GEN h, q = *ptq, qp = muliu(q,p);
1106 7560 : ulong qinv = Fl_inv(umodiu(q,p), p);
1107 7560 : int stable = 1;
1108 7560 : h = Fl_chinese_coprime(*H,Hp,q,p,qinv,qp,shifti(qp,-1));
1109 7560 : if (h) { *H = h; stable = 0; }
1110 7560 : *ptq = qp; return stable;
1111 : }
1112 :
1113 : static int
1114 147611 : ZX_incremental_CRT_raw(GEN *ptH, GEN Hp, GEN q, GEN qp, ulong p)
1115 : {
1116 147611 : GEN H = *ptH, h, qp2 = shifti(qp,-1);
1117 147604 : ulong qinv = Fl_inv(umodiu(q,p), p);
1118 147612 : long i, l = lg(H), lp = lg(Hp);
1119 147612 : int stable = 1;
1120 :
1121 147612 : if (l < lp)
1122 : { /* degree increases */
1123 0 : GEN x = cgetg(lp, t_POL);
1124 0 : for (i=1; i<l; i++) x[i] = H[i];
1125 0 : for ( ; i<lp; i++) gel(x,i) = gen_0;
1126 0 : *ptH = H = x;
1127 0 : stable = 0;
1128 147612 : } else if (l > lp)
1129 : { /* degree decreases */
1130 0 : GEN x = cgetg(l, t_VECSMALL);
1131 0 : for (i=1; i<lp; i++) x[i] = Hp[i];
1132 0 : for ( ; i<l; i++) x[i] = 0;
1133 0 : Hp = x; lp = l;
1134 : }
1135 4933609 : for (i=2; i<lp; i++)
1136 : {
1137 4786108 : h = Fl_chinese_coprime(gel(H,i),Hp[i],q,p,qinv,qp,qp2);
1138 4785997 : if (h) { gel(H,i) = h; stable = 0; }
1139 : }
1140 147501 : (void)ZX_renormalize(H,lp);
1141 147615 : return stable;
1142 : }
1143 :
1144 : int
1145 0 : ZX_incremental_CRT(GEN *ptH, GEN Hp, GEN *ptq, ulong p)
1146 : {
1147 0 : GEN q = *ptq, qp = muliu(q,p);
1148 0 : int stable = ZX_incremental_CRT_raw(ptH, Hp, q, qp, p);
1149 0 : *ptq = qp; return stable;
1150 : }
1151 :
1152 : int
1153 7597 : ZM_incremental_CRT(GEN *pH, GEN Hp, GEN *ptq, ulong p)
1154 : {
1155 7597 : GEN h, H = *pH, q = *ptq, qp = muliu(q, p), qp2 = shifti(qp,-1);
1156 7597 : ulong qinv = Fl_inv(umodiu(q,p), p);
1157 7597 : long i,j, l = lg(H), m = lgcols(H);
1158 7597 : int stable = 1;
1159 26206 : for (j=1; j<l; j++)
1160 202351 : for (i=1; i<m; i++)
1161 : {
1162 183742 : h = Fl_chinese_coprime(gcoeff(H,i,j), coeff(Hp,i,j),q,p,qinv,qp,qp2);
1163 183742 : if (h) { gcoeff(H,i,j) = h; stable = 0; }
1164 : }
1165 7597 : *ptq = qp; return stable;
1166 : }
1167 :
1168 : GEN
1169 679 : ZXM_init_CRT(GEN Hp, long deg, ulong p)
1170 : {
1171 : long i, j, k;
1172 : GEN H;
1173 679 : long m, l = lg(Hp), lim = (long)(p>>1), n;
1174 679 : H = cgetg(l, t_MAT);
1175 679 : if (l==1) return H;
1176 679 : m = lgcols(Hp);
1177 679 : n = deg + 3;
1178 2268 : for (j=1; j<l; j++)
1179 : {
1180 1589 : GEN cp = gel(Hp,j);
1181 1589 : GEN c = cgetg(m, t_COL);
1182 1589 : gel(H,j) = c;
1183 24465 : for (i=1; i<m; i++)
1184 : {
1185 22876 : GEN dp = gel(cp, i);
1186 22876 : long l = lg(dp);
1187 22876 : GEN d = cgetg(n, t_POL);
1188 22876 : gel(c, i) = d;
1189 22876 : d[1] = dp[1] | evalsigne(1);
1190 46459 : for (k=2; k<l; k++)
1191 23583 : gel(d,k) = stoi(Fl_center(dp[k], p, lim));
1192 45493 : for ( ; k<n; k++)
1193 22617 : gel(d,k) = gen_0;
1194 : }
1195 : }
1196 679 : return H;
1197 : }
1198 :
1199 : int
1200 653 : ZXM_incremental_CRT(GEN *pH, GEN Hp, GEN *ptq, ulong p)
1201 : {
1202 653 : GEN v, H = *pH, q = *ptq, qp = muliu(q, p), qp2 = shifti(qp,-1);
1203 653 : ulong qinv = Fl_inv(umodiu(q,p), p);
1204 653 : long i,j,k, l = lg(H), m = lgcols(H), n = lg(gmael(H,1,1));
1205 653 : int stable = 1;
1206 2225 : for (j=1; j<l; j++)
1207 90418 : for (i=1; i<m; i++)
1208 : {
1209 88846 : GEN h = gmael(H,j,i), hp = gmael(Hp,j,i);
1210 88846 : long lh = lg(hp);
1211 246641 : for (k=2; k<lh; k++)
1212 : {
1213 157795 : v = Fl_chinese_coprime(gel(h,k),uel(hp,k),q,p,qinv,qp,qp2);
1214 157795 : if (v) { gel(h,k) = v; stable = 0; }
1215 : }
1216 108763 : for (; k<n; k++)
1217 : {
1218 19917 : v = Fl_chinese_coprime(gel(h,k),0,q,p,qinv,qp,qp2);
1219 19917 : if (v) { gel(h,k) = v; stable = 0; }
1220 : }
1221 : }
1222 653 : *ptq = qp; return stable;
1223 : }
1224 :
1225 : /* record the degrees of Euclidean remainders (make them as large as
1226 : * possible : smaller values correspond to a degenerate sequence) */
1227 : static void
1228 23286 : Flx_resultant_set_dglist(GEN a, GEN b, GEN dglist, ulong p)
1229 : {
1230 : long da,db,dc, ind;
1231 23286 : pari_sp av = avma;
1232 :
1233 23286 : if (lgpol(a)==0 || lgpol(b)==0) return;
1234 22019 : da = degpol(a);
1235 22019 : db = degpol(b);
1236 22019 : if (db > da)
1237 0 : { swapspec(a,b, da,db); }
1238 22019 : else if (!da) return;
1239 22019 : ind = 0;
1240 144374 : while (db)
1241 : {
1242 122360 : GEN c = Flx_rem(a,b, p);
1243 122356 : a = b; b = c; dc = degpol(c);
1244 122355 : if (dc < 0) break;
1245 :
1246 122355 : ind++;
1247 122355 : if (dc > dglist[ind]) dglist[ind] = dc;
1248 122355 : if (gc_needed(av,2))
1249 : {
1250 0 : if (DEBUGMEM>1) pari_warn(warnmem,"Flx_resultant_all");
1251 0 : (void)gc_all(av, 2, &a,&b);
1252 : }
1253 122355 : db = dc; /* = degpol(b) */
1254 : }
1255 22014 : if (ind+1 > lg(dglist)) setlg(dglist,ind+1);
1256 22014 : set_avma(av);
1257 : }
1258 : /* assuming the PRS finishes on a degree 1 polynomial C0 + C1X, with
1259 : * "generic" degree sequence as given by dglist, set *Ci and return
1260 : * resultant(a,b). Modular version of Collins's subresultant */
1261 : static ulong
1262 2085119 : Flx_resultant_all(GEN a, GEN b, long *C0, long *C1, GEN dglist, ulong p)
1263 : {
1264 : long da,db,dc, ind;
1265 2085119 : ulong lb, res, g = 1UL, h = 1UL, ca = 1UL, cb = 1UL;
1266 2085119 : int s = 1;
1267 2085119 : pari_sp av = avma;
1268 :
1269 2085119 : *C0 = 1; *C1 = 0;
1270 2085119 : if (lgpol(a)==0 || lgpol(b)==0) return 0;
1271 2075691 : da = degpol(a);
1272 2075719 : db = degpol(b);
1273 2075696 : if (db > da)
1274 : {
1275 0 : swapspec(a,b, da,db);
1276 0 : if (both_odd(da,db)) s = -s;
1277 : }
1278 2075696 : else if (!da) return 1; /* = a[2] ^ db, since 0 <= db <= da = 0 */
1279 2075696 : ind = 0;
1280 19802961 : while (db)
1281 : { /* sub-resultant algo., applied to ca * a and cb * b, ca,cb scalars,
1282 : * da = deg a, db = deg b */
1283 17731948 : GEN c = Flx_rem(a,b, p);
1284 17618931 : long delta = da - db;
1285 :
1286 17618931 : if (both_odd(da,db)) s = -s;
1287 17619125 : lb = Fl_mul(b[db+2], cb, p);
1288 17636554 : a = b; b = c; dc = degpol(c);
1289 17636294 : ind++;
1290 17636294 : if (dc != dglist[ind]) return gc_ulong(av,0); /* degenerates */
1291 17631345 : if (g == h)
1292 : { /* frequent */
1293 17571501 : ulong cc = Fl_mul(ca, Fl_powu(Fl_div(lb,g,p), delta+1, p), p);
1294 17668094 : ca = cb;
1295 17668094 : cb = cc;
1296 : }
1297 : else
1298 : {
1299 59844 : ulong cc = Fl_mul(ca, Fl_powu(lb, delta+1, p), p);
1300 59844 : ulong ghdelta = Fl_mul(g, Fl_powu(h, delta, p), p);
1301 59844 : ca = cb;
1302 59844 : cb = Fl_div(cc, ghdelta, p);
1303 : }
1304 17728723 : da = db; /* = degpol(a) */
1305 17728723 : db = dc; /* = degpol(b) */
1306 :
1307 17728723 : g = lb;
1308 17728723 : if (delta == 1)
1309 17628543 : h = g; /* frequent */
1310 : else
1311 100180 : h = Fl_mul(h, Fl_powu(Fl_div(g,h,p), delta, p), p);
1312 :
1313 17727435 : if (gc_needed(av,2))
1314 : {
1315 0 : if (DEBUGMEM>1) pari_warn(warnmem,"Flx_resultant_all");
1316 0 : (void)gc_all(av, 2, &a,&b);
1317 : }
1318 : }
1319 2071013 : if (da > 1) return 0; /* Failure */
1320 : /* last nonconstant polynomial has degree 1 */
1321 2071013 : *C0 = Fl_mul(ca, a[2], p);
1322 2070979 : *C1 = Fl_mul(ca, a[3], p);
1323 2070964 : res = Fl_mul(cb, b[2], p);
1324 2070933 : if (s == -1) res = p - res;
1325 2070933 : return gc_ulong(av,res);
1326 : }
1327 :
1328 : /* Q a vector of polynomials representing B in Fp[X][Y], evaluate at X = x,
1329 : * Return 0 in case of degree drop. */
1330 : static GEN
1331 2108526 : FlxY_evalx_drop(GEN Q, ulong x, ulong p)
1332 : {
1333 : GEN z;
1334 2108526 : long i, lb = lg(Q);
1335 2108526 : ulong leadz = Flx_eval(leading_coeff(Q), x, p);
1336 2108378 : long vs=mael(Q,2,1);
1337 2108378 : if (!leadz) return zero_Flx(vs);
1338 :
1339 2097718 : z = cgetg(lb, t_VECSMALL); z[1] = vs;
1340 20069120 : for (i=2; i<lb-1; i++) z[i] = Flx_eval(gel(Q,i), x, p);
1341 2095857 : z[i] = leadz; return z;
1342 : }
1343 :
1344 : GEN
1345 2072 : FpXY_FpXQ_evaly(GEN Q, GEN y, GEN T, GEN p, long vx)
1346 : {
1347 2072 : pari_sp av = avma;
1348 2072 : long i, lb = lg(Q);
1349 : GEN z;
1350 2072 : if (lb == 2) return pol_0(vx);
1351 2072 : z = gel(Q, lb-1);
1352 2072 : if (lb == 3 || !signe(y)) return typ(z)==t_INT? scalar_ZX(z, vx): ZX_copy(z);
1353 :
1354 2072 : if (typ(z) == t_INT) z = scalar_ZX_shallow(z, vx);
1355 48636 : for (i=lb-2; i>=2; i--)
1356 : {
1357 46564 : GEN c = gel(Q,i);
1358 46564 : z = FqX_Fq_mul(z, y, T, p);
1359 46564 : z = typ(c) == t_INT? FqX_Fq_add(z,c,T,p): FqX_add(z,c,T,p);
1360 : }
1361 2072 : return gc_upto(av, z);
1362 : }
1363 :
1364 : static GEN
1365 291737 : ZX_norml1(GEN x)
1366 : {
1367 291737 : long i, l = lg(x);
1368 : GEN s;
1369 :
1370 291737 : if (l == 2) return gen_0;
1371 199183 : s = gel(x, l-1); /* != 0 */
1372 697188 : for (i = l-2; i > 1; i--) {
1373 498015 : GEN xi = gel(x,i);
1374 498015 : if (!signe(xi)) continue;
1375 259400 : s = addii_sign(s,1, xi,1);
1376 : }
1377 199173 : return s;
1378 : }
1379 : /* x >= 0, y != 0, return x + |y| */
1380 : static GEN
1381 25552 : addii_abs(GEN x, GEN y)
1382 : {
1383 25552 : if (!signe(x)) return absi_shallow(y);
1384 16044 : return addii_sign(x,1, y,1);
1385 : }
1386 :
1387 : /* x a ZX, return sum_{i >= k} |x[i]| binomial(i, k) */
1388 : static GEN
1389 31648 : ZX_norml1_1(GEN x, long k)
1390 : {
1391 31648 : long i, d = degpol(x);
1392 : GEN s, C; /* = binomial(i, k) */
1393 :
1394 31648 : if (!d || k > d) return gen_0;
1395 31648 : s = absi_shallow(gel(x, k+2)); /* may be 0 */
1396 31647 : C = gen_1;
1397 68049 : for (i = k+1; i <= d; i++) {
1398 36407 : GEN xi = gel(x,i+2);
1399 36407 : if (k) C = diviuexact(muliu(C, i), i-k);
1400 36409 : if (signe(xi)) s = addii_abs(s, mulii(C, xi));
1401 : }
1402 31642 : return s;
1403 : }
1404 : /* x has non-negative real coefficients */
1405 : static GEN
1406 3283 : RgX_norml1_1(GEN x, long k)
1407 : {
1408 3283 : long i, d = degpol(x);
1409 : GEN s, C; /* = binomial(i, k) */
1410 :
1411 3283 : if (!d || k > d) return gen_0;
1412 3283 : s = gel(x, k+2); /* may be 0 */
1413 3283 : C = gen_1;
1414 9198 : for (i = k+1; i <= d; i++) {
1415 5915 : GEN xi = gel(x,i+2);
1416 5915 : if (k) C = diviuexact(muliu(C, i), i-k);
1417 5915 : if (!gequal0(xi)) s = gadd(s, gmul(C, xi));
1418 : }
1419 3283 : return s;
1420 : }
1421 :
1422 : /* N_2(A)^2 */
1423 : static GEN
1424 8571 : sqrN2(GEN A, long prec)
1425 : {
1426 8571 : pari_sp av = avma;
1427 8571 : long i, l = lg(A);
1428 8571 : GEN a = gen_0;
1429 41935 : for (i = 2; i < l; i++)
1430 : {
1431 33364 : a = gadd(a, gabs(gnorm(gel(A,i)), prec));
1432 33364 : if (gc_needed(av,1))
1433 : {
1434 0 : if(DEBUGMEM>1) pari_warn(warnmem,"RgX_RgXY_ResBound i = %ld",i);
1435 0 : a = gc_upto(av, a);
1436 : }
1437 : }
1438 8571 : return a;
1439 : }
1440 : /* Interpolate at roots of 1 and use Hadamard bound for univariate resultant:
1441 : * bound = N_2(A)^degpol B N_2(B)^degpol(A), where
1442 : * N_2(A) = sqrt(sum (N_1(Ai))^2)
1443 : * Return e such that Res(A, B) < 2^e */
1444 : static GEN
1445 7717 : RgX_RgXY_ResBound(GEN A, GEN B, long prec)
1446 : {
1447 7717 : pari_sp av = avma;
1448 7717 : GEN b = gen_0, bnd;
1449 7717 : long i, lB = lg(B);
1450 30403 : for (i=2; i<lB; i++)
1451 : {
1452 22686 : GEN t = gel(B,i);
1453 22686 : if (typ(t) == t_POL) t = gnorml1(t, prec);
1454 22686 : b = gadd(b, gabs(gsqr(t), prec));
1455 22686 : if (gc_needed(av,1))
1456 : {
1457 0 : if(DEBUGMEM>1) pari_warn(warnmem,"RgX_RgXY_ResBound i = %ld",i);
1458 0 : b = gc_upto(av, b);
1459 : }
1460 : }
1461 7717 : bnd = gsqrt(gmul(gpowgs(sqrN2(A,prec), degpol(B)),
1462 : gpowgs(b, degpol(A))), prec);
1463 7717 : return gc_upto(av, bnd);
1464 : }
1465 : /* A,B in C[X] return RgX_RgXY_ResBound(A, B(x+y)) */
1466 : static GEN
1467 854 : RgX_RgXY_ResBound_1(GEN A, GEN B, long prec)
1468 : {
1469 854 : pari_sp av = avma, av2;
1470 854 : GEN b = gen_0, bnd;
1471 854 : long i, lB = lg(B);
1472 854 : B = shallowcopy(B);
1473 4137 : for (i=2; i<lB; i++) gel(B,i) = gabs(gel(B,i), prec);
1474 854 : av2 = avma;
1475 4137 : for (i=2; i<lB; i++)
1476 : {
1477 3283 : b = gadd(b, gsqr(RgX_norml1_1(B, i-2)));
1478 3283 : if (gc_needed(av2,1))
1479 : {
1480 0 : if(DEBUGMEM>1) pari_warn(warnmem,"RgX_RgXY_ResBound i = %ld",i);
1481 0 : b = gc_upto(av2, b);
1482 : }
1483 : }
1484 854 : bnd = gsqrt(gmul(gpowgs(sqrN2(A,prec), degpol(B)),
1485 : gpowgs(b, degpol(A))), prec);
1486 854 : return gc_upto(av, bnd);
1487 : }
1488 :
1489 : /* log2 N_2(A)^2 */
1490 : static double
1491 176691 : log2N2(GEN A)
1492 : {
1493 176691 : pari_sp av = avma;
1494 176691 : long i, l = lg(A);
1495 176691 : GEN a = gen_0;
1496 1335190 : for (i=2; i < l; i++)
1497 : {
1498 1158495 : a = addii(a, sqri(gel(A,i)));
1499 1158501 : if (gc_needed(av,1))
1500 : {
1501 0 : if(DEBUGMEM>1) pari_warn(warnmem,"ZX_ZXY_ResBound i = %ld",i);
1502 0 : a = gc_upto(av, a);
1503 : }
1504 : }
1505 176695 : return gc_double(av, dbllog2(a));
1506 : }
1507 : /* Interpolate at roots of 1 and use Hadamard bound for univariate resultant:
1508 : * bound = N_2(A)^degpol B N_2(B)^degpol(A), where
1509 : * N_2(A) = sqrt(sum (N_1(Ai))^2)
1510 : * Return e such that Res(A, B) < 2^e */
1511 : ulong
1512 166612 : ZX_ZXY_ResBound(GEN A, GEN B, GEN dB)
1513 : {
1514 166612 : pari_sp av = avma;
1515 166612 : GEN b = gen_0;
1516 166612 : long i, lB = lg(B);
1517 : double logb;
1518 1260638 : for (i=2; i<lB; i++)
1519 : {
1520 1094036 : GEN t = gel(B,i);
1521 1094036 : if (typ(t) == t_POL) t = ZX_norml1(t);
1522 1094031 : b = addii(b, sqri(t));
1523 1094027 : if (gc_needed(av,1))
1524 : {
1525 0 : if(DEBUGMEM>1) pari_warn(warnmem,"ZX_ZXY_ResBound i = %ld",i);
1526 0 : b = gc_upto(av, b);
1527 : }
1528 : }
1529 166602 : logb = dbllog2(b); if (dB) logb -= 2 * dbllog2(dB);
1530 166607 : i = (long)((degpol(B) * log2N2(A) + degpol(A) * logb) / 2);
1531 166612 : return gc_ulong(av, (i <= 0)? 1: 1 + (ulong)i);
1532 : }
1533 : /* A,B ZX. Return ZX_ZXY_ResBound(A(x), B(x+y)) */
1534 : static ulong
1535 10084 : ZX_ZXY_ResBound_1(GEN A, GEN B)
1536 : {
1537 10084 : pari_sp av = avma;
1538 10084 : GEN b = gen_0;
1539 10084 : long i, lB = lg(B);
1540 41735 : for (i=2; i<lB; i++)
1541 : {
1542 31648 : b = addii(b, sqri(ZX_norml1_1(B, i-2)));
1543 31651 : if (gc_needed(av,1))
1544 : {
1545 0 : if(DEBUGMEM>1) pari_warn(warnmem,"ZX_ZXY_ResBound i = %ld",i);
1546 0 : b = gc_upto(av, b);
1547 : }
1548 : }
1549 10087 : i = (long)((degpol(B) * log2N2(A) + degpol(A) * dbllog2(b)) / 2);
1550 10086 : return gc_ulong(av, (i <= 0)? 1: 1 + (ulong)i);
1551 : }
1552 : /* special case B = A' */
1553 : static ulong
1554 1134083 : ZX_discbound(GEN A)
1555 : {
1556 1134083 : pari_sp av = avma;
1557 1134083 : GEN a = gen_0, b = gen_0;
1558 1134083 : long i , lA = lg(A), dA = degpol(A);
1559 : double loga, logb;
1560 6767141 : for (i = 2; i < lA; i++)
1561 : {
1562 5633220 : GEN c = sqri(gel(A,i));
1563 5632895 : a = addii(a, c);
1564 5633033 : if (i > 2) b = addii(b, mulii(c, sqru(i-2)));
1565 5633063 : if (gc_needed(av,1))
1566 : {
1567 0 : if(DEBUGMEM>1) pari_warn(warnmem,"ZX_discbound i = %ld",i);
1568 0 : (void)gc_all(av, 2, &a, &b);
1569 : }
1570 : }
1571 1133921 : loga = dbllog2(a);
1572 1133960 : logb = dbllog2(b); set_avma(av);
1573 1133981 : i = (long)(((dA-1) * loga + dA * logb) / 2);
1574 1133981 : return (i <= 0)? 1: 1 + (ulong)i;
1575 : }
1576 :
1577 : /* return Res(a(Y), b(n,Y)) over Fp. la = leading_coeff(a) [for efficiency] */
1578 : static ulong
1579 5538186 : Flx_FlxY_eval_resultant(GEN a, GEN b, ulong n, ulong p, ulong pi, ulong la)
1580 : {
1581 5538186 : GEN ev = FlxY_evalx_pre(b, n, p, pi);
1582 5538496 : long drop = lg(b) - lg(ev);
1583 5538496 : ulong r = Flx_resultant_pre(a, ev, p, pi);
1584 5537976 : if (drop && la != 1) r = Fl_mul(r, Fl_powu_pre(la, drop, p, pi), p);
1585 5538004 : return r;
1586 : }
1587 : static GEN
1588 284 : FpX_FpXY_eval_resultant(GEN a, GEN b, GEN n, GEN p, GEN la, long db, long vX)
1589 : {
1590 284 : GEN ev = FpXY_evaly(b, n, p, vX);
1591 284 : long drop = db-degpol(ev);
1592 284 : GEN r = FpX_resultant(a, ev, p);
1593 284 : if (drop && !gequal1(la)) r = Fp_mul(r, Fp_powu(la, drop,p),p);
1594 284 : return r;
1595 : }
1596 :
1597 : /* assume dres := deg(Res_X(a,b), Y) <= deg(a,X) * deg(b,Y) < p */
1598 : /* Return a Fly */
1599 : static GEN
1600 368413 : Flx_FlxY_resultant_polint(GEN a, GEN b, ulong p, ulong pi, long dres, long sx)
1601 : {
1602 : long i;
1603 368413 : ulong n, la = Flx_lead(a);
1604 368410 : GEN x = cgetg(dres+2, t_VECSMALL);
1605 368410 : GEN y = cgetg(dres+2, t_VECSMALL);
1606 : /* Evaluate at dres+ 1 points: 0 (if dres even) and +/- n, so that P_n(X) =
1607 : * P_{-n}(-X), where P_i is Lagrange polynomial: P_i(j) = delta_{i,j} */
1608 2956738 : for (i=0,n = 1; i < dres; n++)
1609 : {
1610 2588325 : x[++i] = n; y[i] = Flx_FlxY_eval_resultant(a,b, x[i], p,pi,la);
1611 2588263 : x[++i] = p-n; y[i] = Flx_FlxY_eval_resultant(a,b, x[i], p,pi,la);
1612 : }
1613 368413 : if (i == dres)
1614 : {
1615 362907 : x[++i] = 0; y[i] = Flx_FlxY_eval_resultant(a,b, x[i], p,pi,la);
1616 : }
1617 368410 : return Flv_polint(x,y, p, sx);
1618 : }
1619 :
1620 : static GEN
1621 7497 : FlxX_pseudorem(GEN x, GEN y, ulong p, ulong pi)
1622 : {
1623 7497 : long vx = varn(x), dx, dy, dz, i, lx, dp;
1624 7497 : pari_sp av = avma, av2;
1625 :
1626 7497 : if (!signe(y)) pari_err_INV("FlxX_pseudorem",y);
1627 7497 : (void)new_chunk(2);
1628 7498 : dx=degpol(x); x = RgX_recip_i(x)+2;
1629 7500 : dy=degpol(y); y = RgX_recip_i(y)+2; dz=dx-dy; dp = dz+1;
1630 7499 : av2 = avma;
1631 : for (;;)
1632 : {
1633 61435 : gel(x,0) = Flx_neg(gel(x,0), p); dp--;
1634 230121 : for (i=1; i<=dy; i++)
1635 168438 : gel(x,i) = Flx_add( Flx_mul_pre(gel(y,0), gel(x,i), p, pi),
1636 168660 : Flx_mul_pre(gel(x,0), gel(y,i), p, pi), p );
1637 1116835 : for ( ; i<=dx; i++)
1638 1056080 : gel(x,i) = Flx_mul_pre(gel(y,0), gel(x,i), p, pi);
1639 65327 : do { x++; dx--; } while (dx >= 0 && lg(gel(x,0))==2);
1640 60755 : if (dx < dy) break;
1641 53257 : if (gc_needed(av2,1))
1642 : {
1643 0 : if(DEBUGMEM>1) pari_warn(warnmem,"FlxX_pseudorem dx = %ld >= %ld",dx,dy);
1644 0 : gc_slice(av2,x,dx+1);
1645 : }
1646 : }
1647 7498 : if (dx < 0) return zero_Flx(0);
1648 7498 : lx = dx+3; x -= 2;
1649 7498 : x[0]=evaltyp(t_POL) | _evallg(lx);
1650 7498 : x[1]=evalsigne(1) | evalvarn(vx);
1651 7498 : x = RgX_recip_i(x);
1652 7497 : if (dp)
1653 : { /* multiply by y[0]^dp [beware dummy vars from FpX_FpXY_resultant] */
1654 1959 : GEN t = Flx_powu_pre(gel(y,0), dp, p, pi);
1655 7837 : for (i=2; i<lx; i++) gel(x,i) = Flx_mul_pre(gel(x,i), t, p, pi);
1656 : }
1657 7495 : return gc_GEN(av, x);
1658 : }
1659 :
1660 : /* return a Flx */
1661 : GEN
1662 2508 : FlxX_resultant(GEN u, GEN v, ulong p, long sx)
1663 : {
1664 2508 : pari_sp av = avma, av2;
1665 : long degq, dx, dy, du, dv, dr, signh;
1666 : ulong pi;
1667 : GEN z, g, h, r, p1;
1668 :
1669 2508 : dx = degpol(u); dy = degpol(v); signh = 1;
1670 2508 : if (dx < dy)
1671 : {
1672 7 : swap(u,v); lswap(dx,dy);
1673 7 : if (both_odd(dx, dy)) signh = -signh;
1674 : }
1675 2508 : if (dy < 0) return zero_Flx(sx);
1676 2508 : pi = SMALL_ULONG(p)? 0: get_Fl_red(p);
1677 2508 : if (dy==0) return gc_upto(av, Flx_powu_pre(gel(v,2),dx,p,pi));
1678 :
1679 2508 : g = h = pol1_Flx(sx); av2 = avma;
1680 : for(;;)
1681 : {
1682 7497 : r = FlxX_pseudorem(u,v,p,pi); dr = lg(r);
1683 7501 : if (dr == 2) { set_avma(av); return zero_Flx(sx); }
1684 7501 : du = degpol(u); dv = degpol(v); degq = du-dv;
1685 7501 : u = v; p1 = g; g = leading_coeff(u);
1686 7501 : switch(degq)
1687 : {
1688 0 : case 0: break;
1689 5528 : case 1:
1690 5528 : p1 = Flx_mul_pre(h,p1, p, pi); h = g; break;
1691 1973 : default:
1692 1973 : p1 = Flx_mul_pre(Flx_powu_pre(h,degq,p,pi), p1, p, pi);
1693 1973 : h = Flx_div_pre(Flx_powu_pre(g,degq,p,pi),
1694 1972 : Flx_powu_pre(h,degq-1,p,pi), p, pi);
1695 : }
1696 7497 : if (both_odd(du,dv)) signh = -signh;
1697 7496 : v = FlxY_Flx_div(r, p1, p);
1698 7496 : if (dr==3) break;
1699 4987 : if (gc_needed(av2,1))
1700 : {
1701 0 : if(DEBUGMEM>1) pari_warn(warnmem,"FlxX_resultant, dr = %ld",dr);
1702 0 : (void)gc_all(av2,4, &u, &v, &g, &h);
1703 : }
1704 : }
1705 2509 : z = gel(v,2);
1706 2509 : if (dv > 1) z = Flx_div_pre(Flx_powu_pre(z,dv,p,pi),
1707 0 : Flx_powu_pre(h,dv-1,p,pi), p, pi);
1708 2509 : if (signh < 0) z = Flx_neg(z,p);
1709 2509 : return gc_upto(av, z);
1710 : }
1711 :
1712 : /* Warning:
1713 : * This function switches between valid and invalid variable ordering*/
1714 :
1715 : static GEN
1716 6118 : FlxY_to_FlyX(GEN b, long sv)
1717 : {
1718 6118 : long i, n=-1;
1719 6118 : long sw = b[1]&VARNBITS;
1720 20869 : for(i=2;i<lg(b);i++) n = maxss(n,lgpol(gel(b,i)));
1721 6118 : return Flm_to_FlxX(Flm_transpose(FlxX_to_Flm(b,n)),sv,sw);
1722 : }
1723 :
1724 : /* Return a Fly*/
1725 : GEN
1726 6119 : Flx_FlxY_resultant(GEN a, GEN b, ulong p)
1727 : {
1728 6119 : pari_sp ltop=avma;
1729 6119 : long dres = degpol(a)*degpol(b);
1730 6118 : long sx=a[1], sy=b[1]&VARNBITS;
1731 : GEN z;
1732 6118 : b = FlxY_to_FlyX(b,sx);
1733 6116 : if ((ulong)dres >= p)
1734 2506 : z = FlxX_resultant(Fly_to_FlxY(a, sy), b, p, sx);
1735 : else
1736 : {
1737 3610 : ulong pi = SMALL_ULONG(p)? 0: get_Fl_red(p);
1738 3610 : z = Flx_FlxY_resultant_polint(a, b, p, pi, (ulong)dres, sy);
1739 : }
1740 6120 : return gc_upto(ltop,z);
1741 : }
1742 :
1743 : /* Return a t_POL in variable vc whose coeffs are the coeffs of b in
1744 : * variable v; vc must have higher priority than all variables occuring in b. */
1745 : GEN
1746 146016 : swap_vars(GEN b, long v, long vc)
1747 : {
1748 146016 : long i, n = RgX_degree(b, v);
1749 : GEN c, x;
1750 146014 : if (n < 0) return pol_0(vc);
1751 146014 : c = cgetg(n+3, t_POL); x = c + 2;
1752 146014 : c[1] = evalsigne(1) | evalvarn(vc);
1753 967844 : for (i = 0; i <= n; i++) gel(x,i) = polcoef_i(b, i, v);
1754 146010 : return c;
1755 : }
1756 :
1757 : /* assume varn(b) << varn(a) */
1758 : /* return a FpY*/
1759 : GEN
1760 15 : FpX_FpXY_resultant(GEN a, GEN b, GEN p)
1761 : {
1762 15 : long i,n,dres, db, vY = varn(b), vX = varn(a);
1763 : GEN la,x,y;
1764 :
1765 15 : if (lgefint(p) == 3)
1766 : {
1767 0 : ulong pp = uel(p,2);
1768 0 : b = ZXX_to_FlxX(b, pp, vX);
1769 0 : a = ZX_to_Flx(a, pp);
1770 0 : x = Flx_FlxY_resultant(a, b, pp);
1771 0 : return Flx_to_ZX(x);
1772 : }
1773 15 : db = RgXY_degreex(b);
1774 15 : dres = degpol(a)*degpol(b);
1775 15 : la = leading_coeff(a);
1776 15 : x = cgetg(dres+2, t_VEC);
1777 15 : y = cgetg(dres+2, t_VEC);
1778 : /* Evaluate at dres+ 1 points: 0 (if dres even) and +/- n, so that P_n(X) =
1779 : * P_{-n}(-X), where P_i is Lagrange polynomial: P_i(j) = delta_{i,j} */
1780 157 : for (i=0,n = 1; i < dres; n++)
1781 : {
1782 142 : gel(x,++i) = utoipos(n);
1783 142 : gel(y,i) = FpX_FpXY_eval_resultant(a,b,gel(x,i),p,la,db,vY);
1784 142 : gel(x,++i) = subiu(p,n);
1785 142 : gel(y,i) = FpX_FpXY_eval_resultant(a,b,gel(x,i),p,la,db,vY);
1786 : }
1787 15 : if (i == dres)
1788 : {
1789 0 : gel(x,++i) = gen_0;
1790 0 : gel(y,i) = FpX_FpXY_eval_resultant(a,b, gel(x,i), p,la,db,vY);
1791 : }
1792 15 : return FpV_polint(x,y, p, vY);
1793 : }
1794 :
1795 : GEN
1796 191 : FpX_composedsum(GEN P, GEN Q, GEN p)
1797 : {
1798 191 : pari_sp av = avma;
1799 191 : if (lgefint(p)==3)
1800 : {
1801 0 : ulong pp = p[2];
1802 0 : GEN z = Flx_composedsum(ZX_to_Flx(P, pp), ZX_to_Flx(Q, pp), pp);
1803 0 : return gc_upto(av, Flx_to_ZX(z));
1804 : }
1805 : else
1806 : {
1807 191 : long n = 1+ degpol(P)*degpol(Q);
1808 191 : GEN Pl = FpX_invLaplace(FpX_Newton(P,n,p), p);
1809 191 : GEN Ql = FpX_invLaplace(FpX_Newton(Q,n,p), p);
1810 191 : GEN L = FpX_Laplace(FpXn_mul(Pl, Ql, n, p), p);
1811 191 : GEN lead = Fp_mul(Fp_powu(leading_coeff(P),degpol(Q), p),
1812 191 : Fp_powu(leading_coeff(Q),degpol(P), p), p);
1813 191 : GEN R = FpX_fromNewton(L, p);
1814 191 : return gc_upto(av, FpX_Fp_mul(R, lead, p));
1815 : }
1816 : }
1817 :
1818 : GEN
1819 0 : FpX_composedprod(GEN P, GEN Q, GEN p)
1820 : {
1821 0 : pari_sp av = avma;
1822 0 : if (lgefint(p)==3)
1823 : {
1824 0 : ulong pp = p[2];
1825 0 : GEN z = Flx_composedprod(ZX_to_Flx(P, pp), ZX_to_Flx(Q, pp), pp);
1826 0 : return gc_upto(av, Flx_to_ZX(z));
1827 : }
1828 : else
1829 : {
1830 0 : long n = 1+ degpol(P)*degpol(Q);
1831 0 : GEN L = FpX_convol(FpX_Newton(P,n,p), FpX_Newton(Q,n,p), p);
1832 0 : return gc_upto(av,FpX_fromNewton(L, p));
1833 : }
1834 : }
1835 :
1836 : static GEN
1837 191 : _FpX_composedsum(void *E, GEN a, GEN b)
1838 191 : { return FpX_composedsum(a,b, (GEN)E); }
1839 :
1840 : GEN
1841 1637 : FpXV_composedsum(GEN V, GEN p)
1842 : {
1843 1637 : if (lgefint(p)==3)
1844 : {
1845 0 : ulong pp = p[2];
1846 0 : return Flx_to_ZX(FlxV_composedsum(ZXV_to_FlxV(V, pp), pp));
1847 : }
1848 1637 : return gen_product(V, (void *)p, &_FpX_composedsum);
1849 : }
1850 :
1851 : /* 0, 1, -1, 2, -2, ... */
1852 : #define next_lambda(a) (a>0 ? -a : 1-a)
1853 :
1854 : /* Assume A in Z[Y], B in Q[Y][X], B squarefree in (Q[Y]/(A))[X] and
1855 : * Res_Y(A, B) in Z[X]. Find a small lambda (start from *lambda, use
1856 : * next_lambda successively) such that C(X) = Res_Y(A(Y), B(X + lambda Y))
1857 : * is squarefree, reset *lambda to the chosen value and return C. Set LERS to
1858 : * the Last nonconstant polynomial in the Euclidean Remainder Sequence */
1859 : static GEN
1860 21847 : ZX_ZXY_resultant_LERS(GEN A, GEN B0, long *plambda, GEN *LERS)
1861 : {
1862 : ulong bound, dp;
1863 21847 : pari_sp av = avma, av2 = 0;
1864 21847 : long lambda = *plambda, degA = degpol(A), dres = degA*degpol(B0);
1865 : long stable, checksqfree, i,n, cnt, degB;
1866 21847 : long v, vX = varn(B0), vY = varn(A); /* vY < vX */
1867 : GEN x, y, dglist, B, q, a, b, ev, H, H0, H1, Hp, H0p, H1p, C0, C1;
1868 : forprime_t S;
1869 :
1870 21847 : if (degA == 1)
1871 : {
1872 1260 : GEN a1 = gel(A,3), a0 = gel(A,2);
1873 1260 : B = lambda? RgX_translate(B0, monomial(stoi(lambda), 1, vY)): B0;
1874 1260 : H = gsubst(B, vY, gdiv(gneg(a0),a1));
1875 1260 : if (!equali1(a1)) H = RgX_Rg_mul(H, powiu(a1, poldegree(B,vY)));
1876 1260 : *LERS = mkvec2(scalarpol_shallow(a0,vX), scalarpol_shallow(a1,vX));
1877 1260 : return gc_all(av, 2, &H, LERS);
1878 : }
1879 :
1880 20587 : dglist = Hp = H0p = H1p = C0 = C1 = NULL; /* gcc -Wall */
1881 20587 : C0 = cgetg(dres+2, t_VECSMALL);
1882 20587 : C1 = cgetg(dres+2, t_VECSMALL);
1883 20587 : dglist = cgetg(dres+1, t_VECSMALL);
1884 20587 : x = cgetg(dres+2, t_VECSMALL);
1885 20587 : y = cgetg(dres+2, t_VECSMALL);
1886 20587 : B0 = leafcopy(B0);
1887 20587 : A = leafcopy(A);
1888 20587 : B = B0;
1889 20587 : v = fetch_var_higher(); setvarn(A,v);
1890 : /* make sure p large enough */
1891 21473 : INIT:
1892 : /* always except the first time */
1893 21473 : if (av2) { set_avma(av2); lambda = next_lambda(lambda); }
1894 21473 : if (lambda) B = RgX_translate(B0, monomial(stoi(lambda), 1, vY));
1895 21473 : B = swap_vars(B, vY, v);
1896 : /* B0(lambda v + x, v) */
1897 21473 : if (DEBUGLEVEL>4) err_printf("Trying lambda = %ld\n", lambda);
1898 21473 : av2 = avma;
1899 :
1900 21473 : if (degA <= 3)
1901 : { /* sub-resultant faster for small degrees */
1902 10752 : H = RgX_resultant_all(A,B,&q);
1903 10752 : if (typ(q) != t_POL || degpol(q)!=1) goto INIT;
1904 10031 : H0 = gel(q,2);
1905 10031 : if (typ(H0) == t_POL) setvarn(H0,vX); else H0 = scalarpol(H0,vX);
1906 10031 : H1 = gel(q,3);
1907 10031 : if (typ(H1) == t_POL) setvarn(H1,vX); else H1 = scalarpol(H1,vX);
1908 10031 : if (!ZX_is_squarefree(H)) goto INIT;
1909 9989 : goto END;
1910 : }
1911 :
1912 10721 : H = H0 = H1 = NULL;
1913 10721 : degB = degpol(B);
1914 10721 : bound = ZX_ZXY_ResBound(A, B, NULL);
1915 10720 : if (DEBUGLEVEL>4) err_printf("bound for resultant coeffs: 2^%ld\n",bound);
1916 10720 : dp = 1;
1917 10720 : init_modular_big(&S);
1918 10720 : for(cnt = 0, checksqfree = 1;;)
1919 49204 : {
1920 59924 : ulong p = u_forprime_next(&S);
1921 : GEN Hi;
1922 59924 : a = ZX_to_Flx(A, p);
1923 59926 : b = ZXX_to_FlxX(B, p, varn(A));
1924 59925 : if (degpol(a) < degA || degpol(b) < degB) continue; /* p | lc(A)lc(B) */
1925 59925 : if (checksqfree)
1926 : { /* find degree list for generic Euclidean Remainder Sequence */
1927 10721 : long goal = minss(degpol(a), degpol(b)); /* longest possible */
1928 73203 : for (n=1; n <= goal; n++) dglist[n] = 0;
1929 10721 : setlg(dglist, 1);
1930 23678 : for (n=0; n <= dres; n++)
1931 : {
1932 23286 : ev = FlxY_evalx_drop(b, n, p);
1933 23286 : Flx_resultant_set_dglist(a, ev, dglist, p);
1934 23285 : if (lg(dglist)-1 == goal) break;
1935 : }
1936 : /* last pol in ERS has degree > 1 ? */
1937 10720 : goal = lg(dglist)-1;
1938 10720 : if (degpol(B) == 1) { if (!goal) goto INIT; }
1939 : else
1940 : {
1941 10664 : if (goal <= 1) goto INIT;
1942 10601 : if (dglist[goal] != 0 || dglist[goal-1] != 1) goto INIT;
1943 : }
1944 10657 : if (DEBUGLEVEL>4)
1945 0 : err_printf("Degree list for ERS (trials: %ld) = %Ps\n",n+1,dglist);
1946 : }
1947 :
1948 2145116 : for (i=0,n = 0; i <= dres; n++)
1949 : {
1950 2085249 : ev = FlxY_evalx_drop(b, n, p);
1951 2085109 : x[++i] = n; y[i] = Flx_resultant_all(a, ev, C0+i, C1+i, dglist, p);
1952 2085254 : if (!C1[i]) i--; /* C1(i) = 0. No way to recover C0(i) */
1953 : }
1954 59867 : Hi = Flv_Flm_polint(x, mkvec3(y,C0,C1), p, 0);
1955 59863 : Hp = gel(Hi,1); H0p = gel(Hi,2); H1p = gel(Hi,3);
1956 59863 : if (!H && degpol(Hp) != dres) continue;
1957 59863 : if (dp != 1) Hp = Flx_Fl_mul(Hp, Fl_powu(Fl_inv(dp,p), degA, p), p);
1958 59863 : if (checksqfree) {
1959 10658 : if (!Flx_is_squarefree(Hp, p)) goto INIT;
1960 10598 : if (DEBUGLEVEL>4) err_printf("Final lambda = %ld\n", lambda);
1961 10598 : checksqfree = 0;
1962 : }
1963 :
1964 59803 : if (!H)
1965 : { /* initialize */
1966 10598 : q = utoipos(p); stable = 0;
1967 10598 : H = ZX_init_CRT(Hp, p,vX);
1968 10598 : H0= ZX_init_CRT(H0p, p,vX);
1969 10598 : H1= ZX_init_CRT(H1p, p,vX);
1970 : }
1971 : else
1972 : {
1973 49205 : GEN qp = muliu(q,p);
1974 49203 : stable = ZX_incremental_CRT_raw(&H, Hp, q,qp, p)
1975 49205 : & ZX_incremental_CRT_raw(&H0,H0p, q,qp, p)
1976 49205 : & ZX_incremental_CRT_raw(&H1,H1p, q,qp, p);
1977 49205 : q = qp;
1978 : }
1979 : /* could make it probabilistic for H ? [e.g if stable twice, etc]
1980 : * Probabilistic anyway for H0, H1 */
1981 59803 : if (DEBUGLEVEL>5 && (stable || ++cnt==100))
1982 0 : { cnt=0; err_printf("%ld%%%s ",100*expi(q)/bound,stable?"s":""); }
1983 59802 : if (stable && (ulong)expi(q) >= bound) break; /* DONE */
1984 49204 : if (gc_needed(av,2))
1985 : {
1986 0 : if (DEBUGMEM>1) pari_warn(warnmem,"ZX_ZXY_rnfequation");
1987 0 : (void)gc_all(av2, 4, &H, &q, &H0, &H1);
1988 : }
1989 : }
1990 20587 : END:
1991 20587 : if (DEBUGLEVEL>5) err_printf(" done\n");
1992 20587 : setvarn(H, vX); (void)delete_var();
1993 20587 : *LERS = mkvec2(H0,H1);
1994 20587 : *plambda = lambda; return gc_all(av, 2, &H, LERS);
1995 : }
1996 :
1997 : GEN
1998 59639 : ZX_ZXY_resultant_all(GEN A, GEN B, long *plambda, GEN *LERS)
1999 : {
2000 59639 : if (LERS)
2001 : {
2002 21847 : if (!plambda)
2003 0 : pari_err_BUG("ZX_ZXY_resultant_all [LERS != NULL needs lambda]");
2004 21847 : return ZX_ZXY_resultant_LERS(A, B, plambda, LERS);
2005 : }
2006 37792 : return ZX_ZXY_rnfequation(A, B, plambda);
2007 : }
2008 :
2009 : /* If lambda = NULL, return caract(Mod(A, T)), T,A in Z[X].
2010 : * Otherwise find a small lambda such that caract (Mod(A + lambda X, T)) is
2011 : * squarefree */
2012 : GEN
2013 22594 : ZXQ_charpoly_sqf(GEN A, GEN T, long *lambda, long v)
2014 : {
2015 22594 : pari_sp av = avma;
2016 : GEN R, a;
2017 : long dA;
2018 : int delvar;
2019 :
2020 22594 : if (v < 0) v = 0;
2021 22594 : switch (typ(A))
2022 : {
2023 22594 : case t_POL: dA = degpol(A); if (dA > 0) break;
2024 0 : A = constant_coeff(A);
2025 0 : default:
2026 0 : if (lambda) { A = scalar_ZX_shallow(A,varn(T)); dA = 0; break;}
2027 0 : return gc_upto(av, gpowgs(gsub(pol_x(v), A), degpol(T)));
2028 : }
2029 22594 : delvar = 0;
2030 22594 : if (varncmp(varn(T), 0) <= 0)
2031 : {
2032 3681 : long v0 = fetch_var(); delvar = 1;
2033 3681 : T = leafcopy(T); setvarn(T,v0);
2034 3681 : A = leafcopy(A); setvarn(A,v0);
2035 : }
2036 22594 : R = ZX_ZXY_rnfequation(T, deg1pol_shallow(gen_1, gneg_i(A), 0), lambda);
2037 22594 : if (delvar) (void)delete_var();
2038 22594 : setvarn(R, v); a = leading_coeff(T);
2039 22594 : if (!gequal1(a)) R = gdiv(R, powiu(a, dA));
2040 22594 : return gc_upto(av, R);
2041 : }
2042 :
2043 : /* charpoly(Mod(A,T)), A may be in Q[X], but assume T and result are integral */
2044 : GEN
2045 995251 : ZXQ_charpoly(GEN A, GEN T, long v)
2046 : {
2047 995251 : return (degpol(T) < 16) ? RgXQ_charpoly_i(A,T,v): ZXQ_charpoly_sqf(A,T, NULL, v);
2048 : }
2049 :
2050 : GEN
2051 9772 : QXQ_charpoly(GEN A, GEN T, long v)
2052 : {
2053 9772 : pari_sp av = avma;
2054 9772 : GEN den, B = Q_remove_denom(A, &den);
2055 9772 : GEN P = ZXQ_charpoly(B, T, v);
2056 9772 : return gc_GEN(av, den ? RgX_rescale(P, ginv(den)): P);
2057 : }
2058 :
2059 : static ulong
2060 3863839 : ZX_resultant_prime(GEN a, GEN b, GEN dB, long degA, long degB, ulong p)
2061 : {
2062 3863839 : long dropa = degA - degpol(a), dropb = degB - degpol(b);
2063 : ulong H, dp;
2064 3863721 : if (dropa && dropb) return 0; /* p | lc(A), p | lc(B) */
2065 3863721 : H = Flx_resultant(a, b, p);
2066 3863503 : if (dropa)
2067 : { /* multiply by ((-1)^deg B lc(B))^(deg A - deg a) */
2068 0 : ulong c = b[degB+2]; /* lc(B) */
2069 0 : if (odd(degB)) c = p - c;
2070 0 : c = Fl_powu(c, dropa, p);
2071 0 : if (c != 1) H = Fl_mul(H, c, p);
2072 : }
2073 3863503 : else if (dropb)
2074 : { /* multiply by lc(A)^(deg B - deg b) */
2075 0 : ulong c = a[degA+2]; /* lc(A) */
2076 0 : c = Fl_powu(c, dropb, p);
2077 0 : if (c != 1) H = Fl_mul(H, c, p);
2078 : }
2079 3863503 : dp = dB ? umodiu(dB, p): 1;
2080 3863502 : if (dp != 1) H = Fl_mul(H, Fl_powu(Fl_inv(dp,p), degA, p), p);
2081 3863504 : return H;
2082 : }
2083 :
2084 : /* If B=NULL, assume B=A' */
2085 : static GEN
2086 1494401 : ZX_resultant_slice(GEN A, GEN B, GEN dB, GEN P, GEN *mod)
2087 : {
2088 1494401 : pari_sp av = avma, av2;
2089 1494401 : long degA, degB, i, n = lg(P)-1;
2090 : GEN H, T;
2091 :
2092 1494401 : degA = degpol(A);
2093 1494393 : degB = B? degpol(B): degA - 1;
2094 1494395 : if (n == 1)
2095 : {
2096 810594 : ulong Hp, p = uel(P,1);
2097 810594 : GEN a = ZX_to_Flx(A, p), b = B? ZX_to_Flx(B, p): Flx_deriv(a, p);
2098 810580 : Hp = ZX_resultant_prime(a, b, dB, degA, degB, p);
2099 810579 : set_avma(av); *mod = utoipos(p); return utoi(Hp);
2100 : }
2101 683801 : T = ZV_producttree(P);
2102 683802 : A = ZX_nv_mod_tree(A, P, T);
2103 683799 : if (B) B = ZX_nv_mod_tree(B, P, T);
2104 683799 : H = cgetg(n+1, t_VECSMALL); av2 = avma;
2105 3736767 : for(i=1; i <= n; i++, set_avma(av2))
2106 : {
2107 3052969 : ulong p = P[i];
2108 3052969 : GEN a = gel(A,i), b = B? gel(B,i): Flx_deriv(a, p);
2109 3053263 : H[i] = ZX_resultant_prime(a, b, dB, degA, degB, p);
2110 : }
2111 683798 : H = ZV_chinese_tree(H, P, T, ZV_chinesetree(P,T));
2112 683798 : *mod = gmael(T, lg(T)-1, 1); return gc_all(av, 2, &H, mod);
2113 : }
2114 :
2115 : GEN
2116 1494401 : ZX_resultant_worker(GEN P, GEN A, GEN B, GEN dB)
2117 : {
2118 1494401 : GEN V = cgetg(3, t_VEC);
2119 1494402 : if (typ(B) == t_INT) B = NULL;
2120 1494402 : if (!signe(dB)) dB = NULL;
2121 1494402 : gel(V,1) = ZX_resultant_slice(A, B, dB, P, &gel(V,2));
2122 1494377 : return V;
2123 : }
2124 :
2125 : /* Compute Res(A, B/dB) in Z, assuming A,B in Z[X], dB in Z or NULL (= 1)
2126 : * If B=NULL, take B = A' and assume deg A > 1 and 'bound' is set */
2127 : GEN
2128 1351019 : ZX_resultant_all(GEN A, GEN B, GEN dB, ulong bound)
2129 : {
2130 1351019 : pari_sp av = avma;
2131 : forprime_t S;
2132 : GEN H, worker;
2133 1351019 : if (!B && degpol(A)==2)
2134 : {
2135 114020 : GEN a = gel(A,4), b = gel(A,3), c = gel(A,2);
2136 114020 : H = mulii(a, subii(shifti(mulii(a, c), 2), sqri(b)));
2137 114012 : if (dB) H = diviiexact(H, sqri(dB));
2138 114012 : return gc_INT(av, H);
2139 : }
2140 1236996 : if (B)
2141 : {
2142 155216 : long a = degpol(A), b = degpol(B);
2143 155216 : if (a < 0 || b < 0) return gen_0;
2144 155186 : if (!a) return powiu(gel(A,2), b);
2145 155186 : if (!b) return powiu(gel(B,2), a);
2146 153441 : if (minss(a, b) <= 1)
2147 : {
2148 76738 : H = RgX_resultant_all(A, B, NULL);
2149 76738 : if (dB) H = diviiexact(H, powiu(dB, a));
2150 76738 : return gc_INT(av, H);
2151 : }
2152 76703 : if (!bound) bound = ZX_ZXY_ResBound(A, B, dB);
2153 : }
2154 1158494 : worker = snm_closure(is_entry("_ZX_resultant_worker"),
2155 : mkvec3(A, B? B: gen_0, dB? dB: gen_0));
2156 1158634 : init_modular_big(&S);
2157 1158586 : H = gen_crt("ZX_resultant_all", worker, &S, dB, bound, 0, NULL,
2158 : ZV_chinese_center, Fp_center);
2159 1158643 : return gc_INT(av, H);
2160 : }
2161 :
2162 : /* A0 and B0 in Q[X] */
2163 : GEN
2164 56 : QX_resultant(GEN A0, GEN B0)
2165 : {
2166 : GEN s, a, b, A, B;
2167 56 : pari_sp av = avma;
2168 :
2169 56 : A = Q_primitive_part(A0, &a);
2170 56 : B = Q_primitive_part(B0, &b);
2171 56 : s = ZX_resultant(A, B);
2172 56 : if (!signe(s)) { set_avma(av); return gen_0; }
2173 56 : if (a) s = gmul(s, gpowgs(a,degpol(B)));
2174 56 : if (b) s = gmul(s, gpowgs(b,degpol(A)));
2175 56 : return gc_upto(av, s);
2176 : }
2177 :
2178 : GEN
2179 57309 : ZX_resultant(GEN A, GEN B) { return ZX_resultant_all(A,B,NULL,0); }
2180 :
2181 : GEN
2182 0 : QXQ_intnorm(GEN A, GEN B)
2183 : {
2184 : GEN c, n, R, lB;
2185 0 : long dA = degpol(A), dB = degpol(B);
2186 0 : pari_sp av = avma;
2187 0 : if (dA < 0) return gen_0;
2188 0 : A = Q_primitive_part(A, &c);
2189 0 : if (!c || typ(c) == t_INT) {
2190 0 : n = c;
2191 0 : R = ZX_resultant(B, A);
2192 : } else {
2193 0 : n = gel(c,1);
2194 0 : R = ZX_resultant_all(B, A, gel(c,2), 0);
2195 : }
2196 0 : if (n && !equali1(n)) R = mulii(R, powiu(n, dB));
2197 0 : lB = leading_coeff(B);
2198 0 : if (!equali1(lB)) R = diviiexact(R, powiu(lB, dA));
2199 0 : return gc_INT(av, R);
2200 : }
2201 :
2202 : GEN
2203 19418 : QXQ_norm(GEN A, GEN B)
2204 : {
2205 : GEN c, R, lB;
2206 19418 : long dA = degpol(A), dB = degpol(B);
2207 19418 : pari_sp av = avma;
2208 19418 : if (dA < 0) return gen_0;
2209 19418 : A = Q_primitive_part(A, &c);
2210 19418 : R = ZX_resultant(B, A);
2211 19418 : if (c) R = gmul(R, gpowgs(c, dB));
2212 19418 : lB = leading_coeff(B);
2213 19418 : if (!equali1(lB)) R = gdiv(R, gpowgs(lB, dA));
2214 19418 : return gc_upto(av, R);
2215 : }
2216 :
2217 : /* assume x has integral coefficients */
2218 : GEN
2219 1199175 : ZX_disc_all(GEN x, ulong bound)
2220 : {
2221 1199175 : pari_sp av = avma;
2222 1199175 : long s, d = degpol(x);
2223 : GEN l, R;
2224 :
2225 1199172 : if (d <= 1) return d == 1? gen_1: gen_0;
2226 1195879 : s = (d & 2) ? -1: 1;
2227 1195879 : l = leading_coeff(x);
2228 1195871 : if (!bound) bound = ZX_discbound(x);
2229 1195772 : R = ZX_resultant_all(x, NULL, NULL, bound);
2230 1195919 : if (is_pm1(l))
2231 1016897 : { if (signe(l) < 0) s = -s; }
2232 : else
2233 179023 : R = diviiexact(R,l);
2234 1195920 : if (s == -1) togglesign_safe(&R);
2235 1195916 : return gc_INT(av,R);
2236 : }
2237 :
2238 : GEN
2239 1137331 : ZX_disc(GEN x) { return ZX_disc_all(x,0); }
2240 :
2241 : static GEN
2242 10784 : ZXQX_resultant_prime(GEN a, GEN b, GEN dB, long degA, long degB, GEN T, ulong p)
2243 : {
2244 10784 : long dropa = degA - degpol(a), dropb = degB - degpol(b);
2245 : GEN H, dp;
2246 10785 : if (dropa && dropb) return pol0_Flx(T[1]); /* p | lc(A), p | lc(B) */
2247 10785 : H = FlxqX_saferesultant(a, b, T, p);
2248 10784 : if (!H) return NULL;
2249 10784 : if (dropa)
2250 : { /* multiply by ((-1)^deg B lc(B))^(deg A - deg a) */
2251 0 : GEN c = gel(b,degB+2); /* lc(B) */
2252 0 : if (odd(degB)) c = Flx_neg(c, p);
2253 0 : c = Flxq_powu(c, dropa, T, p);
2254 0 : if (!Flx_equal1(c)) H = Flxq_mul(H, c, T, p);
2255 : }
2256 10784 : else if (dropb)
2257 : { /* multiply by lc(A)^(deg B - deg b) */
2258 0 : GEN c = gel(a,degA+2); /* lc(A) */
2259 0 : c = Flxq_powu(c, dropb, T, p);
2260 0 : if (!Flx_equal1(c)) H = Flxq_mul(H, c, T, p);
2261 : }
2262 10784 : dp = dB ? ZX_to_Flx(dB, p): pol1_Flx(T[1]);
2263 10785 : if (!Flx_equal1(dp))
2264 : {
2265 0 : GEN idp = Flxq_invsafe(dp, T, p);
2266 0 : if (!idp) return NULL;
2267 0 : H = Flxq_mul(H, Flxq_powu(idp, degA, T, p), T, p);
2268 : }
2269 10785 : return H;
2270 : }
2271 :
2272 : /* If B=NULL, assume B=A' */
2273 : static GEN
2274 4687 : ZXQX_resultant_slice(GEN A, GEN B, GEN U, GEN dB, GEN P, GEN *mod)
2275 : {
2276 4687 : pari_sp av = avma;
2277 4687 : long degA, degB, i, n = lg(P)-1;
2278 : GEN H, T;
2279 4687 : long v = varn(U), redo = 0;
2280 :
2281 4687 : degA = degpol(A);
2282 4687 : degB = B? degpol(B): degA - 1;
2283 4687 : if (n == 1)
2284 : {
2285 2953 : ulong p = uel(P,1);
2286 2953 : GEN a = ZXX_to_FlxX(A, p, v), b = B? ZXX_to_FlxX(B, p, v): FlxX_deriv(a, p);
2287 2953 : GEN u = ZX_to_Flx(U, p);
2288 2953 : GEN Hp = ZXQX_resultant_prime(a, b, dB, degA, degB, u, p);
2289 2953 : if (!Hp) { set_avma(av); *mod = gen_1; return pol_0(v); }
2290 2953 : Hp = gc_upto(av, Flx_to_ZX(Hp)); *mod = utoipos(p); return Hp;
2291 : }
2292 1734 : T = ZV_producttree(P);
2293 1734 : A = ZXX_nv_mod_tree(A, P, T, v);
2294 1734 : if (B) B = ZXX_nv_mod_tree(B, P, T, v);
2295 1734 : U = ZX_nv_mod_tree(U, P, T);
2296 1734 : H = cgetg(n+1, t_VEC);
2297 9567 : for(i=1; i <= n; i++)
2298 : {
2299 7833 : ulong p = P[i];
2300 7833 : GEN a = gel(A,i), b = B? gel(B,i): FlxX_deriv(a, p), u = gel(U, i);
2301 7831 : GEN h = ZXQX_resultant_prime(a, b, dB, degA, degB, u, p);
2302 7832 : if (!h)
2303 : {
2304 0 : gel(H,i) = pol_0(v);
2305 0 : P[i] = 1; redo = 1;
2306 : }
2307 : else
2308 7832 : gel(H,i) = h;
2309 : }
2310 1734 : if (redo) T = ZV_producttree(P);
2311 1734 : H = nxV_chinese_center_tree(H, P, T, ZV_chinesetree(P, T));
2312 1734 : *mod = gmael(T, lg(T)-1, 1); return gc_all(av, 2, &H, mod);
2313 : }
2314 :
2315 : GEN
2316 4687 : ZXQX_resultant_worker(GEN P, GEN A, GEN B, GEN T, GEN dB)
2317 : {
2318 4687 : GEN V = cgetg(3, t_VEC);
2319 4687 : if (isintzero(B)) B = NULL;
2320 4687 : if (!signe(dB)) dB = NULL;
2321 4687 : gel(V,1) = ZXQX_resultant_slice(A, B, T, dB, P, &gel(V,2));
2322 4687 : return V;
2323 : }
2324 :
2325 : static ulong
2326 4091 : ZXQX_resultant_bound_i(GEN nf, GEN A, GEN B, GEN (*f)(GEN,GEN,long))
2327 : {
2328 4091 : pari_sp av = avma;
2329 4091 : GEN r, M = nf_L2_bound(nf, NULL, &r);
2330 4091 : long v = nf_get_varn(nf), i, l = lg(r);
2331 4091 : GEN a = cgetg(l, t_COL);
2332 12662 : for (i = 1; i < l; i++)
2333 8571 : gel(a, i) = f(gsubst(A, v, gel(r,i)), gsubst(B, v, gel(r,i)), DEFAULTPREC);
2334 4091 : return gc_ulong(av, (ulong) dbllog2(gmul(M,RgC_fpnorml2(a, DEFAULTPREC))));
2335 : }
2336 : static ulong
2337 3776 : ZXQX_resultant_bound(GEN nf, GEN A, GEN B)
2338 3776 : { return ZXQX_resultant_bound_i(nf, A, B, &RgX_RgXY_ResBound); }
2339 :
2340 : static GEN
2341 56 : _ZXQ_powu(GEN x, ulong u, GEN T)
2342 56 : { return typ(x) == t_INT? powiu(x, u): ZXQ_powu(x, u, T); }
2343 :
2344 : /* Compute Res(A, B/dB) in Z[X]/T, assuming A,B in Z[X,Y], dB in Z or NULL (= 1)
2345 : * If B=NULL, take B = A' and assume deg A > 1 */
2346 : static GEN
2347 3773 : ZXQX_resultant_all(GEN A, GEN B, GEN T, GEN dB, ulong bound)
2348 : {
2349 3773 : pari_sp av = avma;
2350 : forprime_t S;
2351 : GEN H, worker;
2352 3773 : if (B)
2353 : {
2354 63 : long a = degpol(A), b = degpol(B);
2355 63 : if (a < 0 || b < 0) return gen_0;
2356 63 : if (!a) return _ZXQ_powu(gel(A,2), b, T);
2357 63 : if (!b) return _ZXQ_powu(gel(B,2), a, T);
2358 : } else
2359 3710 : if (!bound) B = RgX_deriv(A);
2360 3773 : if (!bound) bound = ZXQX_resultant_bound(nfinit(T, DEFAULTPREC), A, B);
2361 3773 : worker = snm_closure(is_entry("_ZXQX_resultant_worker"),
2362 : mkvec4(A, B? B: gen_0, T, dB? dB: gen_0));
2363 3773 : init_modular_big(&S);
2364 3773 : H = gen_crt("ZXQX_resultant_all", worker, &S, dB, bound, 0, NULL,
2365 : nxV_chinese_center, FpX_center);
2366 3773 : if (DEBUGLEVEL)
2367 0 : err_printf("ZXQX_resultant_all: a priori bound: %lu, a posteriori: %lu\n",
2368 : bound, expi(gsupnorm(H, DEFAULTPREC)));
2369 3773 : return gc_upto(av, H);
2370 : }
2371 :
2372 : GEN
2373 119 : nfX_resultant(GEN nf, GEN x, GEN y)
2374 : {
2375 119 : pari_sp av = avma;
2376 119 : GEN cx, cy, D, T = nf_get_pol(nf);
2377 119 : long dx = degpol(x), dy = degpol(y);
2378 119 : if (dx < 0 || dy < 0) return gen_0;
2379 119 : x = Q_primitive_part(x, &cx); if (cx) cx = gpowgs(cx, dy);
2380 119 : y = Q_primitive_part(y, &cy); if (cy) cy = gpowgs(cy, dx);
2381 119 : if (!dx) D = _ZXQ_powu(gel(x,2), dy, T);
2382 119 : else if (!dy) D = _ZXQ_powu(gel(y,2), dx, T);
2383 : else
2384 : {
2385 63 : ulong bound = ZXQX_resultant_bound(nf, x, y);
2386 63 : D = ZXQX_resultant_all(x, y, T, NULL, bound);
2387 : }
2388 119 : cx = mul_content(cx, cy); if (cx) D = gmul(D, cx);
2389 119 : return gc_upto(av, D);
2390 : }
2391 :
2392 : static GEN
2393 252 : to_ZX(GEN a, long v) { return typ(a)==t_INT? scalarpol(a,v): a; }
2394 :
2395 : static GEN
2396 3710 : ZXQX_disc_all(GEN x, GEN T, ulong bound)
2397 : {
2398 3710 : pari_sp av = avma;
2399 3710 : long s, d = degpol(x), v = varn(T);
2400 : GEN l, R;
2401 :
2402 3710 : if (d <= 1) return d == 1? pol_1(v): pol_0(v);
2403 3710 : s = (d & 2) ? -1: 1;
2404 3710 : l = leading_coeff(x);
2405 3710 : R = ZXQX_resultant_all(x, NULL, T, NULL, bound);
2406 3710 : if (!gequal1(l)) R = QXQ_div(R, to_ZX(l,v), T);
2407 3710 : if (s == -1) R = RgX_neg(R);
2408 3710 : return gc_upto(av, R);
2409 : }
2410 :
2411 : GEN
2412 7 : QX_disc(GEN x)
2413 : {
2414 7 : pari_sp av = avma;
2415 7 : GEN c, d = ZX_disc( Q_primitive_part(x, &c) );
2416 7 : if (c) d = gmul(d, gpowgs(c, 2*degpol(x) - 2));
2417 7 : return gc_upto(av, d);
2418 : }
2419 :
2420 : GEN
2421 3941 : nfX_disc(GEN nf, GEN x)
2422 : {
2423 3941 : pari_sp av = avma;
2424 3941 : GEN c, D, T = nf_get_pol(nf);
2425 : ulong bound;
2426 3941 : long d = degpol(x), v = varn(T);
2427 3941 : if (d <= 1) return d == 1? pol_1(v): pol_0(v);
2428 3710 : x = Q_primitive_part(x, &c);
2429 3710 : bound = ZXQX_resultant_bound(nf, x, RgX_deriv(x));
2430 3710 : D = ZXQX_disc_all(x, T, bound);
2431 3710 : if (c) D = gmul(D, gpowgs(c, 2*d - 2));
2432 3710 : return gc_upto(av, D);
2433 : }
2434 :
2435 : GEN
2436 837295 : QXQ_mul(GEN x, GEN y, GEN T)
2437 : {
2438 837295 : GEN dx, nx = Q_primitive_part(x, &dx);
2439 837293 : GEN dy, ny = Q_primitive_part(y, &dy);
2440 837288 : GEN z = ZXQ_mul(nx, ny, T);
2441 837292 : if (dx || dy)
2442 : {
2443 834492 : GEN d = dx ? dy ? gmul(dx, dy): dx : dy;
2444 834492 : if (!gequal1(d)) z = ZX_Q_mul(z, d);
2445 : }
2446 837295 : return z;
2447 : }
2448 :
2449 : GEN
2450 399481 : QXQ_sqr(GEN x, GEN T)
2451 : {
2452 399481 : GEN dx, nx = Q_primitive_part(x, &dx);
2453 399481 : GEN z = ZXQ_sqr(nx, T);
2454 399481 : if (dx)
2455 397745 : z = ZX_Q_mul(z, gsqr(dx));
2456 399481 : return z;
2457 : }
2458 :
2459 : static GEN
2460 212647 : QXQ_inv_slice(GEN A, GEN B, GEN P, GEN *mod)
2461 : {
2462 212647 : pari_sp av = avma;
2463 212647 : long i, n = lg(P)-1, v = varn(A), redo = 0;
2464 : GEN H, T;
2465 212647 : if (n == 1)
2466 : {
2467 165610 : ulong p = uel(P,1);
2468 165610 : GEN a = ZX_to_Flx(A, p), b = ZX_to_Flx(B, p);
2469 165610 : GEN U = Flxq_invsafe(a, b, p);
2470 165610 : if (!U)
2471 : {
2472 24 : set_avma(av);
2473 24 : *mod = gen_1; return pol_0(v);
2474 : }
2475 165586 : H = gc_GEN(av, Flx_to_ZX(U));
2476 165586 : *mod = utoipos(p); return H;
2477 : }
2478 47037 : T = ZV_producttree(P);
2479 47037 : A = ZX_nv_mod_tree(A, P, T);
2480 47036 : B = ZX_nv_mod_tree(B, P, T);
2481 47037 : H = cgetg(n+1, t_VEC);
2482 237831 : for(i=1; i <= n; i++)
2483 : {
2484 190795 : ulong p = P[i];
2485 190795 : GEN a = gel(A,i), b = gel(B,i);
2486 190795 : GEN U = Flxq_invsafe(a, b, p);
2487 190793 : if (!U)
2488 : {
2489 601 : gel(H,i) = pol_0(v);
2490 601 : P[i] = 1; redo = 1;
2491 : }
2492 : else
2493 190192 : gel(H,i) = U;
2494 : }
2495 47036 : if (redo) T = ZV_producttree(P);
2496 47036 : H = nxV_chinese_center_tree(H, P, T, ZV_chinesetree(P, T));
2497 47037 : *mod = gmael(T, lg(T)-1, 1); return gc_all(av, 2, &H, mod);
2498 : }
2499 :
2500 : GEN
2501 212647 : QXQ_inv_worker(GEN P, GEN A, GEN B)
2502 : {
2503 212647 : GEN V = cgetg(3, t_VEC);
2504 212647 : gel(V,1) = QXQ_inv_slice(A, B, P, &gel(V,2));
2505 212647 : return V;
2506 : }
2507 :
2508 : /* lift(1 / Mod(A,B)). B a ZX, A a scalar or a QX */
2509 : GEN
2510 145976 : QXQ_inv(GEN A, GEN B)
2511 : {
2512 : GEN D, Ap, Bp;
2513 : ulong pp;
2514 145976 : pari_sp av2, av = avma;
2515 : forprime_t S;
2516 145976 : GEN worker, U, H = NULL, mod = gen_1;
2517 : pari_timer ti;
2518 : long k, dA, dB;
2519 145976 : if (is_scalar_t(typ(A))) return scalarpol(ginv(A), varn(B));
2520 : /* A a QX, B a ZX */
2521 145976 : A = Q_primitive_part(A, &D);
2522 145976 : dA = degpol(A); dB= degpol(B);
2523 : /* A, B in Z[X] */
2524 145976 : init_modular_small(&S);
2525 : do {
2526 145976 : pp = u_forprime_next(&S);
2527 145976 : Ap = ZX_to_Flx(A, pp);
2528 145976 : Bp = ZX_to_Flx(B, pp);
2529 145976 : } while (degpol(Ap) != dA || degpol(Bp) != dB);
2530 145976 : if (degpol(Flx_gcd(Ap, Bp, pp)) != 0 && degpol(ZX_gcd(A,B))!=0)
2531 14 : pari_err_INV("QXQ_inv",mkpolmod(A,B));
2532 145962 : worker = snm_closure(is_entry("_QXQ_inv_worker"), mkvec2(A, B));
2533 145962 : av2 = avma;
2534 145962 : for (k = 1; ;k *= 2)
2535 42486 : {
2536 : GEN res, b, N, den;
2537 188448 : gen_inccrt_i("QXQ_inv", worker, NULL, (k+1)>>1, 0, &S, &H, &mod,
2538 : nxV_chinese_center, FpX_center);
2539 188448 : (void)gc_all(av2, 2, &H, &mod);
2540 188448 : b = sqrti(shifti(mod,-1));
2541 188448 : if (DEBUGLEVEL>5) timer_start(&ti);
2542 188448 : U = FpX_ratlift(H, mod, b, b, NULL);
2543 188448 : if (DEBUGLEVEL>5) timer_printf(&ti,"QXQ_inv: ratlift");
2544 194161 : if (!U) continue;
2545 151675 : N = Q_remove_denom(U, &den); if (!den) den = gen_1;
2546 151675 : res = Flx_rem(Flx_Fl_sub(Flx_mul(Ap, ZX_to_Flx(N,pp), pp),
2547 : umodiu(den, pp), pp), Bp, pp);
2548 151675 : if (degpol(res) >= 0) continue;
2549 145962 : res = ZX_Z_sub(ZX_mul(A, N), den);
2550 145962 : res = ZX_is_monic(B) ? ZX_rem(res, B): RgX_pseudorem(res, B);
2551 145962 : if (DEBUGLEVEL>5) timer_printf(&ti,"QXQ_inv: final check");
2552 145962 : if (degpol(res)<0)
2553 : {
2554 145962 : if (D) U = RgX_Rg_div(U, D);
2555 145962 : return gc_GEN(av, U);
2556 : }
2557 : }
2558 : }
2559 :
2560 : static GEN
2561 120511 : QXQ_div_slice(GEN A, GEN B, GEN C, GEN P, GEN *mod)
2562 : {
2563 120511 : pari_sp av = avma;
2564 120511 : long i, n = lg(P)-1, v = varn(A), redo = 0;
2565 : GEN H, T;
2566 120511 : if (n == 1)
2567 : {
2568 44149 : ulong p = uel(P,1);
2569 44149 : GEN a = ZX_to_Flx(A, p), b = ZX_to_Flx(B, p), c = ZX_to_Flx(C, p);
2570 44149 : GEN bi = Flxq_invsafe(b, c, p), U;
2571 44149 : if (!bi)
2572 : {
2573 0 : set_avma(av);
2574 0 : *mod = gen_1; return pol_0(v);
2575 : }
2576 44149 : U = Flxq_mul(a, bi, c, p);
2577 44149 : H = gc_GEN(av, Flx_to_ZX(U));
2578 44149 : *mod = utoipos(p); return H;
2579 : }
2580 76362 : T = ZV_producttree(P);
2581 76362 : A = ZX_nv_mod_tree(A, P, T);
2582 76362 : B = ZX_nv_mod_tree(B, P, T);
2583 76362 : C = ZX_nv_mod_tree(C, P, T);
2584 76362 : H = cgetg(n+1, t_VEC);
2585 337311 : for(i=1; i <= n; i++)
2586 : {
2587 260950 : ulong p = P[i];
2588 260950 : GEN a = gel(A,i), b = gel(B,i), c = gel(C, i);
2589 260950 : GEN bi = Flxq_invsafe(b, c, p);
2590 260953 : if (!bi)
2591 : {
2592 4 : gel(H,i) = pol_0(v);
2593 4 : P[i] = 1; redo = 1;
2594 : }
2595 : else
2596 260949 : gel(H,i) = Flxq_mul(a, bi, c, p);
2597 : }
2598 76361 : if (redo) T = ZV_producttree(P);
2599 76361 : H = nxV_chinese_center_tree(H, P, T, ZV_chinesetree(P, T));
2600 76362 : *mod = gmael(T, lg(T)-1, 1); return gc_all(av, 2, &H, mod);
2601 : }
2602 :
2603 : GEN
2604 120511 : QXQ_div_worker(GEN P, GEN A, GEN B, GEN C)
2605 : {
2606 120511 : GEN V = cgetg(3, t_VEC);
2607 120511 : gel(V,1) = QXQ_div_slice(A, B, C, P, &gel(V,2));
2608 120511 : return V;
2609 : }
2610 :
2611 : /* lift(Mod(A/B, C)). C a ZX, A, B a scalar or a QX */
2612 : GEN
2613 32759 : QXQ_div(GEN A, GEN B, GEN C)
2614 : {
2615 : GEN DA, DB, Ap, Bp, Cp;
2616 : ulong pp;
2617 32759 : pari_sp av2, av = avma;
2618 : forprime_t S;
2619 32759 : GEN worker, U, H = NULL, mod = gen_1;
2620 : pari_timer ti;
2621 : long k, dA, dB, dC;
2622 32759 : if (is_scalar_t(typ(A))) return scalarpol(ginv(A), varn(B));
2623 : /* A a QX, B a ZX */
2624 32759 : A = Q_primitive_part(A, &DA);
2625 32759 : B = Q_primitive_part(B, &DB);
2626 32759 : dA = degpol(A); dB = degpol(B); dC = degpol(C);
2627 : /* A, B in Z[X] */
2628 32759 : init_modular_small(&S);
2629 : do {
2630 32759 : pp = u_forprime_next(&S);
2631 32759 : Ap = ZX_to_Flx(A, pp);
2632 32759 : Bp = ZX_to_Flx(B, pp);
2633 32759 : Cp = ZX_to_Flx(C, pp);
2634 32759 : } while (degpol(Ap) != dA || degpol(Bp) != dB || degpol(Cp) != dC);
2635 32759 : if (degpol(Flx_gcd(Bp, Cp, pp)) != 0 && degpol(ZX_gcd(B,C))!=0)
2636 0 : pari_err_INV("QXQ_div",mkpolmod(B,C));
2637 32758 : worker = snm_closure(is_entry("_QXQ_div_worker"), mkvec3(A, B, C));
2638 32759 : av2 = avma;
2639 32759 : for (k = 1; ;k *= 2)
2640 46576 : {
2641 : GEN res, b, N, den;
2642 79335 : gen_inccrt_i("QXQ_div", worker, NULL, (k+1)>>1, 0, &S, &H, &mod,
2643 : nxV_chinese_center, FpX_center);
2644 79335 : (void)gc_all(av2, 2, &H, &mod);
2645 79335 : b = sqrti(shifti(mod,-1));
2646 79335 : if (DEBUGLEVEL>5) timer_start(&ti);
2647 79335 : U = FpX_ratlift(H, mod, b, b, NULL);
2648 79335 : if (DEBUGLEVEL>5) timer_printf(&ti,"QXQ_div: ratlift");
2649 89945 : if (!U) continue;
2650 43369 : N = Q_remove_denom(U, &den); if (!den) den = gen_1;
2651 43369 : res = Flx_rem(Flx_sub(Flx_mul(Bp, ZX_to_Flx(N,pp), pp),
2652 : Flx_Fl_mul(Ap, umodiu(den, pp), pp), pp), Cp, pp);
2653 43369 : if (degpol(res) >= 0) continue;
2654 32759 : res = ZX_sub(ZX_mul(B, N), ZX_Z_mul(A,den));
2655 32759 : res = ZX_is_monic(C) ? ZX_rem(res, C): RgX_pseudorem(res, C);
2656 32759 : if (DEBUGLEVEL>5) timer_printf(&ti,"QXQ_div: final check");
2657 32759 : if (degpol(res)<0)
2658 : {
2659 32759 : if (DA && DB) U = RgX_Rg_mul(U, gdiv(DA,DB));
2660 27831 : else if (DA) U = RgX_Rg_mul(U, DA);
2661 15743 : else if (DB) U = RgX_Rg_div(U, DB);
2662 32759 : return gc_GEN(av, U);
2663 : }
2664 : }
2665 : }
2666 :
2667 : /************************************************************************
2668 : * *
2669 : * ZXQ_minpoly *
2670 : * *
2671 : ************************************************************************/
2672 :
2673 : static GEN
2674 3523 : ZXQ_minpoly_slice(GEN A, GEN B, long d, GEN P, GEN *mod)
2675 : {
2676 3523 : pari_sp av = avma;
2677 3523 : long i, n = lg(P)-1, v = evalvarn(varn(B));
2678 : GEN H, T;
2679 3523 : if (n == 1)
2680 : {
2681 716 : ulong p = uel(P,1);
2682 716 : GEN a = ZX_to_Flx(A, p), b = ZX_to_Flx(B, p);
2683 716 : GEN Hp = Flxq_minpoly(a, b, p);
2684 716 : if (degpol(Hp) != d) { p = 1; Hp = pol0_Flx(v); }
2685 716 : H = gc_upto(av, Flx_to_ZX(Hp));
2686 716 : *mod = utoipos(p); return H;
2687 : }
2688 2807 : T = ZV_producttree(P);
2689 2807 : A = ZX_nv_mod_tree(A, P, T);
2690 2807 : B = ZX_nv_mod_tree(B, P, T);
2691 2807 : H = cgetg(n+1, t_VEC);
2692 16838 : for(i=1; i <= n; i++)
2693 : {
2694 14031 : ulong p = P[i];
2695 14031 : GEN a = gel(A,i), b = gel(B,i);
2696 14031 : GEN m = Flxq_minpoly(a, b, p);
2697 14031 : if (degpol(m) != d) { P[i] = 1; m = pol0_Flx(v); }
2698 14031 : gel(H, i) = m;
2699 : }
2700 2807 : H = nxV_chinese_center_tree(H, P, T, ZV_chinesetree(P, T));
2701 2807 : *mod = gmael(T, lg(T)-1, 1); return gc_all(av, 2, &H, mod);
2702 : }
2703 :
2704 : GEN
2705 3523 : ZXQ_minpoly_worker(GEN P, GEN A, GEN B, long d)
2706 : {
2707 3523 : GEN V = cgetg(3, t_VEC);
2708 3523 : gel(V,1) = ZXQ_minpoly_slice(A, B, d, P, &gel(V,2));
2709 3523 : return V;
2710 : }
2711 :
2712 : GEN
2713 1701 : ZXQ_minpoly(GEN A, GEN B, long d, ulong bound)
2714 : {
2715 1701 : pari_sp av = avma;
2716 : GEN worker, H, dB;
2717 : forprime_t S;
2718 1701 : B = Q_remove_denom(B, &dB);
2719 1701 : worker = strtoclosure("_ZXQ_minpoly_worker", 3, A, B, stoi(d));
2720 1701 : init_modular_big(&S);
2721 1701 : H = gen_crt("ZXQ_minpoly", worker, &S, dB, bound, 0, NULL,
2722 : nxV_chinese_center, FpX_center_i);
2723 1701 : return gc_GEN(av, H);
2724 : }
2725 :
2726 : /************************************************************************
2727 : * *
2728 : * ZX_ZXY_resultant *
2729 : * *
2730 : ************************************************************************/
2731 :
2732 : static GEN
2733 364807 : ZX_ZXY_resultant_prime(GEN a, GEN b, ulong dp, ulong p,
2734 : long degA, long degB, long dres, long sX)
2735 : {
2736 364807 : long dropa = degA - degpol(a), dropb = degB - degpol(b);
2737 364805 : ulong pi = SMALL_ULONG(p)? 0: get_Fl_red(p);
2738 364803 : GEN Hp = Flx_FlxY_resultant_polint(a, b, p, pi, dres, sX);
2739 364811 : if (dropa && dropb)
2740 0 : Hp = zero_Flx(sX);
2741 : else {
2742 364811 : if (dropa)
2743 : { /* multiply by ((-1)^deg B lc(B))^(deg A - deg a) */
2744 0 : GEN c = gel(b,degB+2); /* lc(B) */
2745 0 : if (odd(degB)) c = Flx_neg(c, p);
2746 0 : if (!Flx_equal1(c)) {
2747 0 : c = Flx_powu_pre(c, dropa, p, pi);
2748 0 : if (!Flx_equal1(c)) Hp = Flx_mul_pre(Hp, c, p, pi);
2749 : }
2750 : }
2751 364811 : else if (dropb)
2752 : { /* multiply by lc(A)^(deg B - deg b) */
2753 0 : ulong c = uel(a, degA+2); /* lc(A) */
2754 0 : c = Fl_powu(c, dropb, p);
2755 0 : if (c != 1) Hp = Flx_Fl_mul_pre(Hp, c, p, pi);
2756 : }
2757 : }
2758 364811 : if (dp != 1) Hp = Flx_Fl_mul_pre(Hp, Fl_powu_pre(Fl_inv(dp,p), degA, p, pi), p, pi);
2759 364811 : return Hp;
2760 : }
2761 :
2762 : static GEN
2763 124940 : ZX_ZXY_resultant_slice(GEN A, GEN B, GEN dB, long degA, long degB, long dres,
2764 : GEN P, GEN *mod, long sX, long vY)
2765 : {
2766 124940 : pari_sp av = avma;
2767 124940 : long i, n = lg(P)-1;
2768 : GEN H, T, D;
2769 124940 : if (n == 1)
2770 : {
2771 40164 : ulong p = uel(P,1);
2772 40164 : ulong dp = dB ? umodiu(dB, p): 1;
2773 40164 : GEN a = ZX_to_Flx(A, p), b = ZXX_to_FlxX(B, p, vY);
2774 40165 : GEN Hp = ZX_ZXY_resultant_prime(a, b, dp, p, degA, degB, dres, sX);
2775 40164 : H = gc_upto(av, Flx_to_ZX(Hp));
2776 40164 : *mod = utoipos(p); return H;
2777 : }
2778 84776 : T = ZV_producttree(P);
2779 84776 : A = ZX_nv_mod_tree(A, P, T);
2780 84776 : B = ZXX_nv_mod_tree(B, P, T, vY);
2781 84776 : D = dB ? Z_ZV_mod_tree(dB, P, T): NULL;
2782 84776 : H = cgetg(n+1, t_VEC);
2783 364102 : for(i=1; i <= n; i++)
2784 : {
2785 279325 : ulong p = P[i];
2786 279325 : GEN a = gel(A,i), b = gel(B,i);
2787 279325 : ulong dp = D ? uel(D, i): 1;
2788 279325 : gel(H,i) = ZX_ZXY_resultant_prime(a, b, dp, p, degA, degB, dres, sX);
2789 : }
2790 84777 : H = nxV_chinese_center_tree(H, P, T, ZV_chinesetree(P, T));
2791 84776 : *mod = gmael(T, lg(T)-1, 1); return gc_all(av, 2, &H, mod);
2792 : }
2793 :
2794 : GEN
2795 124940 : ZX_ZXY_resultant_worker(GEN P, GEN A, GEN B, GEN dB, GEN v)
2796 : {
2797 124940 : GEN V = cgetg(3, t_VEC);
2798 124940 : if (isintzero(dB)) dB = NULL;
2799 124940 : gel(V,1) = ZX_ZXY_resultant_slice(A, B, dB, v[1], v[2], v[3], P, &gel(V,2), v[4], v[5]);
2800 124941 : return V;
2801 : }
2802 :
2803 : GEN
2804 79187 : ZX_ZXY_resultant(GEN A, GEN B)
2805 : {
2806 79187 : pari_sp av = avma;
2807 : forprime_t S;
2808 : ulong bound;
2809 79187 : long v = fetch_var_higher();
2810 79187 : long degA = degpol(A), degB, dres = degA * degpol(B);
2811 79188 : long vX = varn(B), vY = varn(A); /* assume vY has lower priority */
2812 79188 : long sX = evalvarn(vX);
2813 : GEN worker, H, dB;
2814 79188 : B = Q_remove_denom(B, &dB);
2815 79187 : if (!dB) B = leafcopy(B);
2816 79188 : A = leafcopy(A); setvarn(A,v);
2817 79188 : B = swap_vars(B, vY, v); degB = degpol(B);
2818 79187 : bound = ZX_ZXY_ResBound(A, B, dB);
2819 79189 : if (DEBUGLEVEL>4) err_printf("bound for resultant coeffs: 2^%ld\n",bound);
2820 158378 : worker = snm_closure(is_entry("_ZX_ZXY_resultant_worker"),
2821 79189 : mkvec4(A, B, dB? dB: gen_0,
2822 : mkvecsmall5(degA, degB, dres, sX, vY)));
2823 79189 : init_modular_big(&S);
2824 79188 : H = gen_crt("ZX_ZXY_resultant_all", worker, &S, dB, bound, 0, NULL,
2825 : nxV_chinese_center, FpX_center_i);
2826 79189 : setvarn(H, vX); (void)delete_var();
2827 79189 : return gc_GEN(av, H);
2828 : }
2829 :
2830 : static long
2831 40515 : ZX_ZXY_rnfequation_lambda(GEN A, GEN B0, long lambda)
2832 : {
2833 40515 : pari_sp av = avma;
2834 40515 : long degA = degpol(A), degB, dres = degA*degpol(B0);
2835 40515 : long v = fetch_var_higher();
2836 40515 : long vX = varn(B0), vY = varn(A); /* assume vY has lower priority */
2837 40515 : long sX = evalvarn(vX);
2838 : GEN dB, B, a, b, Hp;
2839 : forprime_t S;
2840 :
2841 40515 : B0 = Q_remove_denom(B0, &dB);
2842 40515 : if (!dB) B0 = leafcopy(B0);
2843 40515 : A = leafcopy(A);
2844 40515 : B = B0;
2845 40515 : setvarn(A,v);
2846 45320 : INIT:
2847 45320 : if (lambda) B = RgX_translate(B0, monomial(stoi(lambda), 1, vY));
2848 45320 : B = swap_vars(B, vY, v);
2849 : /* B0(lambda v + x, v) */
2850 45320 : if (DEBUGLEVEL>4) err_printf("Trying lambda = %ld\n", lambda);
2851 :
2852 45320 : degB = degpol(B);
2853 45320 : init_modular_big(&S);
2854 : while (1)
2855 0 : {
2856 45320 : ulong p = u_forprime_next(&S);
2857 45320 : ulong dp = dB ? umodiu(dB, p): 1;
2858 45320 : if (!dp) continue;
2859 45320 : a = ZX_to_Flx(A, p);
2860 45320 : b = ZXX_to_FlxX(B, p, v);
2861 45319 : Hp = ZX_ZXY_resultant_prime(a, b, dp, p, degA, degB, dres, sX);
2862 45320 : if (degpol(Hp) != dres) continue;
2863 45320 : if (dp != 1) Hp = Flx_Fl_mul(Hp, Fl_powu(Fl_inv(dp,p), degA, p), p);
2864 45320 : if (!Flx_is_squarefree(Hp, p)) { lambda = next_lambda(lambda); goto INIT; }
2865 40514 : if (DEBUGLEVEL>4) err_printf("Final lambda = %ld\n", lambda);
2866 40514 : (void)delete_var(); return gc_long(av,lambda);
2867 : }
2868 : }
2869 :
2870 : GEN
2871 60554 : ZX_ZXY_rnfequation(GEN A, GEN B, long *lambda)
2872 : {
2873 60554 : if (lambda)
2874 : {
2875 40515 : *lambda = ZX_ZXY_rnfequation_lambda(A, B, *lambda);
2876 40514 : if (*lambda) B = RgX_translate(B, monomial(stoi(*lambda), 1, varn(A)));
2877 : }
2878 60553 : return ZX_ZXY_resultant(A,B);
2879 : }
2880 :
2881 : static GEN
2882 10346 : ZX_composedsum_slice(GEN A, GEN B, GEN P, GEN *mod)
2883 : {
2884 10346 : pari_sp av = avma;
2885 10346 : long i, n = lg(P)-1;
2886 : GEN H, T;
2887 10346 : if (n == 1)
2888 : {
2889 9844 : ulong p = uel(P,1);
2890 9844 : GEN a = ZX_to_Flx(A, p), b = ZX_to_Flx(B, p);
2891 9847 : GEN Hp = Flx_composedsum(a, b, p);
2892 9846 : H = gc_upto(av, Flx_to_ZX(Hp));
2893 9848 : *mod = utoipos(p); return H;
2894 : }
2895 502 : T = ZV_producttree(P);
2896 502 : A = ZX_nv_mod_tree(A, P, T);
2897 502 : B = ZX_nv_mod_tree(B, P, T);
2898 502 : H = cgetg(n+1, t_VEC);
2899 4526 : for(i=1; i <= n; i++)
2900 : {
2901 4024 : ulong p = P[i];
2902 4024 : GEN a = gel(A,i), b = gel(B,i);
2903 4024 : gel(H,i) = Flx_composedsum(a, b, p);
2904 : }
2905 502 : H = nxV_chinese_center_tree(H, P, T, ZV_chinesetree(P, T));
2906 502 : *mod = gmael(T, lg(T)-1, 1); return gc_all(av, 2, &H, mod);
2907 : }
2908 :
2909 : GEN
2910 10347 : ZX_composedsum_worker(GEN P, GEN A, GEN B)
2911 : {
2912 10347 : GEN V = cgetg(3, t_VEC);
2913 10346 : gel(V,1) = ZX_composedsum_slice(A, B, P, &gel(V,2));
2914 10349 : return V;
2915 : }
2916 :
2917 : static GEN
2918 10085 : ZX_composedsum_i(GEN A, GEN B, GEN lead)
2919 : {
2920 10085 : pari_sp av = avma;
2921 : forprime_t S;
2922 : ulong bound;
2923 : GEN H, worker, mod;
2924 10085 : if (degpol(A) < degpol(B)) swap(A, B);
2925 10084 : if (!lead) lead = mulii(leading_coeff(A),leading_coeff(B));
2926 10084 : bound = ZX_ZXY_ResBound_1(A, B);
2927 10085 : worker = snm_closure(is_entry("_ZX_composedsum_worker"), mkvec2(A,B));
2928 10086 : init_modular_big(&S);
2929 10082 : H = gen_crt("ZX_composedsum", worker, &S, lead, bound, 0, &mod,
2930 : nxV_chinese_center, FpX_center);
2931 10085 : return gc_upto(av, H);
2932 : }
2933 :
2934 : static long
2935 9697 : ZX_compositum_lambda(GEN A, GEN B, GEN lead, long lambda)
2936 : {
2937 9697 : pari_sp av = avma;
2938 : forprime_t S;
2939 : ulong p;
2940 9697 : init_modular_big(&S);
2941 9699 : p = u_forprime_next(&S);
2942 : while (1)
2943 112 : {
2944 : GEN Hp, a;
2945 9811 : if (DEBUGLEVEL>4) err_printf("Trying lambda = %ld\n", lambda);
2946 9811 : if (lead && dvdiu(lead,p)) { p = u_forprime_next(&S); continue; }
2947 9803 : a = ZX_to_Flx(ZX_rescale(A, stoi(-lambda)), p);
2948 9803 : Hp = Flx_composedsum(a, ZX_to_Flx(B, p), p);
2949 9801 : if (!Flx_is_squarefree(Hp, p)) { lambda = next_lambda(lambda); continue; }
2950 9693 : if (DEBUGLEVEL>4) err_printf("Final lambda = %ld\n", lambda);
2951 9693 : return gc_long(av, lambda);
2952 : }
2953 : }
2954 :
2955 : GEN
2956 9701 : ZX_compositum(GEN A, GEN B, long *lambda)
2957 : {
2958 9701 : GEN lead = mulii(leading_coeff(A),leading_coeff(B));
2959 9697 : if (lambda)
2960 : {
2961 9697 : *lambda = ZX_compositum_lambda(A, B, lead, *lambda);
2962 9693 : A = ZX_rescale(A, stoi(-*lambda));
2963 : }
2964 9700 : return ZX_composedsum_i(A, B, lead);
2965 : }
2966 :
2967 : GEN
2968 385 : ZX_composedsum(GEN A, GEN B)
2969 385 : { return ZX_composedsum_i(A, B, NULL); }
2970 :
2971 : static GEN
2972 359 : ZXQX_composedsum_slice(GEN A, GEN B, GEN C, GEN P, GEN *mod)
2973 : {
2974 359 : pari_sp av = avma;
2975 359 : long i, n = lg(P)-1, dC = degpol(C), v = varn(C);
2976 : GEN H, T;
2977 359 : if (n == 1)
2978 : {
2979 181 : ulong p = uel(P,1);
2980 181 : GEN a = ZXX_to_FlxX(A, p, v), b = ZXX_to_FlxX(B, p, v);
2981 181 : GEN c = ZX_to_Flx(C, p);
2982 181 : GEN Hp = FlxX_to_Flm(FlxqX_composedsum(a, b, c, p), dC);
2983 181 : H = gc_upto(av, Flm_to_ZM(Hp));
2984 181 : *mod = utoipos(p); return H;
2985 : }
2986 178 : T = ZV_producttree(P);
2987 178 : A = ZXX_nv_mod_tree(A, P, T, v);
2988 178 : B = ZXX_nv_mod_tree(B, P, T, v);
2989 178 : C = ZX_nv_mod_tree(C, P, T);
2990 178 : H = cgetg(n+1, t_VEC);
2991 660 : for(i=1; i <= n; i++)
2992 : {
2993 482 : ulong p = P[i];
2994 482 : GEN a = gel(A,i), b = gel(B,i), c = gel(C,i);
2995 482 : gel(H,i) = FlxX_to_Flm(FlxqX_composedsum(a, b, c, p), dC);
2996 : }
2997 178 : H = nmV_chinese_center_tree_seq(H, P, T, ZV_chinesetree(P, T));
2998 178 : *mod = gmael(T, lg(T)-1, 1); return gc_all(av, 2, &H, mod);
2999 : }
3000 :
3001 : GEN
3002 359 : ZXQX_composedsum_worker(GEN P, GEN A, GEN B, GEN C)
3003 : {
3004 359 : GEN V = cgetg(3, t_VEC);
3005 359 : gel(V,1) = ZXQX_composedsum_slice(A, B, C, P, &gel(V,2));
3006 359 : return V;
3007 : }
3008 :
3009 : static GEN
3010 315 : ZXQX_composedsum(GEN A, GEN B, GEN T, ulong bound)
3011 : {
3012 315 : pari_sp av = avma;
3013 : forprime_t S;
3014 : GEN H, worker, mod;
3015 315 : GEN lead = mulii(Q_content(leading_coeff(A)), Q_content(leading_coeff(B)));
3016 315 : worker = snm_closure(is_entry("_ZXQX_composedsum_worker")
3017 : , mkvec3(A,B,T));
3018 315 : init_modular_big(&S);
3019 315 : H = gen_crt("ZXQX_composedsum", worker, &S, lead, bound, 0, &mod,
3020 : nmV_chinese_center, FpM_center);
3021 315 : if (DEBUGLEVEL > 4)
3022 0 : err_printf("nfcompositum: a priori bound: %lu, a posteriori: %lu\n",
3023 : bound, expi(gsupnorm(H, DEFAULTPREC)));
3024 315 : return gc_GEN(av, RgM_to_RgXX(H, varn(A), varn(T)));
3025 : }
3026 :
3027 : static long
3028 315 : ZXQX_composedsum_bound(GEN nf, GEN A, GEN B)
3029 315 : { return ZXQX_resultant_bound_i(nf, A, B, &RgX_RgXY_ResBound_1); }
3030 :
3031 : GEN
3032 315 : nf_direct_compositum(GEN nf, GEN A, GEN B)
3033 : {
3034 315 : ulong bnd = ZXQX_composedsum_bound(nf, A, B);
3035 315 : return ZXQX_composedsum(A, B, nf_get_pol(nf), bnd);
3036 : }
3037 :
3038 : /************************************************************************
3039 : * *
3040 : * IRREDUCIBLE POLYNOMIAL / Fp *
3041 : * *
3042 : ************************************************************************/
3043 :
3044 : /* irreducible (unitary) polynomial of degree n over Fp */
3045 : GEN
3046 0 : ffinit_rand(GEN p,long n)
3047 : {
3048 0 : for(;;) {
3049 0 : pari_sp av = avma;
3050 0 : GEN pol = ZX_add(pol_xn(n, 0), random_FpX(n-1,0, p));
3051 0 : if (FpX_is_irred(pol, p)) return pol;
3052 0 : set_avma(av);
3053 : }
3054 : }
3055 :
3056 : /* return an extension of degree 2^l of F_2, assume l > 0
3057 : * Not stack clean. */
3058 : static GEN
3059 598 : ffinit_Artin_Schreier_2(long l)
3060 : {
3061 : GEN Q, T, S;
3062 : long i, v;
3063 :
3064 598 : if (l == 1) return mkvecsmall4(0,1,1,1); /*x^2 + x + 1*/
3065 549 : v = fetch_var_higher();
3066 549 : S = mkvecsmall5(0, 0, 0, 1, 1); /* y(y^2 + y) */
3067 549 : Q = mkpoln(3, pol1_Flx(0), pol1_Flx(0), S); /* x^2 + x + y(y^2+y) */
3068 549 : setvarn(Q, v);
3069 :
3070 : /* x^4+x+1, irred over F_2, minimal polynomial of a root of Q */
3071 549 : T = mkvecsmalln(6,evalvarn(v),1UL,1UL,0UL,0UL,1UL);
3072 : /* Q = x^2 + x + a(y) irred. over K = F2[y] / (T(y))
3073 : * ==> x^2 + x + a(y) b irred. over K for any root b of Q
3074 : * ==> x^2 + x + (b^2+b)b */
3075 3036 : for (i=2; i<l; i++) T = Flx_FlxY_resultant(T, Q, 2); /* minpoly(b) / F2*/
3076 551 : (void)delete_var(); T[1] = 0; return T;
3077 : }
3078 :
3079 : /* return an extension of degree p^l of F_p, assume l > 0
3080 : * Not stack clean. */
3081 : GEN
3082 955 : ffinit_Artin_Schreier(ulong p, long l)
3083 : {
3084 : long i, v;
3085 : GEN Q, R, S, T, xp;
3086 955 : if (p==2) return ffinit_Artin_Schreier_2(l);
3087 357 : xp = polxn_Flx(p,0); /* x^p */
3088 357 : T = Flx_sub(xp, mkvecsmall3(0,1,1),p); /* x^p - x - 1 */
3089 357 : if (l == 1) return T;
3090 :
3091 7 : v = evalvarn(fetch_var_higher());
3092 7 : xp[1] = v;
3093 7 : R = Flx_sub(polxn_Flx(2*p-1,0), polxn_Flx(p,0),p);
3094 7 : S = Flx_sub(xp, polx_Flx(0), p);
3095 7 : Q = FlxX_Flx_sub(Flx_to_FlxX(S, v), R, p); /* x^p - x - (y^(2p-1)-y^p) */
3096 14 : for (i = 2; i <= l; ++i) T = Flx_FlxY_resultant(T, Q, p);
3097 7 : (void)delete_var(); T[1] = 0; return T;
3098 : }
3099 :
3100 : static long
3101 149709 : flinit_check(ulong p, long n, long l)
3102 : {
3103 : ulong q;
3104 149709 : if (!uisprime(n)) return 0;
3105 102496 : q = p % n; if (!q) return 0;
3106 99885 : return ugcd((n-1)/Fl_order(q, n-1, n), l) == 1;
3107 : }
3108 :
3109 : static GEN
3110 31944 : flinit(ulong p, long l)
3111 : {
3112 31944 : ulong n = 1+l;
3113 96782 : while (!flinit_check(p,n,l)) n += l;
3114 31944 : if (DEBUGLEVEL>=4) err_printf("FFInit: using polsubcyclo(%ld, %ld)\n",n,l);
3115 31944 : return ZX_to_Flx(polsubcyclo(n,l,0), p);
3116 : }
3117 :
3118 : static GEN
3119 28989 : ffinit_fact_Flx(ulong p, long n)
3120 : {
3121 28989 : GEN P, F = factoru_pow(n), Fp = gel(F,1), Fe = gel(F,2), Fm = gel(F,3);
3122 28989 : long i, l = lg(Fm);
3123 28989 : P = cgetg(l, t_VEC);
3124 61890 : for (i = 1; i < l; i++)
3125 32899 : gel(P,i) = p==uel(Fp,i) ? ffinit_Artin_Schreier(p, Fe[i])
3126 32899 : : flinit(p, uel(Fm,i));
3127 28991 : return FlxV_composedsum(P, p);
3128 : }
3129 :
3130 : static GEN
3131 52934 : init_Flxq_i(ulong p, long n, long sv)
3132 : {
3133 : GEN P;
3134 52934 : if (!odd(p) && p != 2) pari_err_PRIME("ffinit", utoi(p));
3135 52927 : if (n == 1) return polx_Flx(sv);
3136 52927 : if (flinit_check(p, n+1, n))
3137 : {
3138 23938 : P = const_vecsmall(n+2,1);
3139 23938 : P[1] = sv; return P;
3140 : }
3141 28989 : P = ffinit_fact_Flx(p,n);
3142 28991 : P[1] = sv; return P;
3143 : }
3144 :
3145 : GEN
3146 0 : init_Flxq(ulong p, long n, long v)
3147 : {
3148 0 : pari_sp av = avma;
3149 0 : return gc_upto(av, init_Flxq_i(p, n, v));
3150 : }
3151 :
3152 : /* check if polsubcyclo(n,l,0) is irreducible modulo p */
3153 : static long
3154 8207 : fpinit_check(GEN p, long n, long l)
3155 : {
3156 : ulong q;
3157 8207 : if (!uisprime(n)) return 0;
3158 4842 : q = umodiu(p,n); if (!q) return 0;
3159 4842 : return ugcd((n-1)/Fl_order(q, n-1, n), l) == 1;
3160 : }
3161 :
3162 : /* let k=2 if p%4==1, and k=4 else and assume k*p does not divide l.
3163 : * Return an irreducible polynomial of degree l over F_p.
3164 : * Variant of Adleman and Lenstra "Finding irreducible polynomials over
3165 : * finite fields", ACM, 1986 (5) 350--355.
3166 : * Not stack clean */
3167 : static GEN
3168 1828 : fpinit(GEN p, long l)
3169 : {
3170 1828 : ulong n = 1+l;
3171 6168 : while (!fpinit_check(p,n,l)) n += l;
3172 1828 : if (DEBUGLEVEL>=4) err_printf("FFInit: using polsubcyclo(%ld, %ld)\n",n,l);
3173 1828 : return FpX_red(polsubcyclo(n,l,0),p);
3174 : }
3175 :
3176 : static GEN
3177 1637 : ffinit_fact(GEN p, long n)
3178 : {
3179 1637 : GEN P, F = factoru_pow(n), Fp = gel(F,1), Fe = gel(F,2), Fm = gel(F,3);
3180 1637 : long i, l = lg(Fm);
3181 1637 : P = cgetg(l, t_VEC);
3182 3465 : for (i = 1; i < l; ++i)
3183 3656 : gel(P,i) = absequaliu(p, Fp[i]) ?
3184 0 : Flx_to_ZX(ffinit_Artin_Schreier(Fp[i], Fe[i]))
3185 1828 : : fpinit(p, Fm[i]);
3186 1637 : return FpXV_composedsum(P, p);
3187 : }
3188 :
3189 : static GEN
3190 55239 : init_Fq_i(GEN p, long n, long v)
3191 : {
3192 : GEN P;
3193 55239 : if (n <= 0) pari_err_DOMAIN("ffinit", "degree", "<=", gen_0, stoi(n));
3194 55239 : if (typ(p) != t_INT) pari_err_TYPE("ffinit",p);
3195 55239 : if (cmpiu(p, 2) < 0) pari_err_PRIME("ffinit",p);
3196 55232 : if (v < 0) v = 0;
3197 55232 : if (n == 1) return pol_x(v);
3198 54980 : if (lgefint(p) == 3)
3199 52934 : return Flx_to_ZX(init_Flxq_i(p[2], n, evalvarn(v)));
3200 2046 : if (!mpodd(p)) pari_err_PRIME("ffinit", p);
3201 2039 : if (fpinit_check(p, n+1, n)) return polcyclo(n+1, v);
3202 1637 : P = ffinit_fact(p,n);
3203 1637 : setvarn(P, v); return P;
3204 : }
3205 : GEN
3206 54672 : init_Fq(GEN p, long n, long v)
3207 : {
3208 54672 : pari_sp av = avma;
3209 54672 : return gc_upto(av, init_Fq_i(p, n, v));
3210 : }
3211 : GEN
3212 567 : ffinit(GEN p, long n, long v)
3213 : {
3214 567 : pari_sp av = avma;
3215 567 : return gc_upto(av, FpX_to_mod(init_Fq_i(p, n, v), p));
3216 : }
3217 :
3218 : GEN
3219 3178 : ffnbirred(GEN p, long n)
3220 : {
3221 3178 : pari_sp av = avma;
3222 3178 : GEN s = powiu(p,n), F = factoru(n), D = divisorsu_moebius(gel(F, 1));
3223 3178 : long j, l = lg(D);
3224 6797 : for (j = 2; j < l; j++) /* skip d = 1 */
3225 : {
3226 3619 : long md = D[j]; /* mu(d) * d, d squarefree */
3227 3619 : GEN pd = powiu(p, n / labs(md)); /* p^{n/d} */
3228 3619 : s = md > 0? addii(s, pd): subii(s,pd);
3229 : }
3230 3178 : return gc_INT(av, diviuexact(s, n));
3231 : }
3232 :
3233 : GEN
3234 616 : ffsumnbirred(GEN p, long n)
3235 : {
3236 616 : pari_sp av = avma, av2;
3237 616 : GEN q, t = p, v = vecfactoru_i(1, n);
3238 : long i;
3239 616 : q = cgetg(n+1,t_VEC); gel(q,1) = p;
3240 1764 : for (i=2; i<=n; i++) gel(q,i) = mulii(gel(q,i-1), p);
3241 616 : av2 = avma;
3242 1764 : for (i=2; i<=n; i++)
3243 : {
3244 1148 : GEN s = gel(q,i), F = gel(v,i), D = divisorsu_moebius(gel(F,1));
3245 1148 : long j, l = lg(D);
3246 2534 : for (j = 2; j < l; j++) /* skip 1 */
3247 : {
3248 1386 : long md = D[j];
3249 1386 : GEN pd = gel(q, i / labs(md)); /* p^{i/d} */
3250 1386 : s = md > 0? addii(s, pd): subii(s, pd);
3251 : }
3252 1148 : t = gc_INT(av2, addii(t, diviuexact(s, i)));
3253 : }
3254 616 : return gc_INT(av, t);
3255 : }
3256 :
3257 : GEN
3258 140 : ffnbirred0(GEN p, long n, long flag)
3259 : {
3260 140 : if (typ(p) != t_INT) pari_err_TYPE("ffnbirred", p);
3261 140 : if (n <= 0) pari_err_DOMAIN("ffnbirred", "degree", "<=", gen_0, stoi(n));
3262 140 : switch(flag)
3263 : {
3264 70 : case 0: return ffnbirred(p, n);
3265 70 : case 1: return ffsumnbirred(p, n);
3266 : }
3267 0 : pari_err_FLAG("ffnbirred");
3268 : return NULL; /* LCOV_EXCL_LINE */
3269 : }
3270 :
3271 : static void
3272 2261 : checkmap(GEN m, const char *s)
3273 : {
3274 2261 : if (typ(m)!=t_VEC || lg(m)!=3 || typ(gel(m,1))!=t_FFELT)
3275 0 : pari_err_TYPE(s,m);
3276 2261 : }
3277 :
3278 : GEN
3279 189 : ffembed(GEN a, GEN b)
3280 : {
3281 189 : pari_sp av = avma;
3282 189 : GEN p, Ta, Tb, g, r = NULL;
3283 189 : if (typ(a)!=t_FFELT) pari_err_TYPE("ffembed",a);
3284 189 : if (typ(b)!=t_FFELT) pari_err_TYPE("ffembed",b);
3285 189 : p = FF_p_i(a); g = FF_gen(a);
3286 189 : if (!equalii(p, FF_p_i(b))) pari_err_MODULUS("ffembed",a,b);
3287 189 : Ta = FF_mod(a);
3288 189 : Tb = FF_mod(b);
3289 189 : if (degpol(Tb)%degpol(Ta)!=0)
3290 7 : pari_err_DOMAIN("ffembed",GENtostr_raw(a),"is not a subfield of",b,a);
3291 182 : r = gel(FFX_roots(Ta, b), 1);
3292 182 : return gc_GEN(av, mkvec2(g,r));
3293 : }
3294 :
3295 : GEN
3296 91 : ffextend(GEN a, GEN P, long v)
3297 : {
3298 91 : pari_sp av = avma;
3299 : long n;
3300 : GEN p, T, R, g, m;
3301 91 : if (typ(a)!=t_FFELT) pari_err_TYPE("ffextend",a);
3302 91 : T = a; p = FF_p_i(a);
3303 91 : if (typ(P)!=t_POL || !RgX_is_FpXQX(P,&T,&p)) pari_err_TYPE("ffextend", P);
3304 49 : if (!FF_samefield(a, T)) pari_err_MODULUS("ffextend",a,T);
3305 49 : if (v < 0) v = varn(P);
3306 49 : n = FF_f(T) * degpol(P); R = ffinit(p, n, v); g = ffgen(R, v);
3307 49 : m = ffembed(a, g);
3308 49 : R = FFX_roots(ffmap(m, P),g);
3309 49 : return gc_GEN(av, mkvec2(gel(R,1), m));
3310 : }
3311 :
3312 : GEN
3313 42 : fffrobenius(GEN a, long n)
3314 : {
3315 42 : if (typ(a)!=t_FFELT) pari_err_TYPE("fffrobenius",a);
3316 42 : retmkvec2(FF_gen(a), FF_Frobenius(a, n));
3317 : }
3318 :
3319 : GEN
3320 133 : ffinvmap(GEN m)
3321 : {
3322 133 : pari_sp av = avma;
3323 : long i, l;
3324 133 : GEN T, F, a, g, r, f = NULL;
3325 133 : checkmap(m, "ffinvmap");
3326 133 : a = gel(m,1); r = gel(m,2);
3327 133 : if (typ(r) != t_FFELT)
3328 7 : pari_err_TYPE("ffinvmap", m);
3329 126 : g = FF_gen(a);
3330 126 : T = FF_mod(r);
3331 126 : F = gel(FFX_factor(T, a), 1);
3332 126 : l = lg(F);
3333 490 : for(i=1; i<l; i++)
3334 : {
3335 490 : GEN s = FFX_rem(FF_to_FpXQ_i(r), gel(F, i), a);
3336 490 : if (degpol(s)==0 && gequal(constant_coeff(s),g)) { f = gel(F, i); break; }
3337 : }
3338 126 : if (f==NULL) pari_err_TYPE("ffinvmap", m);
3339 126 : if (degpol(f)==1) f = FF_neg_i(gel(f,2));
3340 126 : return gc_GEN(av, mkvec2(FF_gen(r),f));
3341 : }
3342 :
3343 : static GEN
3344 1260 : ffpartmapimage(const char *s, GEN r)
3345 : {
3346 1260 : GEN a = NULL, p = NULL;
3347 1260 : if (typ(r)==t_POL && degpol(r) >= 1
3348 1260 : && RgX_is_FpXQX(r,&a,&p) && a && typ(a)==t_FFELT) return a;
3349 0 : pari_err_TYPE(s, r);
3350 : return NULL; /* LCOV_EXCL_LINE */
3351 : }
3352 :
3353 : static GEN
3354 2709 : ffeltmap_i(GEN m, GEN x)
3355 : {
3356 2709 : GEN r = gel(m,2);
3357 2709 : if (!FF_samefield(x, gel(m,1)))
3358 84 : pari_err_DOMAIN("ffmap","m","domain does not contain", x, r);
3359 2625 : if (typ(r)==t_FFELT)
3360 1659 : return FF_map(r, x);
3361 : else
3362 966 : return FFX_preimage(x, r, ffpartmapimage("ffmap", r));
3363 : }
3364 :
3365 : static GEN
3366 4459 : ffmap_i(GEN m, GEN x)
3367 : {
3368 : GEN y;
3369 4459 : long i, lx, tx = typ(x);
3370 4459 : switch(tx)
3371 : {
3372 2541 : case t_FFELT:
3373 2541 : return ffeltmap_i(m, x);
3374 1267 : case t_POL: case t_RFRAC: case t_SER:
3375 : case t_VEC: case t_COL: case t_MAT:
3376 1267 : y = cgetg_copy(x, &lx);
3377 1988 : for (i = 1; i < lontyp[tx]; i++) y[i] = x[i];
3378 4564 : for (; i < lx; i++)
3379 : {
3380 3339 : GEN yi = ffmap_i(m, gel(x,i));
3381 3297 : if (!yi) return NULL;
3382 3297 : gel(y,i) = yi;
3383 : }
3384 1225 : return y;
3385 : }
3386 651 : return gcopy(x);
3387 : }
3388 :
3389 : GEN
3390 1036 : ffmap(GEN m, GEN x)
3391 : {
3392 1036 : pari_sp ltop = avma;
3393 : GEN y;
3394 1036 : checkmap(m, "ffmap");
3395 1036 : y = ffmap_i(m, x);
3396 1036 : if (y) return y;
3397 42 : retgc_const(ltop, cgetg(1, t_VEC));
3398 : }
3399 :
3400 : static GEN
3401 252 : ffeltmaprel_i(GEN m, GEN x)
3402 : {
3403 252 : GEN g = gel(m,1), r = gel(m,2);
3404 252 : if (!FF_samefield(x, g))
3405 0 : pari_err_DOMAIN("ffmap","m","domain does not contain", x, r);
3406 252 : if (typ(r)==t_FFELT)
3407 84 : retmkpolmod(FF_map(r, x), pol_x(FF_var(g)));
3408 : else
3409 168 : retmkpolmod(FFX_preimagerel(x, r, ffpartmapimage("ffmap", r)), gcopy(r));
3410 : }
3411 :
3412 : static GEN
3413 252 : ffmaprel_i(GEN m, GEN x)
3414 : {
3415 252 : switch(typ(x))
3416 : {
3417 252 : case t_FFELT:
3418 252 : return ffeltmaprel_i(m, x);
3419 0 : case t_POL: pari_APPLY_pol_normalized(ffmaprel_i(m, gel(x,i)));
3420 0 : case t_SER: pari_APPLY_ser_normalized(ffmaprel_i(m, gel(x,i)));
3421 0 : case t_RFRAC: case t_VEC: case t_COL: case t_MAT:
3422 0 : pari_APPLY_same(ffmaprel_i(m, gel(x,i)));
3423 : }
3424 0 : return gcopy(x);
3425 : }
3426 : GEN
3427 252 : ffmaprel(GEN m, GEN x) { checkmap(m, "ffmaprel"); return ffmaprel_i(m, x); }
3428 :
3429 : static void
3430 84 : err_compo(GEN m, GEN n)
3431 84 : { pari_err_DOMAIN("ffcompomap","m","domain does not contain codomain of",n,m); }
3432 :
3433 : GEN
3434 420 : ffcompomap(GEN m, GEN n)
3435 : {
3436 420 : pari_sp av = avma;
3437 420 : GEN g = gel(n,1), r, m2, n2;
3438 420 : checkmap(m, "ffcompomap");
3439 420 : checkmap(n, "ffcompomap");
3440 420 : m2 = gel(m,2); n2 = gel(n,2);
3441 420 : switch((typ(m2)==t_POL)|((typ(n2)==t_POL)<<1))
3442 : {
3443 84 : case 0:
3444 84 : if (!FF_samefield(gel(m,1),n2)) err_compo(m,n);
3445 42 : r = FF_map(gel(m,2), n2);
3446 42 : break;
3447 84 : case 2:
3448 84 : r = ffmap_i(m, n2);
3449 42 : if (lg(r) == 1) err_compo(m,n);
3450 42 : break;
3451 168 : case 1:
3452 168 : r = ffeltmap_i(m, n2);
3453 126 : if (!r)
3454 : {
3455 : GEN a, A, R, M;
3456 : long dm, dn;
3457 42 : a = ffpartmapimage("ffcompomap",m2);
3458 42 : A = FF_to_FpXQ_i(FF_neg(n2));
3459 42 : setvarn(A, 1);
3460 42 : R = deg1pol(gen_1, A, 0);
3461 42 : setvarn(R, 0);
3462 42 : M = gcopy(m2);
3463 42 : setvarn(M, 1);
3464 42 : r = polresultant0(R, M, 1, 0);
3465 42 : dm = FF_f(gel(m,1)); dn = FF_f(gel(n,1));
3466 42 : if (dm % dn || !FFX_ispower(r, dm/dn, a, &r)) err_compo(m,n);
3467 42 : setvarn(r, varn(FF_mod(g)));
3468 : }
3469 126 : break;
3470 84 : case 3:
3471 : {
3472 : GEN M, R, T, p, a;
3473 84 : a = ffpartmapimage("ffcompomap",n2);
3474 84 : if (!FF_samefield(a, gel(m,1))) err_compo(m,n);
3475 42 : p = FF_p_i(gel(n,1));
3476 42 : T = FF_mod(gel(n,1));
3477 42 : setvarn(T, 1);
3478 42 : R = RgX_to_FpXQX(n2,T,p);
3479 42 : setvarn(R, 0);
3480 42 : M = gcopy(m2);
3481 42 : setvarn(M, 1);
3482 42 : r = polresultant0(R, M, 1, 0);
3483 42 : setvarn(r, varn(n2));
3484 : }
3485 : }
3486 252 : return gc_GEN(av, mkvec2(g,r));
3487 : }
|