Line data Source code
1 : /* Copyright (C) 2000-2003 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 : #include "pari.h"
16 : #include "paripriv.h"
17 :
18 : /*************************************************************************/
19 : /** **/
20 : /** Routines for handling FFELT **/
21 : /** **/
22 : /*************************************************************************/
23 :
24 : /*************************************************************************/
25 : /** **/
26 : /** Low-level constructors **/
27 : /** **/
28 : /*************************************************************************/
29 :
30 : INLINE void
31 30919531 : _getFF(GEN x, GEN *T, GEN *p, ulong *pp)
32 : {
33 30919531 : *T=gel(x,3);
34 30919531 : *p=gel(x,4);
35 30919531 : *pp=(*p)[2];
36 30919531 : }
37 :
38 : INLINE GEN
39 29611893 : _initFF(GEN x, GEN *T, GEN *p, ulong *pp)
40 : {
41 29611893 : _getFF(x,T,p,pp);
42 29611893 : return cgetg(5,t_FFELT);
43 : }
44 :
45 : INLINE void
46 23392481 : _checkFF(GEN x, GEN y, const char *s)
47 23392481 : { if (!FF_samefield(x,y)) pari_err_OP(s,x,y); }
48 :
49 : INLINE GEN
50 29410244 : _mkFF(GEN x, GEN z, GEN r)
51 : {
52 29410244 : z[1]=x[1];
53 29410244 : gel(z,2)=r;
54 29410244 : gel(z,3)=gcopy(gel(x,3));
55 29410244 : gel(z,4)=icopy(gel(x,4));
56 29410244 : return z;
57 : }
58 :
59 : INLINE GEN
60 2883235 : _mkFF_i(GEN x, GEN z, GEN r)
61 : {
62 2883235 : z[1]=x[1];
63 2883235 : gel(z,2)=r;
64 2883235 : gel(z,3)=gel(x,3);
65 2883235 : gel(z,4)=gel(x,4);
66 2883235 : return z;
67 : }
68 :
69 : INLINE GEN
70 2687539 : mkFF_i(GEN x, GEN r)
71 : {
72 2687539 : GEN z = cgetg(5,t_FFELT);
73 2687539 : return _mkFF_i(x,z,r);
74 : }
75 :
76 : /*************************************************************************/
77 : /** **/
78 : /** medium-level constructors **/
79 : /** **/
80 : /*************************************************************************/
81 :
82 : static GEN
83 369244 : Z_to_raw(GEN x, GEN ff)
84 : {
85 : ulong pp;
86 : GEN T, p;
87 369244 : _getFF(ff,&T,&p,&pp);
88 369244 : switch(ff[1])
89 : {
90 10248 : case t_FF_FpXQ:
91 10248 : return scalarpol(x, varn(T));
92 173668 : case t_FF_F2xq:
93 173668 : return Z_to_F2x(x, T[1]);
94 185328 : default:
95 185328 : return Z_to_Flx(x, pp, T[1]);
96 : }
97 : }
98 :
99 : static GEN
100 1895837 : Rg_to_raw(GEN x, GEN ff)
101 : {
102 1895837 : long tx = typ(x);
103 1895837 : switch(tx)
104 : {
105 369244 : case t_INT: case t_FRAC: case t_PADIC: case t_INTMOD:
106 369244 : return Z_to_raw(Rg_to_Fp(x, FF_p_i(ff)), ff);
107 1526593 : case t_FFELT:
108 1526593 : if (!FF_samefield(x,ff))
109 0 : pari_err_MODULUS("Rg_to_raw",x,ff);
110 1526593 : return gel(x,2);
111 : }
112 0 : pari_err_TYPE("Rg_to_raw",x);
113 : return NULL;/* LCOV_EXCL_LINE */
114 : }
115 :
116 : static GEN
117 278144 : FFX_to_raw(GEN x, GEN ff)
118 : {
119 : long i, lx;
120 278144 : GEN y = cgetg_copy(x,&lx);
121 278144 : y[1] = x[1];
122 1615297 : for(i=2; i<lx; i++)
123 1337153 : gel(y, i) = Rg_to_raw(gel(x, i), ff);
124 278144 : switch (ff[1])
125 : {
126 4996 : case t_FF_FpXQ:
127 4996 : return FpXX_renormalize(y, lx);
128 120111 : case t_FF_F2xq:
129 120111 : return F2xX_renormalize(y, lx);
130 153037 : default:
131 153037 : return FlxX_renormalize(y, lx);
132 : }
133 : }
134 :
135 : static GEN
136 596646 : FFC_to_raw(GEN x, GEN ff) { pari_APPLY_same(Rg_to_raw(gel(x, i), ff)) }
137 : static GEN
138 46314 : FFM_to_raw(GEN x, GEN ff) { pari_APPLY_same(FFC_to_raw(gel(x, i), ff)) }
139 :
140 : /* in place */
141 : static GEN
142 568645 : rawFq_to_FF(GEN x, GEN ff)
143 : {
144 568645 : return mkFF_i(ff, typ(x)==t_INT ? scalarpol(x, varn(gel(ff,3))): x);
145 : }
146 :
147 : /* in place */
148 : static GEN
149 94581 : raw_to_FFX(GEN x, GEN ff)
150 : {
151 94581 : long i, lx = lg(x);
152 663226 : for (i=2; i<lx; i++) gel(x,i) = rawFq_to_FF(gel(x,i), ff);
153 94581 : return x;
154 : }
155 :
156 : /* in place */
157 : static GEN
158 106626 : raw_to_FFC(GEN x, GEN ff)
159 : {
160 106626 : long i, lx = lg(x);
161 454890 : for (i=1; i<lx; i++) gel(x,i) = mkFF_i(ff, gel(x,i));
162 106626 : return x;
163 : }
164 :
165 : /* in place */
166 : static GEN
167 4152 : raw_to_FFM(GEN x, GEN ff)
168 : {
169 4152 : long i, lx = lg(x);
170 21738 : for (i=1; i<lx; i++) gel(x,i) = raw_to_FFC(gel(x,i), ff);
171 4152 : return x;
172 : }
173 :
174 : GEN
175 123204 : Fq_to_FF(GEN x, GEN ff)
176 : {
177 : ulong pp;
178 123204 : GEN r, T, p, z = _initFF(ff,&T,&p,&pp);
179 123204 : if (typ(x) == t_INT) switch(ff[1])
180 : {
181 2 : case t_FF_FpXQ: r = scalarpol(x, varn(T)); break;
182 1404 : case t_FF_F2xq: r = Z_to_F2x(x,T[1]); break;
183 4084 : default: r = Z_to_Flx(x,pp,T[1]);
184 : }
185 117714 : else switch(ff[1])
186 : {
187 5958 : case t_FF_FpXQ: r = ZX_copy(x); setvarn(r, varn(T)); break;
188 39882 : case t_FF_F2xq: r = ZX_to_F2x(x); r[1] = T[1]; break;
189 71874 : default: r = ZX_to_Flx(x,pp); r[1] = T[1];
190 : }
191 123204 : return _mkFF_i(ff, z, r);
192 : }
193 :
194 : /*************************************************************************/
195 : /** **/
196 : /** Public functions **/
197 : /** **/
198 : /*************************************************************************/
199 :
200 : /* Return true if x and y are defined in the same field */
201 :
202 : static int
203 26452426 : FF_samechar(GEN x, GEN y)
204 26452426 : { return x[1] == y[1] && equalii(gel(x,4),gel(y,4)); }
205 :
206 : int
207 26452426 : FF_samefield(GEN x, GEN y)
208 26452426 : { return FF_samechar(x, y) && gidentical(gel(x,3),gel(y,3)); }
209 :
210 : int
211 57599 : FF_equal(GEN x, GEN y)
212 57599 : { return FF_samefield(x,y) && gidentical(gel(x,2),gel(y,2)); }
213 :
214 : int
215 7748993 : FF_equal0(GEN x)
216 : {
217 7748993 : return lgpol(gel(x,2))==0;
218 : }
219 :
220 : int
221 20152 : FF_equal1(GEN x)
222 : {
223 20152 : GEN A = gel(x,2);
224 20152 : switch(x[1])
225 : {
226 4500 : case t_FF_FpXQ:
227 4500 : return degpol(A)==0 && gequal1(gel(A,2));
228 15652 : default:
229 15652 : return degpol(A)==0 && A[2]==1;
230 : }
231 : }
232 :
233 : static int
234 0 : Fp_cmp_1(GEN x, GEN p)
235 0 : { pari_sp av = avma; return gc_bool(av, equalii(x, addis(p,-1))); }
236 :
237 : int
238 36 : FF_equalm1(GEN x)
239 : {
240 : ulong pp;
241 36 : GEN T, p, y = gel(x,2);
242 36 : _getFF(x,&T,&p,&pp);
243 36 : switch(x[1])
244 : {
245 0 : case t_FF_FpXQ:
246 0 : return (degpol(y) == 0 && Fp_cmp_1(gel(y,2), p));
247 36 : default:
248 36 : return (degpol(y) == 0 && uel(y,2) == pp-1);
249 : }
250 : }
251 :
252 : GEN
253 5069 : FF_zero(GEN x)
254 : {
255 : ulong pp;
256 5069 : GEN r, T, p, z=_initFF(x,&T,&p,&pp);
257 5069 : switch(x[1])
258 : {
259 407 : case t_FF_FpXQ:
260 407 : r=zeropol(varn(T));
261 407 : break;
262 1124 : case t_FF_F2xq:
263 1124 : r=zero_F2x(T[1]);
264 1124 : break;
265 3538 : default:
266 3538 : r=zero_Flx(T[1]);
267 : }
268 5069 : return _mkFF(x,z,r);
269 : }
270 :
271 : GEN
272 19134 : FF_1(GEN x)
273 : {
274 : ulong pp;
275 19134 : GEN r, T, p, z=_initFF(x,&T,&p,&pp);
276 19134 : switch(x[1])
277 : {
278 813 : case t_FF_FpXQ:
279 813 : r=pol_1(varn(T));
280 813 : break;
281 9705 : case t_FF_F2xq:
282 9705 : r=pol1_F2x(T[1]);
283 9705 : break;
284 8616 : default:
285 8616 : r=pol1_Flx(T[1]);
286 : }
287 19134 : return _mkFF(x,z,r);
288 : }
289 :
290 : GEN
291 656 : FF_gen(GEN x)
292 : {
293 : ulong pp;
294 656 : GEN r, T, p, z=_initFF(x,&T,&p,&pp);
295 656 : switch(x[1])
296 : {
297 139 : case t_FF_FpXQ:
298 139 : r = pol_x(varn(T));
299 139 : if (degpol(T)==1) r = FpX_rem(r, T, p);
300 139 : break;
301 180 : case t_FF_F2xq:
302 180 : r = polx_F2x(T[1]);
303 180 : if (F2x_degree(T)==1) r = F2x_rem(r, T);
304 180 : break;
305 337 : default:
306 337 : r = polx_Flx(T[1]);
307 337 : if (degpol(T)==1) r = Flx_rem(r, T, pp);
308 : }
309 656 : return _mkFF(x,z,r);
310 : }
311 : GEN
312 70026 : FF_q(GEN x)
313 : {
314 : ulong pp;
315 : GEN T, p;
316 70026 : _getFF(x,&T,&p,&pp);
317 70026 : switch(x[1])
318 : {
319 2167 : case t_FF_FpXQ:
320 2167 : return powiu(p, degpol(T));
321 : break;
322 9048 : case t_FF_F2xq:
323 9048 : return int2n(F2x_degree(T));
324 : break;
325 58811 : default:
326 58811 : return powuu(pp,degpol(T));
327 : }
328 : }
329 :
330 : GEN
331 96 : FF_p(GEN x)
332 : {
333 96 : return icopy(gel(x,4));
334 : }
335 :
336 : GEN
337 2258741 : FF_p_i(GEN x)
338 : {
339 2258741 : return gel(x,4);
340 : }
341 :
342 : GEN
343 142998 : FF_mod(GEN x)
344 : {
345 142998 : switch(x[1])
346 : {
347 348 : case t_FF_FpXQ:
348 348 : return ZX_copy(gel(x,3));
349 348 : case t_FF_F2xq:
350 348 : return F2x_to_ZX(gel(x,3));
351 142302 : default:
352 142302 : return Flx_to_ZX(gel(x,3));
353 : }
354 : }
355 :
356 : long
357 120 : FF_var(GEN x)
358 : {
359 120 : switch(x[1])
360 : {
361 48 : case t_FF_FpXQ:
362 48 : return varn(gel(x,3));
363 72 : case t_FF_F2xq:
364 : default:
365 72 : return gel(x,3)[1]>>VARNSHIFT;
366 : }
367 : }
368 :
369 : long
370 318 : FF_f(GEN x)
371 : {
372 318 : switch(x[1])
373 : {
374 126 : case t_FF_F2xq:
375 126 : return F2x_degree(gel(x,3));
376 192 : default:
377 192 : return degpol(gel(x,3));
378 : }
379 : }
380 :
381 : GEN
382 783517 : FF_to_F2xq(GEN x)
383 : {
384 783517 : switch(x[1])
385 : {
386 0 : case t_FF_FpXQ:
387 0 : return ZX_to_F2x(gel(x,2));
388 783517 : case t_FF_F2xq:
389 783517 : return zv_copy(gel(x,2));
390 0 : default:
391 0 : return Flx_to_F2x(gel(x,2));
392 : }
393 : }
394 :
395 : GEN
396 0 : FF_to_F2xq_i(GEN x)
397 : {
398 0 : switch(x[1])
399 : {
400 0 : case t_FF_FpXQ:
401 0 : return ZX_to_F2x(gel(x,2));
402 0 : case t_FF_F2xq:
403 0 : return gel(x,2);
404 0 : default:
405 0 : return Flx_to_F2x(gel(x,2));
406 : }
407 : }
408 :
409 : GEN
410 626556 : FF_to_Flxq(GEN x)
411 : {
412 626556 : switch(x[1])
413 : {
414 0 : case t_FF_FpXQ:
415 0 : return ZX_to_Flx(gel(x,2),itou(gel(x,4)));
416 0 : case t_FF_F2xq:
417 0 : return F2x_to_Flx(gel(x,2));
418 626556 : default:
419 626556 : return zv_copy(gel(x,2));
420 : }
421 : }
422 :
423 : GEN
424 0 : FF_to_Flxq_i(GEN x)
425 : {
426 0 : switch(x[1])
427 : {
428 0 : case t_FF_FpXQ:
429 0 : return ZX_to_Flx(gel(x,2),itou(gel(x,4)));
430 0 : case t_FF_F2xq:
431 0 : return F2x_to_Flx(gel(x,2));
432 0 : default:
433 0 : return gel(x,2);
434 : }
435 : }
436 :
437 : GEN
438 17732 : FF_to_FpXQ(GEN x)
439 : {
440 17732 : switch(x[1])
441 : {
442 14376 : case t_FF_FpXQ:
443 14376 : return ZX_copy(gel(x,2));
444 522 : case t_FF_F2xq:
445 522 : return F2x_to_ZX(gel(x,2));
446 2834 : default:
447 2834 : return Flx_to_ZX(gel(x,2));
448 : }
449 : }
450 :
451 : GEN
452 147546 : FF_to_FpXQ_i(GEN x)
453 : {
454 147546 : switch(x[1])
455 : {
456 1055 : case t_FF_FpXQ:
457 1055 : return gel(x,2);
458 1249 : case t_FF_F2xq:
459 1249 : return F2x_to_ZX(gel(x,2));
460 145242 : default:
461 145242 : return Flx_to_ZX(gel(x,2));
462 : }
463 : }
464 :
465 : GEN
466 11254891 : FF_add(GEN x, GEN y)
467 : {
468 : ulong pp;
469 11254891 : GEN r, T, p, z=_initFF(x,&T,&p,&pp);
470 11254891 : _checkFF(x,y,"+");
471 11254891 : switch(x[1])
472 : {
473 505365 : case t_FF_FpXQ:
474 505365 : r=FpX_add(gel(x,2),gel(y,2),p);
475 505365 : break;
476 964790 : case t_FF_F2xq:
477 964790 : r=F2x_add(gel(x,2),gel(y,2));
478 964790 : break;
479 9784736 : default:
480 9784736 : r=Flx_add(gel(x,2),gel(y,2),pp);
481 : }
482 11254891 : return _mkFF(x,z,r);
483 : }
484 : GEN
485 200266 : FF_sub(GEN x, GEN y)
486 : {
487 : ulong pp;
488 200266 : GEN r, T, p, z=_initFF(x,&T,&p,&pp);
489 200266 : _checkFF(x,y,"+");
490 200266 : switch(x[1])
491 : {
492 1023 : case t_FF_FpXQ:
493 1023 : r=FpX_sub(gel(x,2),gel(y,2),p);
494 1023 : break;
495 156539 : case t_FF_F2xq:
496 156539 : r=F2x_add(gel(x,2),gel(y,2));
497 156539 : break;
498 42704 : default:
499 42704 : r=Flx_sub(gel(x,2),gel(y,2),pp);
500 : }
501 200266 : return _mkFF(x,z,r);
502 : }
503 :
504 : GEN
505 1168347 : FF_Z_add(GEN x, GEN y)
506 : {
507 : ulong pp;
508 1168347 : GEN r, T, p, z=_initFF(x,&T,&p,&pp);
509 1168347 : switch(x[1])
510 : {
511 2400 : case t_FF_FpXQ:
512 : {
513 2400 : pari_sp av=avma;
514 2400 : r=gc_upto(av,FpX_Fp_add(gel(x,2),modii(y,p),p));
515 2400 : break;
516 : }
517 628394 : case t_FF_F2xq:
518 628394 : r=mpodd(y)?F2x_1_add(gel(x,2)):vecsmall_copy(gel(x,2));
519 628394 : break;
520 537553 : default:
521 537553 : r=Flx_Fl_add(gel(x,2),umodiu(y,pp),pp);
522 : }
523 1168347 : return _mkFF(x,z,r);
524 : }
525 :
526 : GEN
527 1146 : FF_Q_add(GEN x, GEN y)
528 : {
529 : ulong pp;
530 1146 : GEN r, T, p, z=_initFF(x,&T,&p,&pp);
531 1146 : switch(x[1])
532 : {
533 1 : case t_FF_FpXQ:
534 : {
535 1 : pari_sp av=avma;
536 1 : r=gc_upto(av,FpX_Fp_add(gel(x,2),Rg_to_Fp(y,p),p));
537 1 : break;
538 : }
539 6 : case t_FF_F2xq:
540 6 : r=Rg_to_Fl(y,pp)?F2x_1_add(gel(x,2)):vecsmall_copy(gel(x,2));
541 6 : break;
542 1139 : default:
543 1139 : r=Flx_Fl_add(gel(x,2),Rg_to_Fl(y,pp),pp);
544 : }
545 1146 : return _mkFF(x,z,r);
546 : }
547 :
548 : GEN
549 69558 : FF_neg(GEN x)
550 : {
551 : ulong pp;
552 69558 : GEN r, T, p, z=_initFF(x,&T,&p,&pp);
553 69558 : switch(x[1])
554 : {
555 4052 : case t_FF_FpXQ:
556 4052 : r=FpX_neg(gel(x,2),p);
557 4052 : break;
558 14673 : case t_FF_F2xq:
559 14673 : r=vecsmall_copy(gel(x,2));
560 14673 : break;
561 50833 : default:
562 50833 : r=Flx_neg(gel(x,2),pp);
563 : }
564 69558 : return _mkFF(x,z,r);
565 : }
566 :
567 : GEN
568 72492 : FF_neg_i(GEN x)
569 : {
570 : ulong pp;
571 72492 : GEN r, T, p, z=_initFF(x,&T,&p,&pp);
572 72492 : switch(x[1])
573 : {
574 1511 : case t_FF_FpXQ:
575 1511 : r=FpX_neg(gel(x,2),p);
576 1511 : break;
577 7326 : case t_FF_F2xq:
578 7326 : r=gel(x,2);
579 7326 : break;
580 63655 : default:
581 63655 : r=Flx_neg(gel(x,2),pp);
582 : }
583 72492 : return _mkFF_i(x,z,r);
584 : }
585 :
586 : GEN
587 1530 : FF_map(GEN m, GEN x)
588 : {
589 : ulong pp;
590 1530 : GEN r, T, p, z=_initFF(m,&T,&p,&pp);
591 1530 : switch(m[1])
592 : {
593 504 : case t_FF_FpXQ:
594 504 : r=FpX_FpXQ_eval(gel(x,2),gel(m,2),T,p);
595 504 : break;
596 600 : case t_FF_F2xq:
597 600 : r=F2x_F2xq_eval(gel(x,2),gel(m,2),T);
598 600 : break;
599 426 : default:
600 426 : r=Flx_Flxq_eval(gel(x,2),gel(m,2),T,pp);
601 : }
602 1530 : return _mkFF(m,z,r);
603 : }
604 :
605 : GEN
606 36 : FF_Frobenius(GEN x, long e)
607 : {
608 : ulong pp;
609 36 : GEN r, T, p, z=_initFF(x,&T,&p,&pp);
610 36 : ulong n = umodsu(e, FF_f(x));
611 36 : pari_sp av = avma;
612 36 : if (n==0) return gcopy(x);
613 36 : switch(x[1])
614 : {
615 12 : case t_FF_FpXQ:
616 12 : r=FpX_Frobenius(T,p);
617 12 : if (n>1) r=FpXQ_autpow(r,n,T,p);
618 12 : break;
619 12 : case t_FF_F2xq:
620 12 : r=F2x_Frobenius(T);
621 12 : if (n>1) r=F2xq_autpow(r,n,T);
622 12 : break;
623 12 : default:
624 12 : r=Flx_Frobenius(T,pp);
625 12 : if (n>1) r=Flxq_autpow(r,n,T,pp);
626 : }
627 36 : r = gc_upto(av, r);
628 36 : return _mkFF(x,z,r);
629 : }
630 :
631 : static GEN
632 882 : FFX_preimage_i(GEN x, GEN y, GEN F, GEN T, GEN p, long pp)
633 : {
634 : GEN r;
635 882 : F = FFX_to_raw(F, y);
636 882 : switch(y[1])
637 : {
638 324 : case t_FF_FpXQ:
639 324 : r = FpXQX_rem(gel(x,2), F, T, p);
640 324 : break;
641 324 : case t_FF_F2xq:
642 324 : r = F2xqX_rem(F2x_to_F2xX(gel(x,2),T[1]), F, T);
643 324 : break;
644 234 : default:
645 234 : r = FlxqX_rem(Flx_to_FlxX(gel(x,2),T[1]), F, T, pp);
646 : }
647 882 : return r;
648 : }
649 :
650 : GEN
651 828 : FFX_preimage(GEN x, GEN F, GEN y)
652 : {
653 : GEN r, T, p, z;
654 : ulong pp;
655 828 : if (FF_equal0(x)) return FF_zero(y);
656 750 : z = _initFF(y,&T,&p,&pp);
657 750 : r = FFX_preimage_i(x, y, F, T, p, pp);
658 750 : if (degpol(r) > 0) return NULL;
659 678 : r = (y[1] == t_FF_FpXQ)? Fq_to_FpXQ(gel(r,2),T, p): gel(r,2);
660 678 : return _mkFF(y,z,r);
661 : }
662 :
663 : GEN
664 144 : FFX_preimagerel(GEN x, GEN F, GEN y)
665 : {
666 144 : pari_sp av = avma;
667 : GEN r, T, p;
668 : ulong pp;
669 144 : if (FF_equal0(x)) return FF_zero(y);
670 132 : _getFF(y,&T,&p,&pp);
671 132 : r = FFX_preimage_i(x, y, F, T, p, pp);
672 132 : return gc_GEN(av, raw_to_FFX(r, y));
673 : }
674 :
675 : GEN
676 11726951 : FF_mul(GEN x, GEN y)
677 : {
678 : ulong pp;
679 11726951 : GEN r, T, p, z=_initFF(x,&T,&p,&pp);
680 11726951 : pari_sp av=avma;
681 11726951 : _checkFF(x,y,"*");
682 11726951 : switch(x[1])
683 : {
684 498037 : case t_FF_FpXQ:
685 498037 : r=FpXQ_mul(gel(x,2),gel(y,2),T,p);
686 498037 : break;
687 694519 : case t_FF_F2xq:
688 694519 : r=F2xq_mul(gel(x,2),gel(y,2),T);
689 694519 : break;
690 10534395 : default:
691 10534395 : r=Flxq_mul(gel(x,2),gel(y,2),T,pp);
692 : }
693 11726951 : return _mkFF(x,z,gc_upto(av, r));
694 : }
695 :
696 : GEN
697 1843381 : FF_Z_mul(GEN x, GEN y)
698 : {
699 : ulong pp;
700 1843381 : GEN r, T, p, A = gel(x,2), z=_initFF(x,&T,&p,&pp);
701 1843381 : switch(x[1])
702 : {
703 49757 : case t_FF_FpXQ: /* modii(y,p) left on stack for efficiency */
704 49757 : r = FpX_Fp_mul(A, modii(y,p),p);
705 49757 : break;
706 951337 : case t_FF_F2xq:
707 951337 : r = mpodd(y)? vecsmall_copy(A): zero_Flx(A[1]);
708 951337 : break;
709 842287 : default:
710 842287 : r = Flx_Fl_mul(A, umodiu(y,pp), pp);
711 : }
712 1843381 : return _mkFF(x,z,r);
713 : }
714 :
715 : GEN
716 2772 : FF_Z_Z_muldiv(GEN x, GEN a, GEN b)
717 : {
718 : ulong pp;
719 2772 : GEN r, T, p, A = gel(x,2), z=_initFF(x,&T,&p,&pp);
720 2772 : switch(x[1])
721 : {
722 80 : case t_FF_FpXQ: /* Fp_div(a,b,p) left on stack for efficiency */
723 80 : r = FpX_Fp_mul(A, Fp_div(a,b,p), p);
724 80 : break;
725 126 : case t_FF_F2xq:
726 126 : if (!mpodd(b)) pari_err_INV("FF_Z_Z_muldiv", b);
727 120 : r = mpodd(a)? vecsmall_copy(A): zero_Flx(A[1]);
728 120 : break;
729 2566 : default:
730 2566 : r = Flx_Fl_mul(A, Fl_div(umodiu(a,pp),umodiu(b,pp),pp),pp);
731 : }
732 2760 : return _mkFF(x,z,r);
733 : }
734 :
735 : GEN
736 2543251 : FF_sqr(GEN x)
737 : {
738 : ulong pp;
739 2543251 : GEN r, T, p, z=_initFF(x,&T,&p,&pp);
740 2543251 : switch(x[1])
741 : {
742 16436 : case t_FF_FpXQ:
743 : {
744 16436 : pari_sp av=avma;
745 16436 : r=gc_upto(av,FpXQ_sqr(gel(x,2),T,p));
746 16436 : break;
747 : }
748 403269 : case t_FF_F2xq:
749 403269 : r=F2xq_sqr(gel(x,2),T);
750 403269 : break;
751 2123546 : default:
752 2123546 : r=Flxq_sqr(gel(x,2),T,pp);
753 : }
754 2543251 : return _mkFF(x,z,r);
755 : }
756 :
757 : GEN
758 185943 : FF_mul2n(GEN x, long n)
759 : {
760 : ulong pp;
761 185943 : GEN r, T, p, A = gel(x,2), z=_initFF(x,&T,&p,&pp);
762 185943 : switch(x[1])
763 : {
764 2552 : case t_FF_FpXQ:
765 : {
766 : GEN p1; /* left on stack for efficiency */
767 2552 : if (n>0) p1=remii(int2n(n),p);
768 31 : else p1=Fp_inv(remii(int2n(-n),p),p);
769 2552 : r = FpX_Fp_mul(A, p1, p);
770 : }
771 2552 : break;
772 77769 : case t_FF_F2xq:
773 77769 : if (n<0) pari_err_INV("FF_mul2n", gen_2);
774 77769 : r = n==0? vecsmall_copy(A): zero_Flx(A[1]);
775 77769 : break;
776 105622 : default:
777 : {
778 : ulong l1;
779 105622 : if (n>0) l1 = umodiu(int2n(n),pp);
780 221 : else l1 = Fl_inv(umodiu(int2n(-n),pp),pp);
781 105622 : r = Flx_Fl_mul(A,l1,pp);
782 : }
783 : }
784 185943 : return _mkFF(x,z,r);
785 : }
786 :
787 : GEN
788 10905 : FF_inv(GEN x)
789 : {
790 : ulong pp;
791 10905 : GEN r, T, p, z=_initFF(x,&T,&p,&pp);
792 10905 : pari_sp av=avma;
793 10905 : switch(x[1])
794 : {
795 2584 : case t_FF_FpXQ:
796 2584 : r=gc_upto(av,FpXQ_inv(gel(x,2),T,p));
797 2584 : break;
798 6793 : case t_FF_F2xq:
799 6793 : r=F2xq_inv(gel(x,2),T);
800 6793 : break;
801 1528 : default:
802 1528 : r=Flxq_inv(gel(x,2),T,pp);
803 : }
804 10905 : return _mkFF(x,z,r);
805 : }
806 :
807 : GEN
808 210211 : FF_div(GEN x, GEN y)
809 : {
810 : ulong pp;
811 210211 : GEN r, T, p, z=_initFF(x,&T,&p,&pp);
812 210211 : pari_sp av=avma;
813 210211 : _checkFF(x,y,"/");
814 210211 : switch(x[1])
815 : {
816 2953 : case t_FF_FpXQ:
817 2953 : r=gc_upto(av,FpXQ_div(gel(x,2),gel(y,2),T,p));
818 2953 : break;
819 88546 : case t_FF_F2xq:
820 88546 : r=gc_upto(av,F2xq_div(gel(x,2),gel(y,2),T));
821 88486 : break;
822 118712 : default:
823 118712 : r=gc_upto(av,Flxq_div(gel(x,2),gel(y,2),T,pp));
824 : }
825 210127 : return _mkFF(x,z,r);
826 : }
827 :
828 : GEN
829 390 : Z_FF_div(GEN n, GEN x)
830 : {
831 : ulong pp;
832 390 : GEN r, T, p, A = gel(x,2), z=_initFF(x,&T,&p,&pp);
833 390 : pari_sp av=avma;
834 390 : switch(x[1])
835 : {
836 14 : case t_FF_FpXQ:
837 14 : r = gc_upto(av,FpX_Fp_mul(FpXQ_inv(A,T,p),modii(n,p),p));
838 14 : break;
839 48 : case t_FF_F2xq:
840 48 : r = F2xq_inv(A,T); /*Check for division by 0*/
841 48 : if(!mpodd(n)) { set_avma(av); r = zero_Flx(A[1]); }
842 48 : break;
843 328 : default:
844 328 : r = gc_upto(av, Flx_Fl_mul(Flxq_inv(A,T,pp),umodiu(n,pp),pp));
845 : }
846 390 : return _mkFF(x,z,r);
847 : }
848 :
849 : GEN
850 1871 : FF_sqrtn(GEN x, GEN n, GEN *zetan)
851 : {
852 : ulong pp;
853 1871 : GEN r, T, p, y=_initFF(x,&T,&p,&pp);
854 1871 : switch (x[1])
855 : {
856 192 : case t_FF_FpXQ:
857 192 : r=FpXQ_sqrtn(gel(x,2),n,T,p,zetan);
858 187 : break;
859 23 : case t_FF_F2xq:
860 23 : r=F2xq_sqrtn(gel(x,2),n,T,zetan);
861 18 : break;
862 1656 : default:
863 1656 : r=Flxq_sqrtn(gel(x,2),n,T,pp,zetan);
864 : }
865 1856 : if (!r) pari_err_SQRTN("FF_sqrtn",x);
866 1856 : (void)_mkFF(x, y, r);
867 1856 : if (zetan)
868 : {
869 140 : GEN z = cgetg(lg(y),t_FFELT);
870 140 : *zetan=_mkFF(x, z, *zetan);
871 : }
872 1856 : return y;
873 : }
874 :
875 : GEN
876 6138 : FF_sqrt(GEN x)
877 : {
878 : ulong pp;
879 6138 : GEN r, T, p, y=_initFF(x,&T,&p,&pp);
880 6138 : switch (x[1])
881 : {
882 22 : case t_FF_FpXQ:
883 22 : r = FpXQ_sqrt(gel(x,2),T,p);
884 22 : break;
885 48 : case t_FF_F2xq:
886 48 : r = F2xq_sqrt(gel(x,2),T);
887 48 : break;
888 6068 : default:
889 6068 : r = Flxq_sqrt(gel(x,2),T,pp);
890 : }
891 6138 : if (!r) pari_err_SQRTN("FF_sqrt",x);
892 138 : return _mkFF(x, y, r);
893 : }
894 :
895 : long
896 6 : FF_issquare(GEN x)
897 : {
898 : GEN T, p;
899 : ulong pp;
900 6 : _getFF(x, &T, &p, &pp);
901 6 : switch(x[1])
902 : {
903 0 : case t_FF_FpXQ:
904 0 : return FpXQ_issquare(gel(x,2), T, p);
905 6 : case t_FF_F2xq:
906 6 : return 1;
907 0 : default: /* case t_FF_Flxq: */
908 0 : return Flxq_issquare(gel(x,2), T, pp);
909 : }
910 : }
911 :
912 : long
913 162 : FF_issquareall(GEN x, GEN *pt)
914 : {
915 162 : if (!pt) return FF_issquare(x);
916 156 : return FF_ispower(x, gen_2, pt);
917 : }
918 :
919 : long
920 180 : FF_ispower(GEN x, GEN K, GEN *pt)
921 : {
922 : ulong pp;
923 : GEN z, r, T, p;
924 180 : pari_sp av = avma;
925 :
926 180 : if (FF_equal0(x)) { if (pt) *pt = gcopy(x); return 1; }
927 180 : _getFF(x, &T, &p, &pp);
928 180 : z = pt? cgetg(5,t_FFELT): NULL;
929 180 : switch(x[1])
930 : {
931 61 : case t_FF_FpXQ:
932 61 : r = FpXQ_sqrtn(gel(x,2),K,T,p,NULL);
933 61 : break;
934 36 : case t_FF_F2xq:
935 36 : r = F2xq_sqrtn(gel(x,2),K,T,NULL);
936 36 : break;
937 83 : default: /* case t_FF_Flxq: */
938 83 : r = Flxq_sqrtn(gel(x,2),K,T,pp,NULL);
939 83 : break;
940 : }
941 180 : if (!r) return gc_long(av,0);
942 108 : if (pt) { *pt = z; (void)_mkFF(x,z,r); }
943 108 : return 1;
944 : }
945 :
946 : GEN
947 109 : FF_pow(GEN x, GEN n)
948 : {
949 : ulong pp;
950 109 : GEN r, T, p, z=_initFF(x,&T,&p,&pp);
951 109 : switch(x[1])
952 : {
953 22 : case t_FF_FpXQ:
954 22 : r = FpXQ_pow(gel(x,2), n, T, p);
955 22 : break;
956 3 : case t_FF_F2xq:
957 3 : r = F2xq_pow(gel(x,2), n, T);
958 3 : break;
959 84 : default:
960 84 : r = Flxq_pow(gel(x,2), n, T, pp);
961 : }
962 109 : return _mkFF(x,z,r);
963 : }
964 :
965 : GEN
966 24 : FF_norm(GEN x)
967 : {
968 : ulong pp;
969 : GEN T,p;
970 24 : _getFF(x,&T,&p,&pp);
971 24 : switch (x[1])
972 : {
973 1 : case t_FF_FpXQ:
974 1 : return FpXQ_norm(gel(x,2),T,p);
975 6 : case t_FF_F2xq:
976 6 : return lgpol(gel(x,2))?gen_1:gen_0;
977 17 : default:
978 17 : return utoi(Flxq_norm(gel(x,2),T,pp));
979 : }
980 : }
981 :
982 : GEN
983 24 : FF_trace(GEN x)
984 : {
985 : ulong pp;
986 : GEN T,p;
987 24 : _getFF(x,&T,&p,&pp);
988 24 : switch(x[1])
989 : {
990 1 : case t_FF_FpXQ:
991 1 : return FpXQ_trace(gel(x,2),T,p);
992 6 : case t_FF_F2xq:
993 6 : return F2xq_trace(gel(x,2),T)?gen_1:gen_0;
994 17 : default:
995 17 : return utoi(Flxq_trace(gel(x,2),T,pp));
996 : }
997 : }
998 :
999 : GEN
1000 24 : FF_conjvec(GEN x)
1001 : {
1002 : ulong pp;
1003 : GEN r,T,p,v;
1004 : long i,l;
1005 : pari_sp av;
1006 24 : _getFF(x,&T,&p,&pp);
1007 24 : av = avma;
1008 24 : switch(x[1])
1009 : {
1010 1 : case t_FF_FpXQ:
1011 1 : v = FpXQ_conjvec(gel(x,2), T, p);
1012 1 : break;
1013 6 : case t_FF_F2xq:
1014 6 : v = F2xq_conjvec(gel(x,2), T);
1015 6 : break;
1016 17 : default:
1017 17 : v = Flxq_conjvec(gel(x,2), T, pp);
1018 : }
1019 24 : l = lg(v); r = cgetg(l, t_COL);
1020 222 : for(i=1; i<l; i++)
1021 198 : gel(r,i) = mkFF_i(x, gel(v,i));
1022 24 : return gc_GEN(av, r);
1023 : }
1024 :
1025 : GEN
1026 24 : FF_charpoly(GEN x)
1027 : {
1028 : ulong pp;
1029 : GEN T,p;
1030 24 : pari_sp av=avma;
1031 24 : _getFF(x,&T,&p,&pp);
1032 24 : switch(x[1])
1033 : {
1034 1 : case t_FF_FpXQ:
1035 1 : return gc_upto(av,FpXQ_charpoly(gel(x,2), T, p));
1036 6 : case t_FF_F2xq:
1037 6 : return gc_upto(av,Flx_to_ZX(Flxq_charpoly(F2x_to_Flx(gel(x,2)),
1038 : F2x_to_Flx(T), 2UL)));
1039 17 : default:
1040 17 : return gc_upto(av,Flx_to_ZX(Flxq_charpoly(gel(x,2), T, pp)));
1041 : }
1042 : }
1043 :
1044 : GEN
1045 132 : FF_minpoly(GEN x)
1046 : {
1047 : ulong pp;
1048 : GEN T,p;
1049 132 : pari_sp av=avma;
1050 132 : _getFF(x,&T,&p,&pp);
1051 132 : switch(x[1])
1052 : {
1053 37 : case t_FF_FpXQ:
1054 37 : return gc_upto(av,FpXQ_minpoly(gel(x,2), T, p));
1055 42 : case t_FF_F2xq:
1056 42 : return gc_upto(av,Flx_to_ZX(Flxq_minpoly(F2x_to_Flx(gel(x,2)),
1057 : F2x_to_Flx(T), 2UL)));
1058 53 : default:
1059 53 : return gc_upto(av,Flx_to_ZX(Flxq_minpoly(gel(x,2), T, pp)));
1060 : }
1061 : }
1062 :
1063 : GEN
1064 162 : FF_log(GEN x, GEN g, GEN ord)
1065 : {
1066 162 : pari_sp av=avma;
1067 : ulong pp;
1068 : GEN r, T, p;
1069 162 : _getFF(x,&T,&p,&pp);
1070 162 : _checkFF(x,g,"log");
1071 162 : switch(x[1])
1072 : {
1073 38 : case t_FF_FpXQ:
1074 38 : if (!ord) ord = factor_pn_1(p,degpol(T));
1075 38 : r = FpXQ_log(gel(x,2), gel(g,2), ord, T, p);
1076 14 : break;
1077 24 : case t_FF_F2xq:
1078 24 : if (!ord) ord = factor_pn_1(gen_2,F2x_degree(T));
1079 24 : r = F2xq_log(gel(x,2), gel(g,2), ord, T);
1080 24 : break;
1081 100 : default:
1082 100 : if (!ord) ord = factor_pn_1(p,degpol(T));
1083 100 : r = Flxq_log(gel(x,2), gel(g,2), ord, T, pp);
1084 : }
1085 138 : return gc_upto(av, r);
1086 : }
1087 :
1088 : GEN
1089 48 : FF_order(GEN x, GEN o)
1090 : {
1091 48 : pari_sp av=avma;
1092 : ulong pp;
1093 : GEN r, T,p;
1094 48 : _getFF(x,&T,&p,&pp);
1095 48 : switch(x[1])
1096 : {
1097 2 : case t_FF_FpXQ:
1098 2 : if (!o) o = factor_pn_1(p,degpol(T));
1099 2 : r = FpXQ_order(gel(x,2), o, T, p);
1100 2 : break;
1101 12 : case t_FF_F2xq:
1102 12 : if (!o) o = factor_pn_1(gen_2,F2x_degree(T));
1103 12 : r = F2xq_order(gel(x,2), o, T);
1104 12 : break;
1105 34 : default:
1106 34 : if (!o) o = factor_pn_1(p,degpol(T));
1107 34 : r = Flxq_order(gel(x,2), o, T, pp);
1108 : }
1109 48 : if (!o) r = gc_INT(av,r);
1110 48 : return r;
1111 : }
1112 :
1113 : GEN
1114 348 : FF_primroot(GEN x, GEN *o)
1115 : {
1116 : ulong pp;
1117 348 : GEN r,T,p,z=_initFF(x,&T,&p,&pp);
1118 348 : switch(x[1])
1119 : {
1120 26 : case t_FF_FpXQ:
1121 26 : r = gener_FpXQ(T, p, o);
1122 26 : break;
1123 138 : case t_FF_F2xq:
1124 138 : r = gener_F2xq(T, o);
1125 138 : break;
1126 184 : default:
1127 184 : r = gener_Flxq(T, pp, o);
1128 : }
1129 348 : return _mkFF(x,z,r);
1130 : }
1131 :
1132 : static GEN
1133 441693 : to_FFE(GEN P, GEN fg)
1134 : {
1135 441693 : if(ell_is_inf(P))
1136 200106 : return ellinf();
1137 : else
1138 241587 : retmkvec2(mkFF_i(fg,gel(P,1)), mkFF_i(fg,gel(P,2)));
1139 : }
1140 :
1141 : static GEN
1142 15528 : to_FFE_vec(GEN x, GEN ff)
1143 : {
1144 15528 : long i, lx = lg(x);
1145 33222 : for (i=1; i<lx; i++) gel(x,i) = to_FFE(gel(x,i), ff);
1146 15528 : return x;
1147 : }
1148 :
1149 : static GEN
1150 1479 : FpXQ_ell_to_a4a6(GEN E, GEN T, GEN p)
1151 : {
1152 : GEN a1, a3, b2, c4, c6;
1153 1479 : a1 = Rg_to_FpXQ(ell_get_a1(E),T,p);
1154 1479 : a3 = Rg_to_FpXQ(ell_get_a3(E),T,p);
1155 1479 : b2 = Rg_to_FpXQ(ell_get_b2(E),T,p);
1156 1479 : c4 = Rg_to_FpXQ(ell_get_c4(E),T,p);
1157 1479 : c6 = Rg_to_FpXQ(ell_get_c6(E),T,p);
1158 1479 : retmkvec3(FpX_neg(FpX_mulu(c4, 27, p), p), FpX_neg(FpX_mulu(c6, 54, p), p),
1159 : mkvec4(Z_to_FpX(utoi(6),p,varn(T)),FpX_mulu(b2,3,p),
1160 : FpX_mulu(a1,3,p),FpX_mulu(a3,108,p)));
1161 : }
1162 :
1163 : static GEN
1164 26964 : F3xq_ell_to_a4a6(GEN E, GEN T)
1165 : {
1166 : GEN a1, a3, b2, b4, b6;
1167 26964 : a1 = Rg_to_Flxq(ell_get_a1(E),T,3);
1168 26964 : a3 = Rg_to_Flxq(ell_get_a3(E),T,3);
1169 26964 : b2 = Rg_to_Flxq(ell_get_b2(E),T,3);
1170 26964 : b4 = Rg_to_Flxq(ell_get_b4(E),T,3);
1171 26964 : b6 = Rg_to_Flxq(ell_get_b6(E),T,3);
1172 26964 : if(lgpol(b2)) /* ordinary case */
1173 : {
1174 18714 : GEN b4b2 = Flxq_div(b4,b2,T,3);
1175 18714 : GEN a6 = Flx_sub(b6,Flxq_mul(b4b2,Flx_add(b4,Flxq_sqr(b4b2,T,3),3),T,3),3);
1176 18714 : retmkvec3(mkvec(b2), a6,
1177 : mkvec4(Fl_to_Flx(1,T[1]),b4b2,Flx_neg(a1,3),Flx_neg(a3,3)));
1178 : }
1179 : else /* super-singular case */
1180 8250 : retmkvec3(Flx_neg(b4, 3), b6,
1181 : mkvec4(Fl_to_Flx(1,T[1]),zero_Flx(T[1]), Flx_neg(a1,3), Flx_neg(a3,3)));
1182 : }
1183 :
1184 : static GEN
1185 63229 : Flxq_ell_to_a4a6(GEN E, GEN T, ulong p)
1186 : {
1187 : GEN a1, a3, b2, c4, c6;
1188 63229 : if(p==3) return F3xq_ell_to_a4a6(E, T);
1189 36265 : a1 = Rg_to_Flxq(ell_get_a1(E),T,p);
1190 36265 : a3 = Rg_to_Flxq(ell_get_a3(E),T,p);
1191 36265 : b2 = Rg_to_Flxq(ell_get_b2(E),T,p);
1192 36265 : c4 = Rg_to_Flxq(ell_get_c4(E),T,p);
1193 36265 : c6 = Rg_to_Flxq(ell_get_c6(E),T,p);
1194 36265 : retmkvec3(Flx_neg(Flx_mulu(c4, 27, p), p), Flx_neg(Flx_mulu(c6, 54, p), p),
1195 : mkvec4(Fl_to_Flx(6%p,T[1]), Flx_triple(b2,p), Flx_triple(a1,p),
1196 : Flx_mulu(a3,108,p)));
1197 : }
1198 :
1199 : static GEN
1200 38777 : F2xq_ell_to_a4a6(GEN E, GEN T)
1201 : {
1202 38777 : long v = T[1];
1203 38777 : GEN a1 = Rg_to_F2xq(ell_get_a1(E),T);
1204 38777 : GEN a2 = Rg_to_F2xq(ell_get_a2(E),T);
1205 38777 : GEN a3 = Rg_to_F2xq(ell_get_a3(E),T);
1206 38777 : GEN a4 = Rg_to_F2xq(ell_get_a4(E),T);
1207 38777 : GEN a6 = Rg_to_F2xq(ell_get_a6(E),T);
1208 38777 : if (lgpol(a1))
1209 : {
1210 6287 : GEN a1i = F2xq_inv(a1,T);
1211 6287 : GEN a1i2 = F2xq_sqr(a1i,T);
1212 6287 : GEN a1i3 = F2xq_mul(a1i,a1i2,T);
1213 6287 : GEN a1i6 = F2xq_sqr(a1i3,T);
1214 6287 : GEN d = F2xq_mul(a3,a1i,T);
1215 6287 : GEN dd = F2xq_mul(d,a1i2,T);
1216 6287 : GEN e = F2xq_mul(F2x_add(a4,F2xq_sqr(d,T)),a1i,T);
1217 6287 : GEN ee = F2xq_mul(e,a1i3,T);
1218 6287 : GEN da2 = F2x_add(a2,d);
1219 6287 : GEN d2 = F2xq_mul(da2,a1i2,T);
1220 6287 : GEN d6 = F2xq_mul(F2x_add(F2x_add(F2xq_mul(F2x_add(F2xq_mul(da2,d,T),a4),d,T),a6),F2xq_sqr(e,T)),a1i6,T);
1221 6287 : retmkvec3(d2, d6, mkvec4(a1i,dd,pol0_F2x(v),ee));
1222 : }
1223 : else
1224 : { /* allow a1 = a3 = 0: singular curve */
1225 32490 : GEN d4 = F2x_add(F2xq_sqr(a2,T),a4);
1226 32490 : GEN d6 = F2x_add(F2xq_mul(a2,a4,T),a6);
1227 32490 : GEN a3i = lgpol(a3)? F2xq_inv(a3,T): a3;
1228 32490 : retmkvec3(mkvec3(a3,d4,a3i), d6, mkvec4(pol1_F2x(v),a2,pol0_F2x(T[1]),pol0_F2x(T[1])));
1229 : }
1230 : }
1231 :
1232 : static GEN
1233 1634 : FqV_to_FpXQV(GEN x, GEN T)
1234 : {
1235 1634 : pari_sp av = avma;
1236 1634 : long v = varn(T), i, s=0, l = lg(x);
1237 1634 : GEN y = shallowcopy(x);
1238 8170 : for(i=1; i<l; i++)
1239 6536 : if (typ(gel(x,i))==t_INT)
1240 : {
1241 0 : gel(y,i) = scalarpol(gel(x,i), v);
1242 0 : s = 1;
1243 : }
1244 1634 : if (!s) { set_avma(av); return x; }
1245 0 : return y;
1246 : }
1247 :
1248 : GEN
1249 88194 : FF_ellcard(GEN E)
1250 : {
1251 : GEN T,p;
1252 : ulong pp;
1253 88194 : GEN e = ellff_get_a4a6(E);
1254 88194 : GEN fg = ellff_get_field(E);
1255 88194 : _getFF(fg,&T,&p,&pp);
1256 88194 : switch(fg[1])
1257 : {
1258 1466 : case t_FF_FpXQ:
1259 1466 : return FpXQ_ellcard(Fq_to_FpXQ(gel(e,1),T,p), Fq_to_FpXQ(gel(e,2),T,p),T,p);
1260 38562 : case t_FF_F2xq:
1261 38562 : return F2xq_ellcard(gel(e,1),gel(e,2),T);
1262 48166 : default:
1263 48166 : return Flxq_ellcard(gel(e,1),gel(e,2),T,pp);
1264 : }
1265 : }
1266 :
1267 : GEN
1268 12 : FF_ellcard_SEA(GEN E, long smallfact)
1269 : {
1270 12 : pari_sp av = avma;
1271 : GEN T,p;
1272 : ulong pp;
1273 12 : GEN e = ellff_get_a4a6(E), a4, a6, card;
1274 12 : GEN fg = ellff_get_field(E);
1275 12 : _getFF(fg,&T,&p,&pp);
1276 12 : switch(fg[1])
1277 : {
1278 0 : case t_FF_FpXQ:
1279 0 : a4 = Fq_to_FpXQ(gel(e,1), T, p);
1280 0 : a6 = Fq_to_FpXQ(gel(e,2), T, p);
1281 0 : card = Fq_ellcard_SEA(a4, a6, powiu(p,degpol(T)), T,p, smallfact);
1282 0 : break;
1283 0 : case t_FF_F2xq:
1284 0 : pari_err_IMPL("SEA for char 2");
1285 12 : default:
1286 12 : a4 = Flx_to_ZX(gel(e,1));
1287 12 : a6 = Flx_to_ZX(gel(e,2));
1288 12 : card = Fq_ellcard_SEA(a4, a6, powuu(pp,degpol(T)), Flx_to_ZX(T), p, smallfact);
1289 : }
1290 12 : return gc_INT(av, card);
1291 : }
1292 :
1293 : /* return G, set m */
1294 : GEN
1295 22386 : FF_ellgroup(GEN E, GEN *pm)
1296 : {
1297 : GEN T, p, e, fg, N;
1298 : ulong pp;
1299 :
1300 22386 : N = ellff_get_card(E);
1301 22386 : e = ellff_get_a4a6(E);
1302 22386 : fg = ellff_get_field(E);
1303 22386 : _getFF(fg,&T,&p,&pp);
1304 22386 : switch(fg[1])
1305 : {
1306 13 : case t_FF_FpXQ:
1307 13 : return FpXQ_ellgroup(Fq_to_FpXQ(gel(e,1),T,p),
1308 13 : Fq_to_FpXQ(gel(e,2),T,p),N,T,p,pm);
1309 3264 : case t_FF_F2xq:
1310 3264 : return F2xq_ellgroup(gel(e,1),gel(e,2),N,T,pm);
1311 19109 : default:
1312 19109 : return Flxq_ellgroup(gel(e,1),gel(e,2),N,T,pp,pm);
1313 : }
1314 : }
1315 :
1316 : GEN
1317 15528 : FF_ellgens(GEN E)
1318 : {
1319 : GEN D, fg, e, m, T, p, F, e3;
1320 : ulong pp;
1321 :
1322 15528 : fg = ellff_get_field(E);
1323 15528 : e = ellff_get_a4a6(E);
1324 15528 : m = ellff_get_m(E);
1325 15528 : D = ellff_get_D(E);
1326 15528 : _getFF(fg,&T,&p,&pp);
1327 15528 : switch(fg[1])
1328 : {
1329 7 : case t_FF_FpXQ:
1330 7 : e3 = FqV_to_FpXQV(gel(e,3),T);
1331 7 : F = FpXQ_ellgens(Fq_to_FpXQ(gel(e,1),T,p),Fq_to_FpXQ(gel(e,2),T,p),e3,D,m,T,p);
1332 7 : break;
1333 3204 : case t_FF_F2xq:
1334 3204 : F = F2xq_ellgens(gel(e,1),gel(e,2),gel(e,3),D,m,T);
1335 3204 : break;
1336 12317 : default:
1337 12317 : F = Flxq_ellgens(gel(e,1),gel(e,2),gel(e,3),D,m,T, pp);
1338 12317 : break;
1339 : }
1340 15528 : return to_FFE_vec(F,fg);
1341 : }
1342 :
1343 : GEN
1344 102 : FF_elldata(GEN E, GEN fg)
1345 : {
1346 : GEN T,p,e;
1347 : ulong pp;
1348 102 : _getFF(fg,&T,&p,&pp);
1349 102 : switch(fg[1])
1350 : {
1351 0 : case t_FF_FpXQ:
1352 0 : e = FpXQ_ell_to_a4a6(E,T,p); break;
1353 48 : case t_FF_F2xq:
1354 48 : e = F2xq_ell_to_a4a6(E,T); break;
1355 54 : default:
1356 54 : e = Flxq_ell_to_a4a6(E,T,pp); break;
1357 : }
1358 102 : return mkvec2(fg,e);
1359 : }
1360 :
1361 : /* allow singular E, set E.j = 0 in this case */
1362 : GEN
1363 103383 : FF_ellinit(GEN E, GEN fg)
1364 : {
1365 103383 : GEN T,p,e, D, c4, y = E;
1366 : ulong pp;
1367 : long i;
1368 103383 : _getFF(fg,&T,&p,&pp);
1369 103383 : switch(fg[1])
1370 : {
1371 1479 : case t_FF_FpXQ:
1372 1479 : e = FpXQ_ell_to_a4a6(E,T,p);
1373 19227 : for(i=1;i<=12;i++) gel(y,i) = mkFF_i(fg,Rg_to_FpXQ(gel(E,i),T,p));
1374 1479 : break;
1375 38729 : case t_FF_F2xq:
1376 38729 : e = F2xq_ell_to_a4a6(E,T);
1377 503477 : for(i=1;i<=12;i++) gel(y,i) = mkFF_i(fg,Rg_to_F2xq(gel(E,i),T));
1378 38729 : break;
1379 63175 : default:
1380 63175 : e = Flxq_ell_to_a4a6(E,T,pp);
1381 821275 : for(i=1;i<=12;i++) gel(y,i) = mkFF_i(fg,Rg_to_Flxq(gel(E,i),T,pp));
1382 63175 : break;
1383 : }
1384 103383 : D = ell_get_disc(y);
1385 103383 : c4 = ell_get_c4(y);
1386 103383 : gel(y,13) = FF_equal0(D)? D: gdiv(gpowgs(c4,3), D);
1387 103383 : gel(y,14) = mkvecsmall(t_ELL_Fq);
1388 103383 : gel(y,15) = mkvec2(fg,e);
1389 103383 : return y;
1390 : }
1391 :
1392 : GEN
1393 23304 : FF_elltwist(GEN E)
1394 : {
1395 23304 : pari_sp av = avma;
1396 23304 : GEN fg = ellff_get_field(E), e = ellff_get_a4a6(E);
1397 : GEN T, p, a, a6, V;
1398 : ulong pp;
1399 23304 : _getFF(fg,&T,&p,&pp);
1400 23304 : switch (fg[1])
1401 : {
1402 720 : case t_FF_FpXQ:
1403 720 : FpXQ_elltwist(gel(e,1), gel(e,2), T, p, &a, &a6);
1404 720 : V = mkvec5(gen_0, gen_0, gen_0, mkFF_i(fg, a), mkFF_i(fg, a6));
1405 720 : break;
1406 3012 : case t_FF_F2xq:
1407 3012 : F2xq_elltwist(gel(e,1), gel(e,2), T, &a, &a6);
1408 6024 : V = typ(a)==t_VECSMALL ?
1409 3012 : mkvec5(gen_1, mkFF_i(fg, a), gen_0, gen_0, mkFF_i(fg, a6)):
1410 0 : mkvec5(gen_0, gen_0, mkFF_i(fg, gel(a,1)), mkFF_i(fg, gel(a,2)), mkFF_i(fg, a6));
1411 3012 : break;
1412 19572 : default:
1413 19572 : Flxq_elltwist(gel(e,1), gel(e,2), T, pp, &a, &a6);
1414 19572 : V = typ(a)==t_VECSMALL ?
1415 19572 : mkvec5(gen_0, gen_0, gen_0, mkFF_i(fg, a), mkFF_i(fg, a6)):
1416 6516 : mkvec5(gen_0, mkFF_i(fg, gel(a,1)), gen_0, gen_0, mkFF_i(fg, a6));
1417 : }
1418 23304 : return gc_GEN(av, V);
1419 : }
1420 :
1421 : static long
1422 0 : F3x_equalm1(GEN x) { return degpol(x)==0 && x[2] == 2; }
1423 : GEN
1424 210555 : FF_ellrandom(GEN E)
1425 : {
1426 210555 : pari_sp av = avma;
1427 210555 : GEN fg = ellff_get_field(E), e = ellff_get_a4a6(E), Q;
1428 : GEN T,p;
1429 : ulong pp;
1430 210555 : _getFF(fg,&T,&p,&pp);
1431 210555 : switch (fg[1])
1432 : {
1433 801 : case t_FF_FpXQ:
1434 801 : Q = random_FpXQE(Fq_to_FpXQ(gel(e,1),T,p), Fq_to_FpXQ(gel(e,2),T,p), T, p);
1435 801 : Q = FpXQE_changepoint(Q, FqV_to_FpXQV(gel(e,3), T) , T, p);
1436 801 : break;
1437 143807 : case t_FF_F2xq:
1438 : {
1439 143807 : long d = F2x_degree(T);
1440 : /* if #E(Fq) = 1 return [0] */
1441 143807 : if (d<=2 && typ(gel(e,1)) == t_VEC)
1442 : { /* over F2 or F4, supersingular */
1443 1182 : GEN v = gel(e,1), A6 = gel(e,2), a3 = gel(v,1), A4 = gel(v,2);
1444 1182 : if (F2x_equal1(a3) &&
1445 120 : ((d==1 && F2x_equal1(A4) && F2x_equal1(A6))
1446 414 : || (d==2 && !lgpol(A4) && F2x_degree(A6)==1))) return ellinf();
1447 : }
1448 143759 : Q = random_F2xqE(gel(e,1), gel(e,2), T);
1449 143759 : Q = F2xqE_changepoint(Q, gel(e,3), T);
1450 143759 : break;
1451 : }
1452 65947 : default:
1453 : /* if #E(Fq) = 1 return [0] */
1454 65947 : if (pp==3 && degpol(T)==1 && typ(gel(e,1))==t_VECSMALL)
1455 : { /* over F3, supersingular */
1456 0 : GEN mb4 = gel(e,1), b6 = gel(e,2);
1457 0 : if (F3x_equalm1(mb4) && F3x_equalm1(b6)) return ellinf();
1458 : }
1459 65947 : Q = random_FlxqE(gel(e,1), gel(e,2), T, pp);
1460 65947 : Q = FlxqE_changepoint(Q, gel(e,3), T, pp);
1461 : }
1462 210507 : return gc_GEN(av, to_FFE(Q, fg));
1463 : }
1464 :
1465 : GEN
1466 213492 : FF_ellmul(GEN E, GEN P, GEN n)
1467 : {
1468 213492 : pari_sp av = avma;
1469 213492 : GEN fg = ellff_get_field(E), e = ellff_get_a4a6(E), Q;
1470 : GEN T,p, Pp, Qp, e3;
1471 : ulong pp;
1472 213492 : _getFF(fg,&T,&p,&pp);
1473 213492 : switch (fg[1])
1474 : {
1475 802 : case t_FF_FpXQ:
1476 802 : e3 = FqV_to_FpXQV(gel(e,3),T);
1477 802 : Pp = FpXQE_changepointinv(RgE_to_FpXQE(P, T, p), e3, T, p);
1478 802 : Qp = FpXQE_mul(Pp, n, gel(e,1), T, p);
1479 802 : Q = FpXQE_changepoint(Qp, gel(e,3), T, p);
1480 802 : break;
1481 157776 : case t_FF_F2xq:
1482 157776 : Pp = F2xqE_changepointinv(RgE_to_F2xqE(P, T), gel(e,3), T);
1483 157776 : Qp = F2xqE_mul(Pp, n, gel(e,1), T);
1484 157776 : Q = F2xqE_changepoint(Qp, gel(e,3), T);
1485 157776 : break;
1486 54914 : default:
1487 54914 : Pp = FlxqE_changepointinv(RgE_to_FlxqE(P, T, pp), gel(e,3), T, pp);
1488 54914 : Qp = FlxqE_mul(Pp, n, gel(e,1), T, pp);
1489 54914 : Q = FlxqE_changepoint(Qp, gel(e,3), T, pp);
1490 : }
1491 213492 : return gc_GEN(av, to_FFE(Q, fg));
1492 : }
1493 :
1494 : GEN
1495 1512 : FF_ellorder(GEN E, GEN P, GEN o)
1496 : {
1497 1512 : pari_sp av = avma;
1498 1512 : GEN fg = ellff_get_field(E), e = ellff_get_a4a6(E);
1499 : GEN r,T,p,Pp,e3;
1500 : ulong pp;
1501 1512 : _getFF(fg,&T,&p,&pp);
1502 1512 : switch (fg[1])
1503 : {
1504 12 : case t_FF_FpXQ:
1505 12 : e3 = FqV_to_FpXQV(gel(e,3),T);
1506 12 : Pp = FpXQE_changepointinv(RgE_to_FpXQE(P,T,p), e3, T, p);
1507 12 : r = FpXQE_order(Pp, o, gel(e,1), T, p);
1508 12 : break;
1509 240 : case t_FF_F2xq:
1510 240 : Pp = F2xqE_changepointinv(RgE_to_F2xqE(P,T), gel(e,3), T);
1511 240 : r = F2xqE_order(Pp, o, gel(e,1), T);
1512 240 : break;
1513 1260 : default:
1514 1260 : Pp = FlxqE_changepointinv(RgE_to_FlxqE(P,T,pp), gel(e,3), T, pp);
1515 1260 : r = FlxqE_order(Pp, o, gel(e,1), T, pp);
1516 : }
1517 1512 : return gc_upto(av, r);
1518 : }
1519 :
1520 : GEN
1521 78 : FF_elllog(GEN E, GEN P, GEN Q, GEN o)
1522 : {
1523 78 : pari_sp av = avma;
1524 78 : GEN fg = ellff_get_field(E), e = ellff_get_a4a6(E);
1525 : GEN r,T,p, Pp,Qp, e3;
1526 : ulong pp;
1527 78 : _getFF(fg,&T,&p,&pp);
1528 78 : switch(fg[1])
1529 : {
1530 0 : case t_FF_FpXQ:
1531 0 : e3 = FqV_to_FpXQV(gel(e,3),T);
1532 0 : Pp = FpXQE_changepointinv(RgE_to_FpXQE(P,T,p), e3, T, p);
1533 0 : Qp = FpXQE_changepointinv(RgE_to_FpXQE(Q,T,p), e3, T, p);
1534 0 : r = FpXQE_log(Pp, Qp, o, gel(e,1), T, p);
1535 0 : break;
1536 36 : case t_FF_F2xq:
1537 36 : Pp = F2xqE_changepointinv(RgE_to_F2xqE(P,T), gel(e,3), T);
1538 36 : Qp = F2xqE_changepointinv(RgE_to_F2xqE(Q,T), gel(e,3), T);
1539 36 : r = F2xqE_log(Pp, Qp, o, gel(e,1), T);
1540 36 : break;
1541 42 : default:
1542 42 : Pp = FlxqE_changepointinv(RgE_to_FlxqE(P,T,pp), gel(e,3), T, pp);
1543 42 : Qp = FlxqE_changepointinv(RgE_to_FlxqE(Q,T,pp), gel(e,3), T, pp);
1544 42 : r = FlxqE_log(Pp, Qp, o, gel(e,1), T, pp);
1545 : }
1546 78 : return gc_upto(av, r);
1547 : }
1548 :
1549 : GEN
1550 4284 : FF_ellweilpairing(GEN E, GEN P, GEN Q, GEN m)
1551 : {
1552 4284 : GEN fg = ellff_get_field(E), e = ellff_get_a4a6(E);
1553 : GEN r,T,p, Pp,Qp, e3;
1554 : ulong pp;
1555 4284 : GEN z=_initFF(fg,&T,&p,&pp);
1556 4284 : pari_sp av = avma;
1557 4284 : switch(fg[1])
1558 : {
1559 6 : case t_FF_FpXQ:
1560 6 : e3 = FqV_to_FpXQV(gel(e,3),T);
1561 6 : Pp = FpXQE_changepointinv(RgE_to_FpXQE(P,T,p), e3, T, p);
1562 6 : Qp = FpXQE_changepointinv(RgE_to_FpXQE(Q,T,p), e3, T, p);
1563 6 : r = FpXQE_weilpairing(Pp, Qp, m, gel(e,1), T, p);
1564 6 : break;
1565 4260 : case t_FF_F2xq:
1566 4260 : Pp = F2xqE_changepointinv(RgE_to_F2xqE(P,T), gel(e,3), T);
1567 4260 : Qp = F2xqE_changepointinv(RgE_to_F2xqE(Q,T), gel(e,3), T);
1568 4260 : r = F2xqE_weilpairing(Pp, Qp, m, gel(e,1), T);
1569 4260 : break;
1570 18 : default:
1571 18 : Pp = FlxqE_changepointinv(RgE_to_FlxqE(P,T,pp), gel(e,3), T, pp);
1572 18 : Qp = FlxqE_changepointinv(RgE_to_FlxqE(Q,T,pp), gel(e,3), T, pp);
1573 18 : r = FlxqE_weilpairing(Pp, Qp, m, gel(e,1), T, pp);
1574 : }
1575 4284 : r = gc_upto(av, r);
1576 4284 : return _mkFF(fg,z,r);
1577 : }
1578 :
1579 : GEN
1580 84 : FF_elltatepairing(GEN E, GEN P, GEN Q, GEN m)
1581 : {
1582 84 : GEN fg = ellff_get_field(E), e = ellff_get_a4a6(E);
1583 : GEN r,T,p, Pp,Qp, e3;
1584 : ulong pp;
1585 84 : GEN z=_initFF(fg,&T,&p,&pp);
1586 84 : pari_sp av = avma;
1587 84 : switch(fg[1])
1588 : {
1589 6 : case t_FF_FpXQ:
1590 6 : e3 = FqV_to_FpXQV(gel(e,3),T);
1591 6 : Pp = FpXQE_changepointinv(RgE_to_FpXQE(P,T,p), e3, T, p);
1592 6 : Qp = FpXQE_changepointinv(RgE_to_FpXQE(Q,T,p), e3, T, p);
1593 6 : r = FpXQE_tatepairing(Pp, Qp, m, gel(e,1), T, p);
1594 6 : break;
1595 24 : case t_FF_F2xq:
1596 24 : Pp = F2xqE_changepointinv(RgE_to_F2xqE(P,T), gel(e,3), T);
1597 24 : Qp = F2xqE_changepointinv(RgE_to_F2xqE(Q,T), gel(e,3), T);
1598 24 : r = F2xqE_tatepairing(Pp, Qp, m, gel(e,1), T);
1599 24 : break;
1600 54 : default:
1601 54 : Pp = FlxqE_changepointinv(RgE_to_FlxqE(P,T,pp), gel(e,3), T, pp);
1602 54 : Qp = FlxqE_changepointinv(RgE_to_FlxqE(Q,T,pp), gel(e,3), T, pp);
1603 54 : r = FlxqE_tatepairing(Pp, Qp, m, gel(e,1), T, pp);
1604 : }
1605 84 : r = gc_upto(av, r);
1606 84 : return _mkFF(fg,z,r);
1607 : }
1608 :
1609 : GEN
1610 88860 : FFX_roots(GEN Pf, GEN ff)
1611 : {
1612 88860 : pari_sp av = avma;
1613 : GEN r,T,p;
1614 : ulong pp;
1615 88860 : GEN P = FFX_to_raw(Pf, ff);
1616 88860 : _getFF(ff,&T,&p,&pp);
1617 88860 : switch(ff[1])
1618 : {
1619 73 : case t_FF_FpXQ:
1620 73 : r = FpXQX_roots(P, T, p);
1621 73 : break;
1622 49296 : case t_FF_F2xq:
1623 49296 : r = F2xqX_roots(P, T);
1624 49296 : break;
1625 39491 : default:
1626 39491 : r = FlxqX_roots(P, T, pp);
1627 : }
1628 88860 : return gc_GEN(av, raw_to_FFC(r, ff));
1629 : }
1630 :
1631 : long
1632 0 : FFX_nbroots(GEN Pf, GEN ff)
1633 : {
1634 0 : pari_sp av = avma;
1635 : GEN T,p;
1636 : ulong pp;
1637 : long r;
1638 0 : GEN P = FFX_to_raw(Pf, ff);
1639 0 : _getFF(ff,&T,&p,&pp);
1640 0 : switch(ff[1])
1641 : {
1642 0 : case t_FF_FpXQ:
1643 0 : r = FpXQX_nbroots(P, T, p);
1644 0 : break;
1645 0 : case t_FF_F2xq:
1646 0 : r = F2xqX_nbroots(P, T);
1647 0 : break;
1648 0 : default:
1649 0 : r = FlxqX_nbroots(P, T, pp);
1650 : }
1651 0 : return gc_long(av, r);
1652 : }
1653 :
1654 : static GEN
1655 120 : raw_to_FFXC(GEN x, GEN ff) { pari_APPLY_type(t_COL, raw_to_FFX(gel(x,i), ff)); }
1656 : static GEN
1657 54 : raw_to_FFXM(GEN x, GEN ff) { pari_APPLY_same(raw_to_FFXC(gel(x,i), ff)); }
1658 : static GEN
1659 384 : raw_to_FFX_fact(GEN F, GEN ff)
1660 : {
1661 : GEN y, u, v;
1662 384 : GEN P = gel(F,1), E = gel(F,2);
1663 384 : long j, l = lg(P);
1664 384 : y = cgetg(3,t_MAT);
1665 384 : u = cgetg(l,t_COL); gel(y,1) = u;
1666 384 : v = cgetg(l,t_COL); gel(y,2) = v;
1667 1794 : for (j=1; j<l; j++)
1668 : {
1669 1410 : gel(u,j) = raw_to_FFX(gel(P,j), ff);
1670 1410 : gel(v,j) = utoi(uel(E,j));
1671 : }
1672 384 : return y;
1673 : }
1674 :
1675 : static GEN
1676 2099 : FFX_zero(GEN ff, long v)
1677 : {
1678 2099 : GEN r = cgetg(3,t_POL);
1679 2099 : r[1] = evalvarn(v);
1680 2099 : gel(r,2) = FF_zero(ff);
1681 2099 : return r;
1682 : }
1683 :
1684 : GEN
1685 0 : FFX_add(GEN Pf, GEN Qf, GEN ff)
1686 : {
1687 0 : pari_sp av = avma;
1688 : GEN r,T,p;
1689 : ulong pp;
1690 0 : GEN P = FFX_to_raw(Pf, ff);
1691 0 : GEN Q = FFX_to_raw(Qf, ff);
1692 0 : _getFF(ff,&T,&p,&pp);
1693 0 : switch(ff[1])
1694 : {
1695 0 : case t_FF_FpXQ:
1696 0 : r = FpXX_add(P, Q, p);
1697 0 : break;
1698 0 : case t_FF_F2xq:
1699 0 : r = F2xX_add(P, Q);
1700 0 : break;
1701 0 : default:
1702 0 : r = FlxX_add(P, Q, pp);
1703 : }
1704 0 : if (!lgpol(r)) { set_avma(av); return FFX_zero(ff, varn(Pf)); }
1705 0 : return gc_GEN(av, raw_to_FFX(r, ff));
1706 : }
1707 :
1708 : static GEN
1709 92554 : FFX_wrap2(GEN Pf, GEN Qf, GEN ff, GEN FpXQX(GEN, GEN, GEN, GEN),
1710 : GEN F2xqX(GEN, GEN, GEN), GEN FlxqX(GEN, GEN, GEN, ulong))
1711 : {
1712 92554 : pari_sp av = avma;
1713 : GEN r,T,p;
1714 : ulong pp;
1715 92554 : GEN P = FFX_to_raw(Pf, ff);
1716 92554 : GEN Q = FFX_to_raw(Qf, ff);
1717 92554 : _getFF(ff,&T,&p,&pp);
1718 92554 : switch(ff[1])
1719 : {
1720 2036 : case t_FF_FpXQ:
1721 2036 : r = FpXQX(P, Q, T, p);
1722 2036 : break;
1723 34575 : case t_FF_F2xq:
1724 34575 : r = F2xqX(P, Q, T);
1725 34575 : break;
1726 55943 : default:
1727 55943 : r = FlxqX(P, Q, T, pp);
1728 : }
1729 92554 : if (!lgpol(r)) { set_avma(av); return FFX_zero(ff, varn(Pf)); }
1730 90455 : return gc_GEN(av, raw_to_FFX(r, ff));
1731 : }
1732 :
1733 : GEN
1734 89846 : FFX_mul(GEN Pf, GEN Qf, GEN ff)
1735 89846 : { return FFX_wrap2(Pf, Qf, ff, FpXQX_mul, F2xqX_mul, FlxqX_mul); }
1736 :
1737 : GEN
1738 2208 : FFX_gcd(GEN Pf, GEN Qf, GEN ff)
1739 2208 : { return FFX_wrap2(Pf, Qf, ff, FpXQX_gcd, F2xqX_gcd, FlxqX_gcd); }
1740 :
1741 : GEN
1742 2130 : FFX_sqr(GEN Pf, GEN ff)
1743 : {
1744 2130 : pari_sp av = avma;
1745 : GEN r,T,p;
1746 : ulong pp;
1747 2130 : GEN P = FFX_to_raw(Pf, ff);
1748 2130 : _getFF(ff,&T,&p,&pp);
1749 2130 : switch(ff[1])
1750 : {
1751 179 : case t_FF_FpXQ:
1752 179 : r = FpXQX_sqr(P, T, p);
1753 179 : break;
1754 911 : case t_FF_F2xq:
1755 911 : r = F2xqX_sqr(P, T);
1756 911 : break;
1757 1040 : default:
1758 1040 : r = FlxqX_sqr(P, T, pp);
1759 : }
1760 2130 : if (!lgpol(r)) { set_avma(av); return FFX_zero(ff, varn(Pf)); }
1761 2130 : return gc_GEN(av, raw_to_FFX(r, ff));
1762 : }
1763 :
1764 : GEN
1765 484 : FFX_rem(GEN Pf, GEN Qf, GEN ff)
1766 484 : { return FFX_wrap2(Pf, Qf, ff, FpXQX_rem, F2xqX_rem, FlxqX_rem); }
1767 :
1768 : GEN
1769 15 : FFX_resultant(GEN Pf, GEN Qf, GEN ff)
1770 : {
1771 15 : pari_sp av = avma;
1772 : GEN r,T,p;
1773 : ulong pp;
1774 15 : GEN P = FFX_to_raw(Pf, ff);
1775 15 : GEN Q = FFX_to_raw(Qf, ff);
1776 15 : GEN z = _initFF(ff,&T,&p,&pp);
1777 15 : switch(ff[1])
1778 : {
1779 5 : case t_FF_FpXQ:
1780 5 : r = FpXQX_resultant(P, Q, T, p);
1781 5 : break;
1782 5 : case t_FF_F2xq:
1783 5 : r = F2xqX_resultant(P, Q, T);
1784 5 : break;
1785 5 : default:
1786 5 : r = FlxqX_resultant(P, Q, T, pp);
1787 : }
1788 15 : return gc_upto(av, _mkFF(ff,z,r));
1789 : }
1790 :
1791 : GEN
1792 28 : FFX_disc(GEN Pf, GEN ff)
1793 : {
1794 28 : pari_sp av = avma;
1795 : GEN r,T,p;
1796 : ulong pp;
1797 28 : GEN P = FFX_to_raw(Pf, ff);
1798 28 : GEN z = _initFF(ff,&T,&p,&pp);
1799 28 : switch(ff[1])
1800 : {
1801 6 : case t_FF_FpXQ:
1802 6 : r = FpXQX_disc(P, T, p);
1803 6 : break;
1804 11 : case t_FF_F2xq:
1805 11 : r = F2xqX_disc(P, T);
1806 11 : break;
1807 11 : default:
1808 11 : r = FlxqX_disc(P, T, pp);
1809 : }
1810 28 : return gc_upto(av, _mkFF(ff,z,r));
1811 : }
1812 :
1813 : static GEN
1814 16 : gc_gcdext(pari_sp av, GEN r, GEN *u, GEN *v)
1815 : {
1816 16 : if (!u && !v) return gc_GEN(av, r);
1817 16 : if (u && v) return gc_all(av, 3, &r, u, v);
1818 0 : else return gc_all(av, 2, &r, u ? u: v);
1819 : }
1820 :
1821 : GEN
1822 16 : FFX_extgcd(GEN Pf, GEN Qf, GEN ff, GEN *pt_Uf, GEN *pt_Vf)
1823 : {
1824 16 : pari_sp av = avma;
1825 : GEN r,T,p;
1826 : ulong pp;
1827 16 : GEN P = FFX_to_raw(Pf, ff);
1828 16 : GEN Q = FFX_to_raw(Qf, ff);
1829 16 : _getFF(ff,&T,&p,&pp);
1830 16 : switch(ff[1])
1831 : {
1832 0 : case t_FF_FpXQ:
1833 0 : r = FpXQX_extgcd(P, Q, T, p, pt_Uf, pt_Vf);
1834 0 : break;
1835 5 : case t_FF_F2xq:
1836 5 : r = F2xqX_extgcd(P, Q, T, pt_Uf, pt_Vf);
1837 5 : break;
1838 11 : default:
1839 11 : r = FlxqX_extgcd(P, Q, T, pp, pt_Uf, pt_Vf);
1840 : }
1841 16 : if (pt_Uf) *pt_Uf = raw_to_FFX(*pt_Uf, ff);
1842 16 : if (pt_Vf) *pt_Vf = raw_to_FFX(*pt_Vf, ff);
1843 16 : return gc_gcdext(av, raw_to_FFX(r, ff), pt_Uf, pt_Vf);
1844 : }
1845 :
1846 : GEN
1847 0 : FFX_halfgcd(GEN Pf, GEN Qf, GEN ff)
1848 : {
1849 0 : pari_sp av = avma;
1850 : GEN r,T,p;
1851 : ulong pp;
1852 0 : GEN P = FFX_to_raw(Pf, ff);
1853 0 : GEN Q = FFX_to_raw(Qf, ff);
1854 0 : _getFF(ff,&T,&p,&pp);
1855 0 : switch(ff[1])
1856 : {
1857 0 : case t_FF_FpXQ:
1858 0 : r = FpXQX_halfgcd(P, Q, T, p);
1859 0 : break;
1860 0 : case t_FF_F2xq:
1861 0 : r = F2xqX_halfgcd(P, Q, T);
1862 0 : break;
1863 0 : default:
1864 0 : r = FlxqX_halfgcd(P, Q, T, pp);
1865 : }
1866 0 : return gc_GEN(av, raw_to_FFXM(r, ff));
1867 : }
1868 :
1869 : GEN
1870 18 : FFX_halfgcd_all(GEN Pf, GEN Qf, GEN ff, GEN *a, GEN *b)
1871 : {
1872 18 : pari_sp av = avma;
1873 : GEN r,T,p,R;
1874 : ulong pp;
1875 18 : GEN P = FFX_to_raw(Pf, ff);
1876 18 : GEN Q = FFX_to_raw(Qf, ff);
1877 18 : _getFF(ff,&T,&p,&pp);
1878 18 : switch(ff[1])
1879 : {
1880 6 : case t_FF_FpXQ:
1881 6 : r = FpXQX_halfgcd_all(P, Q, T, p, a, b);
1882 6 : break;
1883 6 : case t_FF_F2xq:
1884 6 : r = F2xqX_halfgcd_all(P, Q, T, a, b);
1885 6 : break;
1886 6 : default:
1887 6 : r = FlxqX_halfgcd_all(P, Q, T, pp, a, b);
1888 : }
1889 18 : if (*a) *a = raw_to_FFX(*a, ff);
1890 18 : if (*b) *b = raw_to_FFX(*b, ff);
1891 18 : R = raw_to_FFXM(r, ff);
1892 18 : return !a && b ? gc_all(av, 2, &R, b): gc_all(av, 1+!!a+!!b, &R, a, b);
1893 : }
1894 :
1895 : GEN
1896 6 : FFXQ_sqr(GEN Pf, GEN Qf, GEN ff)
1897 6 : { return FFX_wrap2(Pf, Qf, ff, FpXQXQ_sqr, F2xqXQ_sqr, FlxqXQ_sqr); }
1898 :
1899 : GEN
1900 10 : FFXQ_inv(GEN Pf, GEN Qf, GEN ff)
1901 10 : { return FFX_wrap2(Pf, Qf, ff, FpXQXQ_inv, F2xqXQ_inv, FlxqXQ_inv); }
1902 :
1903 : GEN
1904 88 : FFXQ_mul(GEN Pf, GEN Qf, GEN Sf, GEN ff)
1905 : {
1906 88 : pari_sp av = avma;
1907 : GEN r,T,p;
1908 : ulong pp;
1909 88 : GEN P = FFX_to_raw(Pf, ff);
1910 88 : GEN Q = FFX_to_raw(Qf, ff);
1911 88 : GEN S = FFX_to_raw(Sf, ff);
1912 88 : _getFF(ff,&T,&p,&pp);
1913 88 : switch(ff[1])
1914 : {
1915 24 : case t_FF_FpXQ:
1916 24 : r = FpXQXQ_mul(P, Q, S, T, p);
1917 24 : break;
1918 47 : case t_FF_F2xq:
1919 47 : r = F2xqXQ_mul(P, Q, S, T);
1920 47 : break;
1921 17 : default:
1922 17 : r = FlxqXQ_mul(P, Q, S, T, pp);
1923 : }
1924 88 : if (!lgpol(r)) { set_avma(av); return FFX_zero(ff, varn(Pf)); }
1925 88 : return gc_GEN(av, raw_to_FFX(r, ff));
1926 : }
1927 :
1928 : GEN
1929 66 : FFXQ_minpoly(GEN Pf, GEN Qf, GEN ff)
1930 : {
1931 66 : pari_sp av = avma;
1932 : GEN r,T,p;
1933 : ulong pp;
1934 66 : GEN P = FFX_to_raw(Pf, ff);
1935 66 : GEN Q = FFX_to_raw(Qf, ff);
1936 66 : _getFF(ff,&T,&p,&pp);
1937 66 : switch(ff[1])
1938 : {
1939 24 : case t_FF_FpXQ:
1940 24 : r = FpXQXQ_minpoly(P, Q, T, p);
1941 24 : break;
1942 24 : case t_FF_F2xq:
1943 24 : r = FlxX_to_F2xX(FlxqXQ_minpoly(F2xX_to_FlxX(P), F2xX_to_FlxX(Q), F2x_to_Flx(T), 2UL));
1944 24 : break;
1945 18 : default:
1946 18 : r = FlxqXQ_minpoly(P, Q, T, pp);
1947 : }
1948 66 : return gc_GEN(av, raw_to_FFX(r, ff));
1949 : }
1950 :
1951 : long
1952 144 : FFX_ispower(GEN Pf, long k, GEN ff, GEN *pt_r)
1953 : {
1954 144 : pari_sp av = avma;
1955 : GEN P,T,p;
1956 : ulong pp;
1957 : long s;
1958 144 : if (degpol(Pf) % k) return 0;
1959 144 : P = FFX_to_raw(Pf, ff);
1960 144 : _getFF(ff,&T,&p,&pp);
1961 144 : switch(ff[1])
1962 : {
1963 48 : case t_FF_FpXQ:
1964 48 : s = FpXQX_ispower(P, k, T, p, pt_r);
1965 48 : break;
1966 48 : case t_FF_F2xq:
1967 48 : s = F2xqX_ispower(P, k, T, pt_r);
1968 48 : break;
1969 48 : default:
1970 48 : s = FlxqX_ispower(P, k, T, pp, pt_r);
1971 : }
1972 144 : if (s==0) return gc_long(av,0);
1973 126 : if (pt_r)
1974 126 : *pt_r = gc_GEN(av, raw_to_FFX(*pt_r, ff));
1975 0 : else set_avma(av);
1976 126 : return 1;
1977 : }
1978 :
1979 : GEN
1980 378 : FFX_factor(GEN Pf, GEN ff)
1981 : {
1982 378 : pari_sp av = avma;
1983 : GEN r,T,p;
1984 : ulong pp;
1985 378 : GEN P = FFX_to_raw(Pf, ff);
1986 378 : _getFF(ff,&T,&p,&pp);
1987 378 : switch(ff[1])
1988 : {
1989 104 : case t_FF_FpXQ:
1990 104 : r = FpXQX_factor(P, T, p);
1991 104 : break;
1992 114 : case t_FF_F2xq:
1993 114 : r = F2xqX_factor(P, T);
1994 114 : break;
1995 160 : default:
1996 160 : r = FlxqX_factor(P, T, pp);
1997 : }
1998 378 : return gc_GEN(av, raw_to_FFX_fact(r, ff));
1999 : }
2000 :
2001 : GEN
2002 6 : FFX_factor_squarefree(GEN Pf, GEN ff)
2003 : {
2004 6 : pari_sp av = avma;
2005 : GEN r,T,p;
2006 : ulong pp;
2007 6 : GEN P = FFX_to_raw(Pf, ff);
2008 6 : _getFF(ff,&T,&p,&pp);
2009 6 : switch(ff[1])
2010 : {
2011 6 : case t_FF_FpXQ:
2012 6 : r = FpXQX_factor_squarefree(P, T, p);
2013 6 : break;
2014 0 : case t_FF_F2xq:
2015 0 : r = F2xqX_factor_squarefree(P, T);
2016 0 : break;
2017 0 : default:
2018 0 : r = FlxqX_factor_squarefree(P, T, pp);
2019 : }
2020 6 : return gc_GEN(av, raw_to_FFXC(r, ff));
2021 : }
2022 :
2023 : GEN
2024 6 : FFX_ddf(GEN Pf, GEN ff)
2025 : {
2026 6 : pari_sp av = avma;
2027 : GEN r,T,p;
2028 : ulong pp;
2029 6 : GEN P = FFX_to_raw(Pf, ff);
2030 6 : _getFF(ff,&T,&p,&pp);
2031 6 : switch(ff[1])
2032 : {
2033 6 : case t_FF_FpXQ:
2034 6 : r = FpXQX_ddf(P, T, p);
2035 6 : break;
2036 0 : case t_FF_F2xq:
2037 0 : r = F2xqX_ddf(P, T);
2038 0 : break;
2039 0 : default:
2040 0 : r = FlxqX_ddf(P, T, pp);
2041 : }
2042 6 : return gc_GEN(av, raw_to_FFX_fact(r, ff));
2043 : }
2044 :
2045 : GEN
2046 108 : FFX_degfact(GEN Pf, GEN ff)
2047 : {
2048 108 : pari_sp av = avma;
2049 : GEN r,T,p;
2050 : ulong pp;
2051 108 : GEN P = FFX_to_raw(Pf, ff);
2052 108 : _getFF(ff, &T, &p, &pp);
2053 108 : switch(ff[1])
2054 : {
2055 36 : case t_FF_FpXQ:
2056 36 : r = FpXQX_degfact(P, T, p);
2057 36 : break;
2058 36 : case t_FF_F2xq:
2059 36 : r = F2xqX_degfact(P, T);
2060 36 : break;
2061 36 : default:
2062 36 : r = FlxqX_degfact(P, T, pp);
2063 : }
2064 108 : return gc_GEN(av, r);
2065 : }
2066 :
2067 : GEN
2068 0 : FqX_to_FFX(GEN x, GEN ff)
2069 0 : { pari_APPLY_pol_normalized(Fq_to_FF(gel(x,i), ff)); }
2070 :
2071 : GEN
2072 3198 : FqC_to_FFC(GEN x, GEN ff)
2073 12372 : { pari_APPLY_type(t_COL, Fq_to_FF(gel(x,i), ff)) }
2074 :
2075 : GEN
2076 0 : FqV_to_FFV(GEN x, GEN ff)
2077 0 : { pari_APPLY_type(t_VEC, Fq_to_FF(gel(x,i), ff)) }
2078 :
2079 : GEN
2080 1470 : FqM_to_FFM(GEN x, GEN ff)
2081 4668 : { pari_APPLY_same(FqC_to_FFC(gel(x,i), ff)) }
2082 :
2083 : GEN
2084 2500 : ffgen(GEN T, long v)
2085 : {
2086 2500 : GEN A, p = NULL, ff = cgetg(5,t_FFELT);
2087 : long d;
2088 2500 : switch(typ(T))
2089 : {
2090 224 : case t_FFELT:
2091 224 : if (v < 0 || FF_var(T)==v) return FF_gen(T);
2092 0 : p = FF_p_i(T); T = FF_mod(T); d = degpol(T);
2093 0 : break;
2094 328 : case t_POL:
2095 328 : d = degpol(T); p = NULL;
2096 328 : if (d < 1 || !RgX_is_FpX(T, &p) || !p) pari_err_TYPE("ffgen",T);
2097 328 : T = RgX_to_FpX(T, p);
2098 : /* testing for irreducibility is too costly */
2099 328 : if (!FpX_is_squarefree(T,p)) pari_err_IRREDPOL("ffgen",T);
2100 322 : break;
2101 1306 : case t_INT:
2102 1306 : d = ispseudoprimepower(T,&p);
2103 1306 : if (!d) pari_err_PRIME("ffgen",T);
2104 1306 : T = init_Fq(p, d, v);
2105 1306 : break;
2106 642 : case t_VEC: case t_COL:
2107 642 : if (lg(T) == 3) {
2108 642 : p = gel(T,1);
2109 642 : A = gel(T,2);
2110 642 : if (typ(p) == t_INT && typ(A) == t_INT)
2111 : {
2112 642 : d = itos(A);
2113 642 : T = init_Fq(p, d, v);
2114 642 : break;
2115 : }
2116 : }
2117 : default:
2118 0 : pari_err_TYPE("ffgen",T);
2119 : return NULL;/* LCOV_EXCL_LINE */
2120 : }
2121 2270 : if (v < 0) v = varn(T);
2122 2270 : if (lgefint(p)==3)
2123 : {
2124 1950 : ulong pp = p[2];
2125 1950 : long sv = evalvarn(v);
2126 1950 : if (pp==2)
2127 : {
2128 641 : ff[1] = t_FF_F2xq;
2129 641 : T = ZX_to_F2x(T); T[1] = sv;
2130 641 : A = polx_F2x(sv); if (d == 1) A = F2x_rem(A, T);
2131 641 : p = gen_2;
2132 : }
2133 : else
2134 : {
2135 1309 : ff[1] = t_FF_Flxq;
2136 1309 : T = ZX_to_Flx(T,pp); T[1] = sv;
2137 1309 : A = polx_Flx(sv); if (d == 1) A = Flx_rem(A, T, pp);
2138 1309 : p = icopy(p);
2139 : }
2140 : }
2141 : else
2142 : {
2143 320 : ff[1] = t_FF_FpXQ;
2144 320 : setvarn(T,v);
2145 320 : A = pol_x(v); if (d == 1) A = FpX_rem(A, T, p);
2146 320 : p = icopy(p);
2147 : }
2148 2270 : gel(ff,2) = A;
2149 2270 : gel(ff,3) = T;
2150 2270 : gel(ff,4) = p; return ff;
2151 : }
2152 :
2153 : GEN
2154 8196 : p_to_FF(GEN p, long v)
2155 : {
2156 : GEN A, T;
2157 8196 : GEN ff = cgetg(5,t_FFELT);
2158 8196 : if (lgefint(p)==3)
2159 : {
2160 8195 : ulong pp = p[2];
2161 8195 : long sv = evalvarn(v);
2162 8195 : if (pp==2)
2163 : {
2164 210 : ff[1] = t_FF_F2xq;
2165 210 : T = polx_F2x(sv);
2166 210 : A = pol1_F2x(sv);
2167 210 : p = gen_2;
2168 : }
2169 : else
2170 : {
2171 7985 : ff[1] = t_FF_Flxq;
2172 7985 : T = polx_Flx(sv);
2173 7985 : A = pol1_Flx(sv);
2174 7985 : p = icopy(p);
2175 : }
2176 : }
2177 : else
2178 : {
2179 1 : ff[1] = t_FF_FpXQ;
2180 1 : T = pol_x(v);
2181 1 : A = pol_1(v);
2182 1 : p = icopy(p);
2183 : }
2184 8196 : gel(ff,2) = A;
2185 8196 : gel(ff,3) = T;
2186 8196 : gel(ff,4) = p; return ff;
2187 : }
2188 : GEN
2189 7050 : Tp_to_FF(GEN T, GEN p)
2190 : {
2191 : GEN A, ff;
2192 : long v;
2193 7050 : if (!T) return p_to_FF(p,0);
2194 6984 : ff = cgetg(5,t_FFELT);
2195 6984 : v = varn(T);
2196 6984 : if (lgefint(p)==3)
2197 : {
2198 6796 : ulong pp = p[2];
2199 6796 : long sv = evalvarn(v);
2200 6796 : if (pp==2)
2201 : {
2202 1650 : ff[1] = t_FF_F2xq;
2203 1650 : T = ZX_to_F2x(T);
2204 1650 : A = pol1_F2x(sv);
2205 1650 : p = gen_2;
2206 : }
2207 : else
2208 : {
2209 5146 : ff[1] = t_FF_Flxq;
2210 5146 : T = ZX_to_Flx(T, pp);
2211 5146 : A = pol1_Flx(sv);
2212 5146 : p = icopy(p);
2213 : }
2214 : }
2215 : else
2216 : {
2217 188 : ff[1] = t_FF_FpXQ;
2218 188 : T = ZX_copy(T);
2219 188 : A = pol_1(v);
2220 188 : p = icopy(p);
2221 : }
2222 6984 : gel(ff,2) = A;
2223 6984 : gel(ff,3) = T;
2224 6984 : gel(ff,4) = p; return ff;
2225 : }
2226 :
2227 : GEN
2228 54 : fforder(GEN x, GEN o)
2229 : {
2230 54 : if (typ(x)!=t_FFELT) pari_err_TYPE("fforder",x);
2231 54 : if (FF_equal0(x)) pari_err_DOMAIN("fforder", "x", "=", gen_0, x);
2232 48 : return FF_order(x,o);
2233 : }
2234 :
2235 : GEN
2236 348 : ffprimroot(GEN x, GEN *o)
2237 : {
2238 348 : if (typ(x)!=t_FFELT) pari_err_TYPE("ffprimroot",x);
2239 348 : return FF_primroot(x,o);
2240 : }
2241 :
2242 : GEN
2243 162 : fflog(GEN x, GEN g, GEN o)
2244 : {
2245 162 : if (typ(x)!=t_FFELT) pari_err_TYPE("fflog",x);
2246 162 : if (typ(g)!=t_FFELT) pari_err_TYPE("fflog",g);
2247 162 : return FF_log(x,g,o);
2248 : }
2249 :
2250 : GEN
2251 158133 : ffrandom(GEN ff)
2252 : {
2253 : ulong pp;
2254 158133 : GEN r, T, p, z = _initFF(ff,&T,&p,&pp);
2255 158133 : switch(ff[1])
2256 : {
2257 8689 : case t_FF_FpXQ:
2258 8689 : r = random_FpX(degpol(T), varn(T), p);
2259 8689 : break;
2260 100517 : case t_FF_F2xq:
2261 100517 : r = random_F2x(F2x_degree(T), T[1]);
2262 100517 : break;
2263 48927 : default:
2264 48927 : r = random_Flx(degpol(T), T[1], pp);
2265 : }
2266 158133 : return _mkFF(ff,z,r);
2267 : }
2268 :
2269 : int
2270 0 : Rg_is_FF(GEN c, GEN *ff)
2271 : {
2272 0 : switch(typ(c))
2273 : {
2274 0 : case t_FFELT:
2275 0 : if (!*ff) *ff = c;
2276 0 : else if (!FF_samefield(*ff, c)) return 0;
2277 0 : break;
2278 0 : default:
2279 0 : return 0;
2280 : }
2281 0 : return 1;
2282 : }
2283 :
2284 : int
2285 0 : RgC_is_FFC(GEN x, GEN *ff)
2286 : {
2287 0 : long i, lx = lg(x);
2288 0 : for (i=lx-1; i>0; i--)
2289 0 : if (!Rg_is_FF(gel(x,i), ff)) return 0;
2290 0 : return (*ff != NULL);
2291 : }
2292 :
2293 : int
2294 0 : RgM_is_FFM(GEN x, GEN *ff)
2295 : {
2296 0 : long j, lx = lg(x);
2297 0 : for (j=lx-1; j>0; j--)
2298 0 : if (!RgC_is_FFC(gel(x,j), ff)) return 0;
2299 0 : return (*ff != NULL);
2300 : }
2301 :
2302 : static GEN
2303 2987 : FqC_to_FpXQC(GEN x, GEN T, GEN p)
2304 48100 : { pari_APPLY_same(Fq_to_FpXQ(gel(x,i), T, p)); }
2305 :
2306 : static GEN
2307 410 : FqM_to_FpXQM(GEN x, GEN T, GEN p)
2308 3331 : { pari_APPLY_same(FqC_to_FpXQC(gel(x,i), T, p)); }
2309 :
2310 : /* for functions t_MAT -> t_MAT */
2311 : static GEN
2312 264 : FFM_wrap(GEN M, GEN ff, GEN (*Fq)(GEN,GEN,GEN),
2313 : GEN (*Flxq)(GEN,GEN,ulong),
2314 : GEN (*F2xq)(GEN,GEN))
2315 : {
2316 264 : pari_sp av = avma;
2317 : ulong pp;
2318 : GEN T, p;
2319 264 : _getFF(ff,&T,&p,&pp); M = FFM_to_raw(M, ff);
2320 264 : switch(ff[1])
2321 : {
2322 102 : case t_FF_FpXQ: M = Fq(M,T,p); if (M) M = FqM_to_FpXQM(M,T,p);
2323 102 : break;
2324 60 : case t_FF_F2xq: M = F2xq(M,T); break;
2325 102 : default: M = Flxq(M,T,pp); break;
2326 : }
2327 264 : if (!M) return gc_NULL(av);
2328 234 : return gc_GEN(av, raw_to_FFM(M, ff));
2329 : }
2330 :
2331 : /* for functions (t_MAT, t_MAT) -> t_MAT */
2332 : static GEN
2333 3996 : FFM_FFM_wrap(GEN M, GEN N, GEN ff,
2334 : GEN (*Fq)(GEN, GEN, GEN, GEN),
2335 : GEN (*Flxq)(GEN, GEN, GEN, ulong),
2336 : GEN (*F2xq)(GEN, GEN, GEN))
2337 : {
2338 3996 : pari_sp av = avma;
2339 : ulong pp;
2340 : GEN T, p;
2341 3996 : int is_sqr = M==N;
2342 3996 : _getFF(ff, &T, &p, &pp);
2343 3996 : M = FFM_to_raw(M, ff);
2344 3996 : N = is_sqr? M: FFM_to_raw(N, ff);
2345 3996 : switch(ff[1])
2346 : {
2347 350 : case t_FF_FpXQ: M = Fq(M, N, T, p); if (M) M = FqM_to_FpXQM(M, T, p);
2348 350 : break;
2349 1278 : case t_FF_F2xq: M = F2xq(M, N, T); break;
2350 2368 : default: M = Flxq(M, N, T, pp); break;
2351 : }
2352 3996 : if (!M) return gc_NULL(av);
2353 3918 : return gc_GEN(av, raw_to_FFM(M, ff));
2354 : }
2355 :
2356 : /* for functions (t_MAT, t_COL) -> t_COL */
2357 : static GEN
2358 180 : FFM_FFC_wrap(GEN M, GEN C, GEN ff,
2359 : GEN (*Fq)(GEN, GEN, GEN, GEN),
2360 : GEN (*Flxq)(GEN, GEN, GEN, ulong),
2361 : GEN (*F2xq)(GEN, GEN, GEN))
2362 : {
2363 180 : pari_sp av = avma;
2364 : ulong pp;
2365 : GEN T, p;
2366 180 : _getFF(ff, &T, &p, &pp);
2367 180 : M = FFM_to_raw(M, ff);
2368 180 : C = FFC_to_raw(C, ff);
2369 180 : switch(ff[1])
2370 : {
2371 60 : case t_FF_FpXQ: C = Fq(M, C, T, p); if (C) C = FqC_to_FpXQC(C, T, p);
2372 60 : break;
2373 60 : case t_FF_F2xq: C = F2xq(M, C, T); break;
2374 60 : default: C = Flxq(M, C, T, pp); break;
2375 : }
2376 180 : if (!C) return gc_NULL(av);
2377 138 : return gc_GEN(av, raw_to_FFC(C, ff));
2378 : }
2379 :
2380 : GEN
2381 66 : FFM_ker(GEN M, GEN ff)
2382 66 : { return FFM_wrap(M,ff, &FqM_ker,&FlxqM_ker,&F2xqM_ker); }
2383 : GEN
2384 42 : FFM_image(GEN M, GEN ff)
2385 42 : { return FFM_wrap(M,ff, &FqM_image,&FlxqM_image,&F2xqM_image); }
2386 : GEN
2387 126 : FFM_inv(GEN M, GEN ff)
2388 126 : { return FFM_wrap(M,ff, &FqM_inv,&FlxqM_inv,&F2xqM_inv); }
2389 : GEN
2390 30 : FFM_suppl(GEN M, GEN ff)
2391 30 : { return FFM_wrap(M,ff, &FqM_suppl,&FlxqM_suppl,&F2xqM_suppl); }
2392 :
2393 : GEN
2394 72 : FFM_deplin(GEN M, GEN ff)
2395 : {
2396 72 : pari_sp av = avma;
2397 : ulong pp;
2398 : GEN C, T, p;
2399 72 : _getFF(ff, &T, &p, &pp); M = FFM_to_raw(M, ff);
2400 72 : switch(ff[1]) {
2401 30 : case t_FF_FpXQ: C = FqM_deplin(M, T, p);
2402 30 : if (C) C = FqC_to_FpXQC(C, T, p); break;
2403 12 : case t_FF_F2xq: C = F2xqM_deplin(M, T); break;
2404 30 : default: C = FlxqM_deplin(M, T, pp); break;
2405 : }
2406 72 : if (!C) return gc_NULL(av);
2407 42 : return gc_GEN(av, raw_to_FFC(C, ff));
2408 : }
2409 :
2410 : GEN
2411 18 : FFM_indexrank(GEN M, GEN ff)
2412 : {
2413 18 : pari_sp av = avma;
2414 : ulong pp;
2415 : GEN R, T, p;
2416 18 : _getFF(ff,&T,&p,&pp); M = FFM_to_raw(M, ff);
2417 18 : switch(ff[1]) {
2418 6 : case t_FF_FpXQ: R = FqM_indexrank(M,T,p); break;
2419 6 : case t_FF_F2xq: R = F2xqM_indexrank(M,T); break;
2420 6 : default: R = FlxqM_indexrank(M,T,pp); break;
2421 : }
2422 18 : return gc_upto(av, R);
2423 : }
2424 :
2425 : long
2426 60 : FFM_rank(GEN M, GEN ff)
2427 : {
2428 60 : pari_sp av = avma;
2429 : long r;
2430 : ulong pp;
2431 : GEN T, p;
2432 60 : _getFF(ff,&T,&p,&pp); M = FFM_to_raw(M, ff);
2433 60 : switch(ff[1])
2434 : {
2435 24 : case t_FF_FpXQ: r = FqM_rank(M,T,p); break;
2436 6 : case t_FF_F2xq: r = F2xqM_rank(M,T); break;
2437 30 : default: r = FlxqM_rank(M,T,pp); break;
2438 : }
2439 60 : return gc_long(av,r);
2440 : }
2441 : GEN
2442 54 : FFM_det(GEN M, GEN ff)
2443 : {
2444 54 : pari_sp av = avma;
2445 : ulong pp;
2446 : GEN d, T, p;
2447 54 : _getFF(ff,&T,&p,&pp); M = FFM_to_raw(M, ff);
2448 54 : switch(ff[1])
2449 : {
2450 24 : case t_FF_FpXQ: d = FqM_det(M,T,p); break;
2451 6 : case t_FF_F2xq: d = F2xqM_det(M,T); break;
2452 24 : default: d = FlxqM_det(M,T,pp); break;
2453 : }
2454 54 : return gc_GEN(av, mkFF_i(ff, d));
2455 : }
2456 :
2457 : GEN
2458 48 : FFM_FFC_gauss(GEN M, GEN C, GEN ff)
2459 : {
2460 48 : return FFM_FFC_wrap(M, C, ff, FqM_FqC_gauss,
2461 : FlxqM_FlxqC_gauss, F2xqM_F2xqC_gauss);
2462 : }
2463 :
2464 : GEN
2465 54 : FFM_gauss(GEN M, GEN N, GEN ff)
2466 : {
2467 54 : return FFM_FFM_wrap(M, N, ff, FqM_gauss,
2468 : FlxqM_gauss, F2xqM_gauss);
2469 : }
2470 :
2471 : GEN
2472 54 : FFM_FFC_invimage(GEN M, GEN C, GEN ff)
2473 : {
2474 54 : return FFM_FFC_wrap(M, C, ff, FqM_FqC_invimage,
2475 : FlxqM_FlxqC_invimage, F2xqM_F2xqC_invimage);
2476 : }
2477 :
2478 : GEN
2479 90 : FFM_invimage(GEN M, GEN N, GEN ff)
2480 : {
2481 90 : return FFM_FFM_wrap(M, N, ff, FqM_invimage,
2482 : FlxqM_invimage, F2xqM_invimage);
2483 : }
2484 :
2485 : GEN
2486 78 : FFM_FFC_mul(GEN M, GEN C, GEN ff)
2487 : {
2488 78 : return FFM_FFC_wrap(M, C, ff, FqM_FqC_mul,
2489 : FlxqM_FlxqC_mul, F2xqM_F2xqC_mul);
2490 : }
2491 :
2492 : GEN
2493 3852 : FFM_mul(GEN M, GEN N, GEN ff)
2494 : {
2495 3852 : return FFM_FFM_wrap(M, N, ff, FqM_mul, FlxqM_mul, F2xqM_mul);
2496 : }
2497 :
2498 : GEN
2499 12 : FFV_roots_to_pol(GEN V, GEN ff, long v)
2500 : {
2501 12 : pari_sp av = avma;
2502 : ulong pp;
2503 : GEN P, T, p;
2504 12 : long w = fetch_var_higher();
2505 12 : _getFF(ff,&T,&p,&pp); V = FFC_to_raw(V, ff);
2506 12 : switch(ff[1])
2507 : {
2508 0 : case t_FF_FpXQ: P = FqV_roots_to_pol(V, T, p, w); break;
2509 6 : case t_FF_F2xq: P = F2xqV_roots_to_pol(V, T, w); break;
2510 6 : default: P = FlxqV_roots_to_pol(V, T, pp, w); break;
2511 : }
2512 12 : if (!P) return gc_NULL(av);
2513 12 : P = raw_to_FFX(P, ff);
2514 12 : setvarn(P, v);
2515 12 : delete_var();
2516 12 : return gc_GEN(av, P);
2517 : }
|