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