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