Code coverage tests

This page documents the degree to which the PARI/GP source code is tested by our public test suite, distributed with the source distribution in directory src/test/. This is measured by the gcov utility; we then process gcov output using the lcov frond-end.

We test a few variants depending on Configure flags on the pari.math.u-bordeaux.fr machine (x86_64 architecture), and agregate them in the final report:

The target is 90% coverage for all mathematical modules (given that branches depending on DEBUGLEVEL or DEBUGMEM are not covered). This script is run to produce the results below.

LCOV - code coverage report
Current view: top level - basemath - FF.c (source / functions) Hit Total Coverage
Test: PARI/GP v2.10.0 lcov report (development 21741-70cf009) Lines: 1042 1102 94.6 %
Date: 2018-01-21 06:18:30 Functions: 125 131 95.4 %
Legend: Lines: hit not hit

          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. It is distributed in the hope that it will be useful, but WITHOUT
       8             : ANY WARRANTY WHATSOEVER.
       9             : 
      10             : Check the License for details. You should have received a copy of it, along
      11             : with the package; see the file 'COPYING'. If not, write to the Free Software
      12             : Foundation, Inc., 51 Franklin Street, Fifth Floor, Boston, MA 02110-1301 USA. */
      13             : 
      14             : #include "pari.h"
      15             : #include "paripriv.h"
      16             : 
      17             : /*************************************************************************/
      18             : /**                                                                     **/
      19             : /**                   Routines for handling FFELT                       **/
      20             : /**                                                                     **/
      21             : /*************************************************************************/
      22             : 
      23             : /*************************************************************************/
      24             : /**                                                                     **/
      25             : /**                   Low-level constructors                            **/
      26             : /**                                                                     **/
      27             : /*************************************************************************/
      28             : 
      29             : INLINE void
      30    26991879 : _getFF(GEN x, GEN *T, GEN *p, ulong *pp)
      31             : {
      32    26991879 :   *T=gel(x,3);
      33    26991879 :   *p=gel(x,4);
      34    26991879 :   *pp=(*p)[2];
      35    26991879 : }
      36             : 
      37             : INLINE GEN
      38    25532241 : _initFF(GEN x, GEN *T, GEN *p, ulong *pp)
      39             : {
      40    25532241 :   _getFF(x,T,p,pp);
      41    25532246 :   return cgetg(5,t_FFELT);
      42             : }
      43             : 
      44             : INLINE void
      45    18732039 : _checkFF(GEN x, GEN y, const char *s)
      46    18732039 : { if (!FF_samefield(x,y)) pari_err_OP(s,x,y); }
      47             : 
      48             : INLINE GEN
      49    25447053 : _mkFF(GEN x, GEN z, GEN r)
      50             : {
      51    25447053 :   z[1]=x[1];
      52    25447053 :   gel(z,2)=r;
      53    25447053 :   gel(z,3)=gcopy(gel(x,3));
      54    25447057 :   gel(z,4)=icopy(gel(x,4));
      55    25447057 :   return z;
      56             : }
      57             : 
      58             : INLINE GEN
      59     2968947 : _mkFF_i(GEN x, GEN z, GEN r)
      60             : {
      61     2968947 :   z[1]=x[1];
      62     2968947 :   gel(z,2)=r;
      63     2968947 :   gel(z,3)=gel(x,3);
      64     2968947 :   gel(z,4)=gel(x,4);
      65     2968947 :   return z;
      66             : }
      67             : 
      68             : INLINE GEN
      69     2883750 : mkFF_i(GEN x, GEN r)
      70             : {
      71     2883750 :   GEN z = cgetg(5,t_FFELT);
      72     2883751 :   return _mkFF_i(x,z,r);
      73             : }
      74             : 
      75             : /*************************************************************************/
      76             : /**                                                                     **/
      77             : /**                   medium-level constructors                         **/
      78             : /**                                                                     **/
      79             : /*************************************************************************/
      80             : 
      81             : static GEN
      82      430192 : Z_to_raw(GEN x, GEN ff)
      83             : {
      84             :   ulong pp;
      85             :   GEN T, p;
      86      430192 :   _getFF(ff,&T,&p,&pp);
      87      430192 :   switch(ff[1])
      88             :   {
      89             :   case t_FF_FpXQ:
      90       10745 :     return scalarpol(x, varn(T));
      91             :   case t_FF_F2xq:
      92      201495 :     return Z_to_F2x(x,varn(T));
      93             :   default:
      94      217952 :     return Z_to_Flx(x,pp,T[1]);
      95             :   }
      96             : }
      97             : 
      98             : static GEN
      99     1754892 : Rg_to_raw(GEN x, GEN ff)
     100             : {
     101     1754892 :   long tx = typ(x);
     102     1754892 :   switch(tx)
     103             :   {
     104             :   case t_INT: case t_FRAC: case t_PADIC: case t_INTMOD:
     105      430192 :     return Z_to_raw(Rg_to_Fp(x, FF_p_i(ff)), ff);
     106             :   case t_FFELT:
     107     1324700 :     if (!FF_samefield(x,ff))
     108           0 :       pari_err_MODULUS("Rg_to_raw",x,ff);
     109     1324700 :     return gel(x,2);
     110             :   }
     111           0 :   pari_err_TYPE("Rg_to_raw",x);
     112             :   return NULL; /* LCOV_EXCL_LINE */
     113             : }
     114             : 
     115             : static GEN
     116      316045 : FFX_to_raw(GEN x, GEN ff)
     117             : {
     118             :   long i, lx;
     119      316045 :   GEN y = cgetg_copy(x,&lx);
     120      316045 :   y[1] = x[1];
     121     1842569 :   for(i=2; i<lx; i++)
     122     1526524 :     gel(y, i) = Rg_to_raw(gel(x, i), ff);
     123      316045 :   switch (ff[1])
     124             :   {
     125             :   case t_FF_FpXQ:
     126        4742 :     return FpXX_renormalize(y, lx);
     127             :   case t_FF_F2xq:
     128      136898 :     return F2xX_renormalize(y, lx);
     129             :   default:
     130      174405 :     return FlxX_renormalize(y, lx);
     131             :   }
     132             : }
     133             : 
     134             : static GEN
     135       13804 : FFC_to_raw(GEN x, GEN ff)
     136       13804 : { pari_APPLY_same(Rg_to_raw(gel(x, i), ff)) }
     137             : 
     138             : static GEN
     139        1470 : FFM_to_raw(GEN x, GEN ff)
     140        1470 : { pari_APPLY_same(FFC_to_raw(gel(x, i), ff)) }
     141             : 
     142             : /* in place */
     143             : static GEN
     144      632927 : rawFq_to_FF(GEN x, GEN ff)
     145             : {
     146      632927 :   return mkFF_i(ff, typ(x)==t_INT ? scalarpol(x, varn(gel(ff,3))): x);
     147             : }
     148             : 
     149             : /* in place */
     150             : static GEN
     151      105893 : raw_to_FFX(GEN x, GEN ff)
     152             : {
     153      105893 :   long i, lx = lg(x);
     154      105893 :   for (i=2; i<lx; i++) gel(x,i) = rawFq_to_FF(gel(x,i), ff);
     155      105893 :   return x;
     156             : }
     157             : 
     158             : /* in place */
     159             : static GEN
     160      108458 : raw_to_FFC(GEN x, GEN ff)
     161             : {
     162      108458 :   long i, lx = lg(x);
     163      108458 :   for (i=1; i<lx; i++) gel(x,i) = mkFF_i(ff, gel(x,i));
     164      108458 :   return x;
     165             : }
     166             : 
     167             : /* in place */
     168             : static GEN
     169         567 : raw_to_FFM(GEN x, GEN ff)
     170             : {
     171         567 :   long i, lx = lg(x);
     172         567 :   for (i=1; i<lx; i++) gel(x,i) = raw_to_FFC(gel(x,i), ff);
     173         567 :   return x;
     174             : }
     175             : 
     176             : GEN
     177         161 : Fq_to_FF(GEN x, GEN ff)
     178             : {
     179             :   ulong pp;
     180         161 :   GEN r, T, p, z=_initFF(ff,&T,&p,&pp);
     181         161 :   int is_int = typ(x)==t_INT;
     182         161 :   switch(ff[1])
     183             :   {
     184             :   case t_FF_FpXQ:
     185           0 :     r= is_int ? scalarpol(x, varn(T)): x;
     186           0 :     break;
     187             :   case t_FF_F2xq:
     188         126 :     r= is_int ? Z_to_F2x(x,T[1]): ZX_to_F2x(x);
     189         126 :     break;
     190             :   default:
     191          35 :     r= is_int ? Z_to_Flx(x,pp,T[1]): ZX_to_Flx(x,pp);
     192             :   }
     193         161 :   setvarn(r, varn(T)); /* paranoia */
     194         161 :   return _mkFF_i(ff,z,r);
     195             : }
     196             : 
     197             : /*************************************************************************/
     198             : /**                                                                     **/
     199             : /**                   Public functions                                  **/
     200             : /**                                                                     **/
     201             : /*************************************************************************/
     202             : 
     203             : /* Return true if x and y are defined in the same field */
     204             : 
     205             : static int
     206    21394472 : FF_samechar(GEN x, GEN y)
     207             : {
     208    21394472 :   return x[1] == y[1] && equalii(gel(x,4),gel(y,4));
     209             : }
     210             : 
     211             : int
     212    21394472 : FF_samefield(GEN x, GEN y)
     213             : {
     214    21394472 :   return FF_samechar(x, y) && gidentical(gel(x,3),gel(y,3));
     215             : }
     216             : 
     217             : int
     218       39683 : FF_equal(GEN x, GEN y)
     219       39683 : { return FF_samefield(x,y) && gidentical(gel(x,2),gel(y,2)); }
     220             : 
     221             : int
     222     8812732 : FF_equal0(GEN x)
     223             : {
     224     8812732 :   return lgpol(gel(x,2))==0;
     225             : }
     226             : 
     227             : int
     228       14868 : FF_equal1(GEN x)
     229             : {
     230       14868 :   GEN A = gel(x,2);
     231       14868 :   switch(x[1])
     232             :   {
     233             :   case t_FF_FpXQ:
     234        1441 :     return degpol(A)==0 && gequal1(gel(A,2));
     235             :   default:
     236       13427 :     return degpol(A)==0 && A[2]==1;
     237             :   }
     238             : }
     239             : 
     240             : static int
     241           0 : Fp_cmp_1(GEN x, GEN p)
     242             : {
     243           0 :   pari_sp av = avma;
     244           0 :   int b = equalii(x, addis(p,-1));
     245           0 :   avma = av; return b;
     246             : }
     247             : 
     248             : int
     249          42 : FF_equalm1(GEN x)
     250             : {
     251             :   ulong pp;
     252          42 :   GEN T, p, y = gel(x,2);
     253          42 :   _getFF(x,&T,&p,&pp);
     254          42 :   switch(x[1])
     255             :   {
     256             :   case t_FF_FpXQ:
     257           0 :     return (degpol(y) == 0 && Fp_cmp_1(gel(y,2), p));
     258             :   default:
     259          42 :     return (degpol(y) == 0 && uel(y,2) == pp-1);
     260             :   }
     261             : }
     262             : 
     263             : GEN
     264        5961 : FF_zero(GEN x)
     265             : {
     266             :   ulong pp;
     267        5961 :   GEN r, T, p, z=_initFF(x,&T,&p,&pp);
     268        5961 :   switch(x[1])
     269             :   {
     270             :   case t_FF_FpXQ:
     271         329 :     r=zeropol(varn(T));
     272         329 :     break;
     273             :   case t_FF_F2xq:
     274        1495 :     r=zero_F2x(T[1]);
     275        1495 :     break;
     276             :   default:
     277        4137 :     r=zero_Flx(T[1]);
     278             :   }
     279        5961 :   return _mkFF(x,z,r);
     280             : }
     281             : 
     282             : GEN
     283       16072 : FF_1(GEN x)
     284             : {
     285             :   ulong pp;
     286       16072 :   GEN r, T, p, z=_initFF(x,&T,&p,&pp);
     287       16072 :   switch(x[1])
     288             :   {
     289             :   case t_FF_FpXQ:
     290         181 :     r=pol_1(varn(T));
     291         181 :     break;
     292             :   case t_FF_F2xq:
     293        9086 :     r=pol1_F2x(T[1]);
     294        9086 :     break;
     295             :   default:
     296        6805 :     r=pol1_Flx(T[1]);
     297             :   }
     298       16072 :   return _mkFF(x,z,r);
     299             : }
     300             : 
     301             : GEN
     302         469 : FF_gen(GEN x)
     303             : {
     304             :   ulong pp;
     305         469 :   GEN r, T, p, z=_initFF(x,&T,&p,&pp);
     306         469 :   switch(x[1])
     307             :   {
     308             :   case t_FF_FpXQ:
     309         154 :     r = pol_x(varn(T));
     310         154 :     if (degpol(T)==1) r = FpX_rem(r, T, p);
     311         154 :     break;
     312             :   case t_FF_F2xq:
     313         154 :     r = polx_F2x(T[1]);
     314         154 :     if (F2x_degree(T)==1) r = F2x_rem(r, T);
     315         154 :     break;
     316             :   default:
     317         161 :     r = polx_Flx(T[1]);
     318         161 :     if (degpol(T)==1) r = Flx_rem(r, T, pp);
     319             :   }
     320         469 :   return _mkFF(x,z,r);
     321             : }
     322             : GEN
     323       54481 : FF_q(GEN x)
     324             : {
     325             :   ulong pp;
     326             :   GEN T, p;
     327       54481 :   _getFF(x,&T,&p,&pp);
     328       54481 :   switch(x[1])
     329             :   {
     330             :   case t_FF_FpXQ:
     331        1681 :     return powiu(p, degpol(T));
     332             :     break;
     333             :   case t_FF_F2xq:
     334        7028 :     return int2n(F2x_degree(T));
     335             :     break;
     336             :   default:
     337       45772 :     return powuu(pp,degpol(T));
     338             :   }
     339             : }
     340             : 
     341             : GEN
     342          56 : FF_p(GEN x)
     343             : {
     344          56 :   return icopy(gel(x,4));
     345             : }
     346             : 
     347             : GEN
     348     2072080 : FF_p_i(GEN x)
     349             : {
     350     2072080 :   return gel(x,4);
     351             : }
     352             : 
     353             : GEN
     354      165214 : FF_mod(GEN x)
     355             : {
     356      165214 :   switch(x[1])
     357             :   {
     358             :   case t_FF_FpXQ:
     359         267 :     return ZX_copy(gel(x,3));
     360             :   case t_FF_F2xq:
     361         273 :     return F2x_to_ZX(gel(x,3));
     362             :   default:
     363      164674 :     return Flx_to_ZX(gel(x,3));
     364             :   }
     365             : }
     366             : 
     367             : long
     368         364 : FF_f(GEN x)
     369             : {
     370         364 :   switch(x[1])
     371             :   {
     372             :   case t_FF_F2xq:
     373         147 :     return F2x_degree(gel(x,3));
     374             :   default:
     375         217 :     return degpol(gel(x,3));
     376             :   }
     377             : }
     378             : 
     379             : GEN
     380      961646 : FF_to_F2xq(GEN x)
     381             : {
     382      961646 :   switch(x[1])
     383             :   {
     384             :   case t_FF_FpXQ:
     385           0 :     return ZX_to_F2x(gel(x,2));
     386             :   case t_FF_F2xq:
     387      961646 :     return zv_copy(gel(x,2));
     388             :   default:
     389           0 :     return Flx_to_F2x(gel(x,2));
     390             :   }
     391             : }
     392             : 
     393             : GEN
     394           0 : FF_to_F2xq_i(GEN x)
     395             : {
     396           0 :   switch(x[1])
     397             :   {
     398             :   case t_FF_FpXQ:
     399           0 :     return ZX_to_F2x(gel(x,2));
     400             :   case t_FF_F2xq:
     401           0 :     return gel(x,2);
     402             :   default:
     403           0 :     return Flx_to_F2x(gel(x,2));
     404             :   }
     405             : }
     406             : 
     407             : GEN
     408      786758 : FF_to_Flxq(GEN x)
     409             : {
     410      786758 :   switch(x[1])
     411             :   {
     412             :   case t_FF_FpXQ:
     413           0 :     return ZX_to_Flx(gel(x,2),itou(gel(x,4)));
     414             :   case t_FF_F2xq:
     415           0 :     return F2x_to_Flx(gel(x,2));
     416             :   default:
     417      786758 :     return zv_copy(gel(x,2));
     418             :   }
     419             : }
     420             : 
     421             : GEN
     422           0 : FF_to_Flxq_i(GEN x)
     423             : {
     424           0 :   switch(x[1])
     425             :   {
     426             :   case t_FF_FpXQ:
     427           0 :     return ZX_to_Flx(gel(x,2),itou(gel(x,4)));
     428             :   case t_FF_F2xq:
     429           0 :     return F2x_to_Flx(gel(x,2));
     430             :   default:
     431           0 :     return gel(x,2);
     432             :   }
     433             : }
     434             : 
     435             : GEN
     436       19663 : FF_to_FpXQ(GEN x)
     437             : {
     438       19663 :   switch(x[1])
     439             :   {
     440             :   case t_FF_FpXQ:
     441       18180 :     return ZX_copy(gel(x,2));
     442             :   case t_FF_F2xq:
     443          63 :     return F2x_to_ZX(gel(x,2));
     444             :   default:
     445        1420 :     return Flx_to_ZX(gel(x,2));
     446             :   }
     447             : }
     448             : 
     449             : GEN
     450      170149 : FF_to_FpXQ_i(GEN x)
     451             : {
     452      170149 :   switch(x[1])
     453             :   {
     454             :   case t_FF_FpXQ:
     455         764 :     return gel(x,2);
     456             :   case t_FF_F2xq:
     457        1169 :     return F2x_to_ZX(gel(x,2));
     458             :   default:
     459      168216 :     return Flx_to_ZX(gel(x,2));
     460             :   }
     461             : }
     462             : 
     463             : GEN
     464     8855204 : FF_add(GEN x, GEN y)
     465             : {
     466             :   ulong pp;
     467     8855204 :   GEN r, T, p, z=_initFF(x,&T,&p,&pp);
     468     8855203 :   _checkFF(x,y,"+");
     469     8855204 :   switch(x[1])
     470             :   {
     471             :   case t_FF_FpXQ:
     472      296646 :     r=FpX_add(gel(x,2),gel(y,2),p);
     473      296646 :     break;
     474             :   case t_FF_F2xq:
     475     1010184 :     r=F2x_add(gel(x,2),gel(y,2));
     476     1010183 :     break;
     477             :   default:
     478     7548374 :     r=Flx_add(gel(x,2),gel(y,2),pp);
     479             :   }
     480     8855203 :   return _mkFF(x,z,r);
     481             : }
     482             : GEN
     483      281382 : FF_sub(GEN x, GEN y)
     484             : {
     485             :   ulong pp;
     486      281382 :   GEN r, T, p, z=_initFF(x,&T,&p,&pp);
     487      281383 :   _checkFF(x,y,"+");
     488      281383 :   switch(x[1])
     489             :   {
     490             :   case t_FF_FpXQ:
     491        2682 :     r=FpX_sub(gel(x,2),gel(y,2),p);
     492        2682 :     break;
     493             :   case t_FF_F2xq:
     494      154790 :     r=F2x_add(gel(x,2),gel(y,2));
     495      154790 :     break;
     496             :   default:
     497      123911 :     r=Flx_sub(gel(x,2),gel(y,2),pp);
     498             :   }
     499      281383 :   return _mkFF(x,z,r);
     500             : }
     501             : 
     502             : GEN
     503     1654722 : FF_Z_add(GEN x, GEN y)
     504             : {
     505             :   ulong pp;
     506     1654722 :   GEN r, T, p, z=_initFF(x,&T,&p,&pp);
     507     1654722 :   switch(x[1])
     508             :   {
     509             :   case t_FF_FpXQ:
     510             :     {
     511       10291 :       pari_sp av=avma;
     512       10291 :       r=gerepileupto(av,FpX_Fp_add(gel(x,2),modii(y,p),p));
     513       10291 :       break;
     514             :     }
     515             :   case t_FF_F2xq:
     516      791300 :     r=mpodd(y)?F2x_1_add(gel(x,2)):vecsmall_copy(gel(x,2));
     517      791298 :     break;
     518             :   default:
     519      853131 :     r=Flx_Fl_add(gel(x,2),umodiu(y,pp),pp);
     520             :   }
     521     1654720 :   return _mkFF(x,z,r);
     522             : }
     523             : 
     524             : GEN
     525        1344 : FF_Q_add(GEN x, GEN y)
     526             : {
     527             :   ulong pp;
     528        1344 :   GEN r, T, p, z=_initFF(x,&T,&p,&pp);
     529        1344 :   switch(x[1])
     530             :   {
     531             :   case t_FF_FpXQ:
     532             :     {
     533           1 :       pari_sp av=avma;
     534           1 :       r=gerepileupto(av,FpX_Fp_add(gel(x,2),Rg_to_Fp(y,p),p));
     535           1 :       break;
     536             :     }
     537             :   case t_FF_F2xq:
     538           7 :     r=Rg_to_Fl(y,pp)?F2x_1_add(gel(x,2)):vecsmall_copy(gel(x,2));
     539           7 :     break;
     540             :   default:
     541        1336 :     r=Flx_Fl_add(gel(x,2),Rg_to_Fl(y,pp),pp);
     542             :   }
     543        1344 :   return _mkFF(x,z,r);
     544             : }
     545             : 
     546             : GEN
     547       34329 : FF_neg(GEN x)
     548             : {
     549             :   ulong pp;
     550       34329 :   GEN r, T, p, z=_initFF(x,&T,&p,&pp);
     551       34329 :   switch(x[1])
     552             :   {
     553             :   case t_FF_FpXQ:
     554        1354 :     r=FpX_neg(gel(x,2),p);
     555        1354 :     break;
     556             :   case t_FF_F2xq:
     557       12804 :     r=vecsmall_copy(gel(x,2));
     558       12804 :     break;
     559             :   default:
     560       20171 :     r=Flx_neg(gel(x,2),pp);
     561             :   }
     562       34329 :   return _mkFF(x,z,r);
     563             : }
     564             : 
     565             : GEN
     566       85036 : FF_neg_i(GEN x)
     567             : {
     568             :   ulong pp;
     569       85036 :   GEN r, T, p, z=_initFF(x,&T,&p,&pp);
     570       85036 :   switch(x[1])
     571             :   {
     572             :   case t_FF_FpXQ:
     573        1748 :     r=FpX_neg(gel(x,2),p);
     574        1748 :     break;
     575             :   case t_FF_F2xq:
     576        8477 :     r=gel(x,2);
     577        8477 :     break;
     578             :   default:
     579       74811 :     r=Flx_neg(gel(x,2),pp);
     580             :   }
     581       85036 :   return _mkFF_i(x,z,r);
     582             : }
     583             : 
     584             : GEN
     585        1687 : FF_map(GEN m, GEN x)
     586             : {
     587             :   ulong pp;
     588        1687 :   GEN r, T, p, z=_initFF(m,&T,&p,&pp);
     589        1687 :   switch(m[1])
     590             :   {
     591             :   case t_FF_FpXQ:
     592         560 :     r=FpX_FpXQ_eval(gel(x,2),gel(m,2),T,p);
     593         560 :     break;
     594             :   case t_FF_F2xq:
     595         672 :     r=F2x_F2xq_eval(gel(x,2),gel(m,2),T);
     596         672 :     break;
     597             :   default:
     598         455 :     r=Flx_Flxq_eval(gel(x,2),gel(m,2),T,pp);
     599             :   }
     600        1687 :   return _mkFF(m,z,r);
     601             : }
     602             : 
     603             : GEN
     604          42 : FF_Frobenius(GEN x, long e)
     605             : {
     606             :   ulong pp;
     607          42 :   GEN r, T, p, z=_initFF(x,&T,&p,&pp);
     608          42 :   ulong n = umodsu(e, FF_f(x));
     609          42 :   if (n==0) return gcopy(x);
     610          42 :   switch(x[1])
     611             :   {
     612             :   case t_FF_FpXQ:
     613          14 :     r=FpXQ_pow(gel(x,2),p,T,p);
     614          14 :     if (n>1) r=FpXQ_autpow(r,n,T,p);
     615          14 :     break;
     616             :   case t_FF_F2xq:
     617          14 :     r=F2xq_sqr(gel(x,2),T);
     618          14 :     if (n>1) r=F2xq_autpow(r,n,T);
     619          14 :     break;
     620             :   default:
     621          14 :     r=Flxq_powu(gel(x,2),pp,T,pp);
     622          14 :     if (n>1) r=Flxq_autpow(r,n,T,pp);
     623             :   }
     624          42 :   return _mkFF(x,z,r);
     625             : }
     626             : 
     627             : GEN
     628         966 : FFX_preimage(GEN x, GEN F, GEN y)
     629             : {
     630             :   GEN r, T, p, z;
     631             :   ulong pp;
     632         966 :   if (FF_equal0(x)) return FF_zero(y);
     633         875 :   z=_initFF(y,&T,&p,&pp);
     634         875 :   F = FFX_to_raw(F, y);
     635         875 :   switch(y[1])
     636             :   {
     637             :   case t_FF_FpXQ:
     638         322 :     r = FpXQX_rem(gel(x,2), F, T, p);
     639         322 :     break;
     640             :   case t_FF_F2xq:
     641         322 :     r = F2xqX_rem(F2x_to_F2xX(gel(x,2),T[1]), F, T);
     642         322 :     break;
     643             :   default:
     644         231 :     r = FlxqX_rem(Flx_to_FlxX(gel(x,2),T[1]), F, T, pp);
     645             :   }
     646         875 :   if (degpol(r) > 0) return NULL;
     647         791 :   r = (y[1] == t_FF_FpXQ)? Fq_to_FpXQ(gel(r,2),T, p): gel(r,2);
     648         791 :   return _mkFF(y,z,r);
     649             : }
     650             : 
     651             : GEN
     652     9474240 : FF_mul(GEN x, GEN y)
     653             : {
     654             :   ulong pp;
     655     9474240 :   GEN r, T, p, z=_initFF(x,&T,&p,&pp);
     656     9474241 :   pari_sp av=avma;
     657     9474241 :   _checkFF(x,y,"*");
     658     9474240 :   switch(x[1])
     659             :   {
     660             :   case t_FF_FpXQ:
     661      301708 :     r=FpXQ_mul(gel(x,2),gel(y,2),T,p);
     662      301708 :     break;
     663             :   case t_FF_F2xq:
     664      730736 :     r=F2xq_mul(gel(x,2),gel(y,2),T);
     665      730737 :     break;
     666             :   default:
     667     8441796 :     r=Flxq_mul(gel(x,2),gel(y,2),T,pp);
     668             :   }
     669     9474241 :   return _mkFF(x,z,gerepileupto(av, r));
     670             : }
     671             : 
     672             : GEN
     673     1843643 : FF_Z_mul(GEN x, GEN y)
     674             : {
     675             :   ulong pp;
     676     1843643 :   GEN r, T, p, A = gel(x,2), z=_initFF(x,&T,&p,&pp);
     677     1843642 :   switch(x[1])
     678             :   {
     679             :   case t_FF_FpXQ: /* modii(y,p) left on stack for efficiency */
     680       20831 :     r = FpX_Fp_mul(A, modii(y,p),p);
     681       20831 :     break;
     682             :   case t_FF_F2xq:
     683      984532 :     r = mpodd(y)? vecsmall_copy(A): zero_Flx(A[1]);
     684      984533 :     break;
     685             :   default:
     686      838279 :     r = Flx_Fl_mul(A, umodiu(y,pp), pp);
     687             :   }
     688     1843643 :   return _mkFF(x,z,r);
     689             : }
     690             : 
     691             : GEN
     692        4823 : FF_Z_Z_muldiv(GEN x, GEN a, GEN b)
     693             : {
     694             :   ulong pp;
     695        4823 :   GEN r, T, p, A = gel(x,2), z=_initFF(x,&T,&p,&pp);
     696        4823 :   switch(x[1])
     697             :   {
     698             :   case t_FF_FpXQ: /* Fp_div(a,b,p) left on stack for efficiency */
     699          93 :     r = FpX_Fp_mul(A, Fp_div(a,b,p), p);
     700          93 :     break;
     701             :   case t_FF_F2xq:
     702         763 :     if (!mpodd(b)) pari_err_INV("FF_Z_Z_muldiv", b);
     703         756 :     r = mpodd(a)? vecsmall_copy(A): zero_Flx(A[1]);
     704         756 :     break;
     705             :   default:
     706        3967 :     r = Flx_Fl_mul(A, Fl_div(umodiu(a,pp),umodiu(b,pp),pp),pp);
     707             :   }
     708        4809 :   return _mkFF(x,z,r);
     709             : }
     710             : 
     711             : GEN
     712     2777187 : FF_sqr(GEN x)
     713             : {
     714             :   ulong pp;
     715     2777187 :   GEN r, T, p, z=_initFF(x,&T,&p,&pp);
     716     2777188 :   switch(x[1])
     717             :   {
     718             :   case t_FF_FpXQ:
     719             :     {
     720       14088 :       pari_sp av=avma;
     721       14088 :       r=gerepileupto(av,FpXQ_sqr(gel(x,2),T,p));
     722       14088 :       break;
     723             :     }
     724             :   case t_FF_F2xq:
     725      421810 :     r=F2xq_sqr(gel(x,2),T);
     726      421808 :     break;
     727             :   default:
     728     2341290 :     r=Flxq_sqr(gel(x,2),T,pp);
     729             :   }
     730     2777186 :   return _mkFF(x,z,r);
     731             : }
     732             : 
     733             : GEN
     734      215724 : FF_mul2n(GEN x, long n)
     735             : {
     736             :   ulong pp;
     737      215724 :   GEN r, T, p, A = gel(x,2), z=_initFF(x,&T,&p,&pp);
     738      215724 :   switch(x[1])
     739             :   {
     740             :   case t_FF_FpXQ:
     741             :     {
     742             :       GEN p1; /* left on stack for efficiency */
     743        2944 :       if (n>0) p1=remii(int2n(n),p);
     744          22 :       else p1=Fp_inv(remii(int2n(-n),p),p);
     745        2944 :       r = FpX_Fp_mul(A, p1, p);
     746             :     }
     747        2944 :     break;
     748             :   case t_FF_F2xq:
     749       90487 :     if (n<0) pari_err_INV("FF_mul2n", gen_2);
     750       90487 :     r = n==0? vecsmall_copy(A): zero_Flx(A[1]);
     751       90487 :     break;
     752             :   default:
     753             :     {
     754             :       ulong l1;
     755      122293 :       if (n>0) l1 = umodiu(int2n(n),pp);
     756         384 :       else l1 = Fl_inv(umodiu(int2n(-n),pp),pp);
     757      122293 :       r = Flx_Fl_mul(A,l1,pp);
     758             :     }
     759             :   }
     760      215724 :   return _mkFF(x,z,r);
     761             : }
     762             : 
     763             : GEN
     764        7754 : FF_inv(GEN x)
     765             : {
     766             :   ulong pp;
     767        7754 :   GEN r, T, p, z=_initFF(x,&T,&p,&pp);
     768        7753 :   pari_sp av=avma;
     769        7753 :   switch(x[1])
     770             :   {
     771             :   case t_FF_FpXQ:
     772         428 :     r=gerepileupto(av,FpXQ_inv(gel(x,2),T,p));
     773         428 :     break;
     774             :   case t_FF_F2xq:
     775        6822 :     r=F2xq_inv(gel(x,2),T);
     776        6823 :     break;
     777             :   default:
     778         503 :     r=Flxq_inv(gel(x,2),T,pp);
     779             :   }
     780        7754 :   return _mkFF(x,z,r);
     781             : }
     782             : 
     783             : GEN
     784      121037 : FF_div(GEN x, GEN y)
     785             : {
     786             :   ulong pp;
     787      121037 :   GEN r, T, p, z=_initFF(x,&T,&p,&pp);
     788      121037 :   pari_sp av=avma;
     789      121037 :   _checkFF(x,y,"/");
     790      121037 :   switch(x[1])
     791             :   {
     792             :   case t_FF_FpXQ:
     793        1676 :     r=gerepileupto(av,FpXQ_div(gel(x,2),gel(y,2),T,p));
     794        1676 :     break;
     795             :   case t_FF_F2xq:
     796       55706 :     r=gerepileupto(av,F2xq_div(gel(x,2),gel(y,2),T));
     797       55636 :     break;
     798             :   default:
     799       63655 :     r=gerepileupto(av,Flxq_div(gel(x,2),gel(y,2),T,pp));
     800             :   }
     801      120939 :   return _mkFF(x,z,r);
     802             : }
     803             : 
     804             : GEN
     805         259 : Z_FF_div(GEN n, GEN x)
     806             : {
     807             :   ulong pp;
     808         259 :   GEN r, T, p, A = gel(x,2), z=_initFF(x,&T,&p,&pp);
     809         259 :   pari_sp av=avma;
     810         259 :   switch(x[1])
     811             :   {
     812             :   case t_FF_FpXQ:
     813           8 :     r = gerepileupto(av,FpX_Fp_mul(FpXQ_inv(A,T,p),modii(n,p),p));
     814           8 :     break;
     815             :   case t_FF_F2xq:
     816          56 :     r = F2xq_inv(A,T); /*Check for division by 0*/
     817          56 :     if(!mpodd(n)) { avma = av; r = zero_Flx(A[1]); }
     818          56 :     break;
     819             :   default:
     820         195 :     r = gerepileupto(av, Flx_Fl_mul(Flxq_inv(A,T,pp),umodiu(n,pp),pp));
     821             :   }
     822         259 :   return _mkFF(x,z,r);
     823             : }
     824             : 
     825             : GEN
     826         217 : FF_sqrtn(GEN x, GEN n, GEN *zetan)
     827             : {
     828             :   ulong pp;
     829         217 :   GEN r, T, p, y=_initFF(x,&T,&p,&pp);
     830         217 :   switch (x[1])
     831             :   {
     832             :   case t_FF_FpXQ:
     833          66 :     r=FpXQ_sqrtn(gel(x,2),n,T,p,zetan);
     834          59 :     break;
     835             :   case t_FF_F2xq:
     836          28 :     r=F2xq_sqrtn(gel(x,2),n,T,zetan);
     837          21 :     break;
     838             :   default:
     839         123 :     r=Flxq_sqrtn(gel(x,2),n,T,pp,zetan);
     840             :   }
     841         196 :   if (!r) pari_err_SQRTN("FF_sqrtn",x);
     842         196 :   (void)_mkFF(x, y, r);
     843         196 :   if (zetan)
     844             :   {
     845         140 :     GEN z = cgetg(lg(y),t_FFELT);
     846         140 :     *zetan=_mkFF(x, z, *zetan);
     847             :   }
     848         196 :   return y;
     849             : }
     850             : 
     851             : GEN
     852          70 : FF_sqrt(GEN x)
     853             : {
     854             :   ulong pp;
     855          70 :   GEN r, T, p, y=_initFF(x,&T,&p,&pp);
     856          70 :   switch (x[1])
     857             :   {
     858             :   case t_FF_FpXQ:
     859           1 :     r = FpXQ_sqrt(gel(x,2),T,p);
     860           1 :     break;
     861             :   case t_FF_F2xq:
     862          49 :     r = F2xq_sqrt(gel(x,2),T);
     863          49 :     break;
     864             :   default:
     865          20 :     r = Flxq_sqrt(gel(x,2),T,pp);
     866             :   }
     867          70 :   if (!r) pari_err_SQRTN("FF_sqrt",x);
     868          70 :   return _mkFF(x, y, r);
     869             : }
     870             : 
     871             : long
     872           7 : FF_issquare(GEN x)
     873             : {
     874             :   GEN T, p;
     875             :   ulong pp;
     876           7 :   _getFF(x, &T, &p, &pp);
     877           7 :   switch(x[1])
     878             :   {
     879             :   case t_FF_FpXQ:
     880           0 :     return FpXQ_issquare(gel(x,2), T, p);
     881             :   case t_FF_F2xq:
     882           7 :     return 1;
     883             :   default: /* case t_FF_Flxq: */
     884           0 :     return Flxq_issquare(gel(x,2), T, pp);
     885             :   }
     886             : }
     887             : 
     888             : long
     889         182 : FF_issquareall(GEN x, GEN *pt)
     890             : {
     891         182 :   if (!pt) return FF_issquare(x);
     892         175 :   return FF_ispower(x, gen_2, pt);
     893             : }
     894             : 
     895             : long
     896         203 : FF_ispower(GEN x, GEN K, GEN *pt)
     897             : {
     898             :   ulong pp;
     899             :   GEN r, T, p;
     900         203 :   pari_sp av = avma;
     901             : 
     902         203 :   if (FF_equal0(x)) { if (pt) *pt = gcopy(x); return 1; }
     903         203 :   _getFF(x, &T, &p, &pp);
     904         203 :   if (pt) *pt = cgetg(5,t_FFELT);
     905         203 :   switch(x[1])
     906             :   {
     907             :   case t_FF_FpXQ:
     908          71 :     r = FpXQ_sqrtn(gel(x,2),K,T,p,NULL);
     909          71 :     break;
     910             :   case t_FF_F2xq:
     911          42 :     r = F2xq_sqrtn(gel(x,2),K,T,NULL);
     912          42 :     break;
     913             :   default: /* case t_FF_Flxq: */
     914          90 :     r = Flxq_sqrtn(gel(x,2),K,T,pp,NULL);
     915          90 :     break;
     916             :   }
     917         203 :   if (!r) { avma = av; return 0; }
     918         112 :   if (pt) { (void)_mkFF(x,*pt,r); }
     919         112 :   return 1;
     920             : }
     921             : 
     922             : GEN
     923          79 : FF_pow(GEN x, GEN n)
     924             : {
     925             :   ulong pp;
     926          79 :   GEN r, T, p, z=_initFF(x,&T,&p,&pp);
     927          79 :   switch(x[1])
     928             :    {
     929             :    case t_FF_FpXQ:
     930          17 :      r = FpXQ_pow(gel(x,2), n, T, p);
     931          17 :      break;
     932             :    case t_FF_F2xq:
     933           3 :      r = F2xq_pow(gel(x,2), n, T);
     934           3 :      break;
     935             :    default:
     936          59 :      r = Flxq_pow(gel(x,2), n, T, pp);
     937             :    }
     938          79 :   return _mkFF(x,z,r);
     939             : }
     940             : 
     941             : GEN
     942          28 : FF_norm(GEN x)
     943             : {
     944             :   ulong pp;
     945             :   GEN T,p;
     946          28 :   _getFF(x,&T,&p,&pp);
     947          28 :   switch (x[1])
     948             :   {
     949             :   case t_FF_FpXQ:
     950           1 :     return FpXQ_norm(gel(x,2),T,p);
     951             :   case t_FF_F2xq:
     952           7 :     return lgpol(gel(x,2))?gen_1:gen_0;
     953             :   default:
     954          20 :     return utoi(Flxq_norm(gel(x,2),T,pp));
     955             :   }
     956             : }
     957             : 
     958             : GEN
     959          28 : FF_trace(GEN x)
     960             : {
     961             :   ulong pp;
     962             :   GEN T,p;
     963          28 :   _getFF(x,&T,&p,&pp);
     964          28 :   switch(x[1])
     965             :   {
     966             :   case t_FF_FpXQ:
     967           1 :     return FpXQ_trace(gel(x,2),T,p);
     968             :   case t_FF_F2xq:
     969           7 :     return F2xq_trace(gel(x,2),T)?gen_1:gen_0;
     970             :   default:
     971          20 :     return utoi(Flxq_trace(gel(x,2),T,pp));
     972             :   }
     973             : }
     974             : 
     975             : GEN
     976          28 : FF_conjvec(GEN x)
     977             : {
     978             :   ulong pp;
     979             :   GEN r,T,p,v;
     980             :   long i,l;
     981             :   pari_sp av;
     982          28 :   _getFF(x,&T,&p,&pp);
     983          28 :   av = avma;
     984          28 :   switch(x[1])
     985             :   {
     986             :   case t_FF_FpXQ:
     987           1 :     v = FpXQ_conjvec(gel(x,2), T, p);
     988           1 :     break;
     989             :   case t_FF_F2xq:
     990           7 :     v = F2xq_conjvec(gel(x,2), T);
     991           7 :     break;
     992             :   default:
     993          20 :     v = Flxq_conjvec(gel(x,2), T, pp);
     994             :   }
     995          28 :   l = lg(v); r = cgetg(l, t_COL);
     996         252 :   for(i=1; i<l; i++)
     997         224 :     gel(r,i) = mkFF_i(x, gel(v,i));
     998          28 :   return gerepilecopy(av, r);
     999             : }
    1000             : 
    1001             : GEN
    1002          28 : FF_charpoly(GEN x)
    1003             : {
    1004             :   ulong pp;
    1005             :   GEN T,p;
    1006          28 :   pari_sp av=avma;
    1007          28 :   _getFF(x,&T,&p,&pp);
    1008          28 :   switch(x[1])
    1009             :   {
    1010             :   case t_FF_FpXQ:
    1011           1 :     return gerepileupto(av,FpXQ_charpoly(gel(x,2), T, p));
    1012             :   case t_FF_F2xq:
    1013           7 :     return gerepileupto(av,Flx_to_ZX(Flxq_charpoly(F2x_to_Flx(gel(x,2)),
    1014             :                                                    F2x_to_Flx(T), 2UL)));
    1015             :   default:
    1016          20 :     return gerepileupto(av,Flx_to_ZX(Flxq_charpoly(gel(x,2), T, pp)));
    1017             :   }
    1018             : }
    1019             : 
    1020             : GEN
    1021         154 : FF_minpoly(GEN x)
    1022             : {
    1023             :   ulong pp;
    1024             :   GEN T,p;
    1025         154 :   pari_sp av=avma;
    1026         154 :   _getFF(x,&T,&p,&pp);
    1027         154 :   switch(x[1])
    1028             :   {
    1029             :   case t_FF_FpXQ:
    1030          43 :     return gerepileupto(av,FpXQ_minpoly(gel(x,2), T, p));
    1031             :   case t_FF_F2xq:
    1032          49 :     return gerepileupto(av,Flx_to_ZX(Flxq_minpoly(F2x_to_Flx(gel(x,2)),
    1033             :                                                   F2x_to_Flx(T), 2UL)));
    1034             :   default:
    1035          62 :     return gerepileupto(av,Flx_to_ZX(Flxq_minpoly(gel(x,2), T, pp)));
    1036             :   }
    1037             : }
    1038             : 
    1039             : GEN
    1040         175 : FF_log(GEN x, GEN g, GEN ord)
    1041             : {
    1042         175 :   pari_sp av=avma;
    1043             :   ulong pp;
    1044             :   GEN r, T, p;
    1045         175 :   _getFF(x,&T,&p,&pp);
    1046         175 :   _checkFF(x,g,"log");
    1047         175 :   switch(x[1])
    1048             :   {
    1049             :   case t_FF_FpXQ:
    1050          44 :     if (!ord) ord = factor_pn_1(p,degpol(T));
    1051          44 :     r = FpXQ_log(gel(x,2), gel(g,2), ord, T, p);
    1052          16 :     break;
    1053             :   case t_FF_F2xq:
    1054          28 :     if (!ord) ord = factor_pn_1(gen_2,F2x_degree(T));
    1055          28 :     r = F2xq_log(gel(x,2), gel(g,2), ord, T);
    1056          28 :     break;
    1057             :   default:
    1058         103 :     if (!ord) ord = factor_pn_1(p,degpol(T));
    1059         103 :     r = Flxq_log(gel(x,2), gel(g,2), ord, T, pp);
    1060             :   }
    1061         147 :   return gerepileupto(av, r);
    1062             : }
    1063             : 
    1064             : GEN
    1065          28 : FF_order(GEN x, GEN o)
    1066             : {
    1067          28 :   pari_sp av=avma;
    1068             :   ulong pp;
    1069             :   GEN r, T,p;
    1070          28 :   _getFF(x,&T,&p,&pp);
    1071          28 :   switch(x[1])
    1072             :   {
    1073             :   case t_FF_FpXQ:
    1074           1 :     if (!o) o = factor_pn_1(p,degpol(T));
    1075           1 :     r = FpXQ_order(gel(x,2), o, T, p);
    1076           1 :     break;
    1077             :   case t_FF_F2xq:
    1078           7 :     if (!o) o = factor_pn_1(gen_2,F2x_degree(T));
    1079           7 :     r = F2xq_order(gel(x,2), o, T);
    1080           7 :     break;
    1081             :   default:
    1082          20 :     if (!o) o = factor_pn_1(p,degpol(T));
    1083          20 :     r = Flxq_order(gel(x,2), o, T, pp);
    1084             :   }
    1085          28 :   if (!o) r = gerepileuptoint(av,r);
    1086          28 :   return r;
    1087             : }
    1088             : 
    1089             : GEN
    1090         385 : FF_primroot(GEN x, GEN *o)
    1091             : {
    1092             :   ulong pp;
    1093         385 :   GEN r,T,p,z=_initFF(x,&T,&p,&pp);
    1094         385 :   switch(x[1])
    1095             :   {
    1096             :   case t_FF_FpXQ:
    1097          30 :     r = gener_FpXQ(T, p, o);
    1098          30 :     break;
    1099             :   case t_FF_F2xq:
    1100         154 :     r = gener_F2xq(T, o);
    1101         154 :     break;
    1102             :   default:
    1103         201 :     r = gener_Flxq(T, pp, o);
    1104             :   }
    1105         385 :   return _mkFF(x,z,r);
    1106             : }
    1107             : 
    1108             : static GEN
    1109      510755 : to_FFE(GEN P, GEN fg)
    1110             : {
    1111      510755 :   if(ell_is_inf(P))
    1112      231833 :     return ellinf();
    1113             :   else
    1114      278922 :     retmkvec2(mkFF_i(fg,gel(P,1)), mkFF_i(fg,gel(P,2)));
    1115             : }
    1116             : 
    1117             : static GEN
    1118       16772 : to_FFE_vec(GEN x, GEN ff)
    1119             : {
    1120       16772 :   long i, lx = lg(x);
    1121       16772 :   for (i=1; i<lx; i++) gel(x,i) = to_FFE(gel(x,i), ff);
    1122       16772 :   return x;
    1123             : }
    1124             : 
    1125             : static GEN
    1126        1710 : FpXQ_ell_to_a4a6(GEN E, GEN T, GEN p)
    1127             : {
    1128             :   GEN a1, a3, b2, c4, c6;
    1129        1710 :   a1 = Rg_to_FpXQ(ell_get_a1(E),T,p);
    1130        1710 :   a3 = Rg_to_FpXQ(ell_get_a3(E),T,p);
    1131        1710 :   b2 = Rg_to_FpXQ(ell_get_b2(E),T,p);
    1132        1710 :   c4 = Rg_to_FpXQ(ell_get_c4(E),T,p);
    1133        1710 :   c6 = Rg_to_FpXQ(ell_get_c6(E),T,p);
    1134        1710 :   retmkvec3(FpX_neg(FpX_mulu(c4, 27, p), p), FpX_neg(FpX_mulu(c6, 54, p), p),
    1135             :             mkvec4(Z_to_FpX(utoi(6),p,varn(T)),FpX_mulu(b2,3,p),
    1136             :                    FpX_mulu(a1,3,p),FpX_mulu(a3,108,p)));
    1137             : }
    1138             : 
    1139             : static GEN
    1140       23541 : F3xq_ell_to_a4a6(GEN E, GEN T)
    1141             : {
    1142             :   GEN a1, a3, b2, b4, b6;
    1143       23541 :   a1 = Rg_to_Flxq(ell_get_a1(E),T,3);
    1144       23541 :   a3 = Rg_to_Flxq(ell_get_a3(E),T,3);
    1145       23541 :   b2 = Rg_to_Flxq(ell_get_b2(E),T,3);
    1146       23541 :   b4 = Rg_to_Flxq(ell_get_b4(E),T,3);
    1147       23541 :   b6 = Rg_to_Flxq(ell_get_b6(E),T,3);
    1148       23541 :   if(lgpol(b2)) /* ordinary case */
    1149             :   {
    1150       17108 :     GEN b4b2 = Flxq_div(b4,b2,T,3);
    1151       17108 :     GEN a6 = Flx_sub(b6,Flxq_mul(b4b2,Flx_add(b4,Flxq_sqr(b4b2,T,3),3),T,3),3);
    1152       17108 :     retmkvec3(mkvec(b2), a6,
    1153             :        mkvec4(Fl_to_Flx(1,T[1]),b4b2,Flx_neg(a1,3),Flx_neg(a3,3)));
    1154             :   }
    1155             :   else /* super-singular case */
    1156        6433 :     retmkvec3(Flx_neg(b4, 3), b6,
    1157             :        mkvec4(Fl_to_Flx(1,T[1]),zero_Flx(T[1]), Flx_neg(a1,3), Flx_neg(a3,3)));
    1158             : }
    1159             : 
    1160             : static GEN
    1161       65602 : Flxq_ell_to_a4a6(GEN E, GEN T, ulong p)
    1162             : {
    1163             :   GEN a1, a3, b2, c4, c6;
    1164       65602 :   if(p==3) return F3xq_ell_to_a4a6(E, T);
    1165       42061 :   a1 = Rg_to_Flxq(ell_get_a1(E),T,p);
    1166       42061 :   a3 = Rg_to_Flxq(ell_get_a3(E),T,p);
    1167       42061 :   b2 = Rg_to_Flxq(ell_get_b2(E),T,p);
    1168       42061 :   c4 = Rg_to_Flxq(ell_get_c4(E),T,p);
    1169       42061 :   c6 = Rg_to_Flxq(ell_get_c6(E),T,p);
    1170       42061 :   retmkvec3(Flx_neg(Flx_mulu(c4, 27, p), p), Flx_neg(Flx_mulu(c6, 54, p), p),
    1171             :             mkvec4(Fl_to_Flx(6%p,T[1]),Flx_mulu(b2,3,p),
    1172             :                    Flx_mulu(a1,3,p),Flx_mulu(a3,108,p)));
    1173             : }
    1174             : 
    1175             : static GEN
    1176       45365 : F2xq_ell_to_a4a6(GEN E, GEN T)
    1177             : {
    1178       45365 :   long v = T[1];
    1179       45365 :   GEN a1 = Rg_to_F2xq(ell_get_a1(E),T);
    1180       45365 :   GEN a2 = Rg_to_F2xq(ell_get_a2(E),T);
    1181       45365 :   GEN a3 = Rg_to_F2xq(ell_get_a3(E),T);
    1182       45365 :   GEN a4 = Rg_to_F2xq(ell_get_a4(E),T);
    1183       45365 :   GEN a6 = Rg_to_F2xq(ell_get_a6(E),T);
    1184       45365 :   if (lgpol(a1))
    1185             :   {
    1186        7481 :     GEN a1i = F2xq_inv(a1,T);
    1187        7481 :     GEN a1i2 = F2xq_sqr(a1i,T);
    1188        7480 :     GEN a1i3 = F2xq_mul(a1i,a1i2,T);
    1189        7481 :     GEN a1i6 = F2xq_sqr(a1i3,T);
    1190        7481 :     GEN d  = F2xq_mul(a3,a1i,T);
    1191        7481 :     GEN dd = F2xq_mul(d,a1i2,T);
    1192        7480 :     GEN e  = F2xq_mul(F2x_add(a4,F2xq_sqr(d,T)),a1i,T);
    1193        7480 :     GEN ee = F2xq_mul(e,a1i3,T);
    1194        7481 :     GEN da2 = F2x_add(a2,d);
    1195        7480 :     GEN d2 = F2xq_mul(da2,a1i2,T);
    1196        7480 :     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);
    1197        7480 :     retmkvec3(d2, d6, mkvec4(a1i,dd,pol0_F2x(v),ee));
    1198             :   }
    1199       37884 :   else if (lgpol(a3))
    1200             :   {
    1201       37870 :     GEN d4 = F2x_add(F2xq_sqr(a2,T),a4);
    1202       37870 :     GEN d6 = F2x_add(F2xq_mul(a2,a4,T),a6);
    1203       37870 :     GEN a3i = F2xq_inv(a3,T);
    1204       37870 :     retmkvec3(mkvec3(a3,d4,a3i), d6, mkvec4(pol1_F2x(v),a2,pol0_F2x(T[1]),pol0_F2x(T[1])));
    1205             :   }
    1206          14 :   return NULL;
    1207             : }
    1208             : 
    1209             : static GEN
    1210        1743 : FqV_to_FpXQV(GEN x, GEN T)
    1211             : {
    1212        1743 :   pari_sp av = avma;
    1213        1743 :   long v = varn(T), i, s=0, l = lg(x);
    1214        1743 :   GEN y = shallowcopy(x);
    1215        8715 :   for(i=1; i<l; i++)
    1216        6972 :     if (typ(gel(x,i))==t_INT)
    1217             :     {
    1218           0 :       gel(y,i) = scalarpol(gel(x,i), v);
    1219           0 :       s = 1;
    1220             :     }
    1221        1743 :   if (!s) { avma = av; return x; }
    1222           0 :   return y;
    1223             : }
    1224             : 
    1225             : GEN
    1226       94252 : FF_ellcard(GEN E)
    1227             : {
    1228             :   GEN T,p;
    1229             :   ulong pp;
    1230       94252 :   GEN e = ellff_get_a4a6(E);
    1231       94252 :   GEN fg = ellff_get_field(E);
    1232       94252 :   _getFF(fg,&T,&p,&pp);
    1233       94253 :   switch(fg[1])
    1234             :   {
    1235             :   case t_FF_FpXQ:
    1236        1695 :     return FpXQ_ellcard(Fq_to_FpXQ(gel(e,1),T,p), Fq_to_FpXQ(gel(e,2),T,p),T,p);
    1237             :   case t_FF_F2xq:
    1238       45141 :     return F2xq_ellcard(gel(e,1),gel(e,2),T);
    1239             :   default:
    1240       47417 :     return Flxq_ellcard(gel(e,1),gel(e,2),T,pp);
    1241             :   }
    1242             : }
    1243             : 
    1244             : GEN
    1245          14 : FF_ellcard_SEA(GEN E, long smallfact)
    1246             : {
    1247          14 :   pari_sp av = avma;
    1248             :   GEN T,p;
    1249             :   ulong pp;
    1250          14 :   GEN e = ellff_get_a4a6(E), a4, a6, card;
    1251          14 :   GEN fg = ellff_get_field(E);
    1252          14 :   _getFF(fg,&T,&p,&pp);
    1253          14 :   switch(fg[1])
    1254             :   {
    1255             :   case t_FF_FpXQ:
    1256           0 :     a4 = Fq_to_FpXQ(gel(e,1), T, p);
    1257           0 :     a6 = Fq_to_FpXQ(gel(e,2), T, p);
    1258           0 :     card = Fq_ellcard_SEA(a4, a6, powiu(p,degpol(T)), T,p, smallfact);
    1259           0 :     break;
    1260             :   case t_FF_F2xq:
    1261           0 :     pari_err_IMPL("SEA for char 2");
    1262             :   default:
    1263          14 :     a4 = Flx_to_ZX(gel(e,1));
    1264          14 :     a6 = Flx_to_ZX(gel(e,2));
    1265          14 :     card = Fq_ellcard_SEA(a4, a6, powuu(pp,degpol(T)), Flx_to_ZX(T), p, smallfact);
    1266             :   }
    1267          14 :   return gerepileuptoint(av, card);
    1268             : }
    1269             : 
    1270             : /* return [G,m] */
    1271             : GEN
    1272       18081 : FF_ellgroup(GEN E)
    1273             : {
    1274             :   GEN T, p, e, fg, N, G, m;
    1275             :   ulong pp;
    1276             : 
    1277       18081 :   N = ellff_get_card(E);
    1278       18081 :   e = ellff_get_a4a6(E);
    1279       18081 :   fg = ellff_get_field(E);
    1280       18081 :   _getFF(fg,&T,&p,&pp);
    1281       18081 :   switch(fg[1])
    1282             :   {
    1283             :   case t_FF_FpXQ:
    1284          14 :     G = FpXQ_ellgroup(Fq_to_FpXQ(gel(e,1),T,p),Fq_to_FpXQ(gel(e,2),T,p),N,T,p,&m);
    1285          14 :     break;
    1286             :   case t_FF_F2xq:
    1287        3668 :     G = F2xq_ellgroup(gel(e,1),gel(e,2),N,T,&m); break;
    1288             :   default:
    1289       14399 :     G = Flxq_ellgroup(gel(e,1),gel(e,2),N,T,pp,&m); break;
    1290             :   }
    1291       18081 :   return mkvec2(G,m);
    1292             : }
    1293             : 
    1294             : GEN
    1295       16772 : FF_ellgens(GEN E)
    1296             : {
    1297             :   GEN fg, e, Gm, G, m, T, p, F, e3;
    1298             :   ulong pp;
    1299             : 
    1300       16772 :   fg = ellff_get_field(E);
    1301       16772 :   e = ellff_get_a4a6(E);
    1302       16772 :   Gm = ellff_get_group(E); G = gel(Gm,1); m = gel(Gm,2);
    1303       16772 :   _getFF(fg,&T,&p,&pp);
    1304       16772 :   switch(fg[1])
    1305             :   {
    1306             :   case t_FF_FpXQ:
    1307           7 :     e3 = FqV_to_FpXQV(gel(e,3),T);
    1308           7 :     F = FpXQ_ellgens(Fq_to_FpXQ(gel(e,1),T,p),Fq_to_FpXQ(gel(e,2),T,p),e3,G,m,T,p);
    1309           7 :     break;
    1310             :   case t_FF_F2xq:
    1311        3598 :     F = F2xq_ellgens(gel(e,1),gel(e,2),gel(e,3),G,m,T);
    1312        3598 :     break;
    1313             :   default:
    1314       13167 :     F = Flxq_ellgens(gel(e,1),gel(e,2),gel(e,3),G,m,T, pp);
    1315       13167 :     break;
    1316             :   }
    1317       16772 :   return to_FFE_vec(F,fg);
    1318             : }
    1319             : 
    1320             : GEN
    1321         119 : FF_elldata(GEN E, GEN fg)
    1322             : {
    1323             :   GEN T,p,e;
    1324             :   ulong pp;
    1325         119 :   _getFF(fg,&T,&p,&pp);
    1326         119 :   switch(fg[1])
    1327             :   {
    1328             :   case t_FF_FpXQ:
    1329           0 :     e = FpXQ_ell_to_a4a6(E,T,p); break;
    1330             :   case t_FF_F2xq:
    1331          56 :     e = F2xq_ell_to_a4a6(E,T); break;
    1332             :   default:
    1333          63 :     e = Flxq_ell_to_a4a6(E,T,pp); break;
    1334             :   }
    1335         119 :   return mkvec2(fg,e);
    1336             : }
    1337             : 
    1338             : GEN
    1339      112557 : FF_ellinit(GEN E, GEN fg)
    1340             : {
    1341             :   GEN T,p,e;
    1342             :   ulong pp;
    1343             :   long i;
    1344      112557 :   GEN y=E;
    1345      112557 :   _getFF(fg,&T,&p,&pp);
    1346      112558 :   switch(fg[1])
    1347             :   {
    1348             :   case t_FF_FpXQ:
    1349        1710 :     e = FpXQ_ell_to_a4a6(E,T,p);
    1350       22230 :     for(i=1;i<=12;i++)
    1351       20520 :       gel(y,i) = mkFF_i(fg,Rg_to_FpXQ(gel(E,i),T,p));
    1352        1710 :     if (FF_equal0(gel(y,12))) return NULL;
    1353        1710 :     gel(y,i) = mkFF_i(fg,Rg_to_FpXQ(gel(E,i),T,p));
    1354        1710 :     break;
    1355             :   case t_FF_F2xq:
    1356       45309 :     e = F2xq_ell_to_a4a6(E,T);
    1357       45309 :     if (!e) return NULL;
    1358      588835 :     for(i=1;i<=12;i++)
    1359      543540 :       gel(y,i) = mkFF_i(fg,Rg_to_F2xq(gel(E,i),T));
    1360       45295 :     if (FF_equal0(gel(y,12))) return NULL;
    1361       45295 :     gel(y,i) = mkFF_i(fg,Rg_to_F2xq(gel(E,i),T));
    1362       45295 :     break;
    1363             :   default:
    1364       65539 :     e = Flxq_ell_to_a4a6(E,T,pp);
    1365      852007 :     for(i=1;i<=12;i++)
    1366      786468 :       gel(y,i) = mkFF_i(fg,Rg_to_Flxq(gel(E,i),T,pp));
    1367       65539 :     if (FF_equal0(gel(y,12))) return NULL;
    1368       65539 :     gel(y,i) = mkFF_i(fg,Rg_to_Flxq(gel(E,i),T,pp));
    1369       65539 :     break;
    1370             :   }
    1371      112544 :   gel(y,14) = mkvecsmall(t_ELL_Fq);
    1372      112544 :   gel(y,15) = mkvec2(fg,e);
    1373      112544 :   return y;
    1374             : }
    1375             : 
    1376             : GEN
    1377       27188 : FF_elltwist(GEN E)
    1378             : {
    1379       27188 :   pari_sp av = avma;
    1380       27188 :   GEN fg = ellff_get_field(E), e = ellff_get_a4a6(E);
    1381             :   GEN T, p, a, a6, V;
    1382             :   ulong pp;
    1383       27188 :   _getFF(fg,&T,&p,&pp);
    1384       27188 :   switch (fg[1])
    1385             :   {
    1386             :   case t_FF_FpXQ:
    1387         840 :     FpXQ_elltwist(gel(e,1), gel(e,2), T, p, &a, &a6);
    1388         840 :     V = mkvec5(gen_0, gen_0, gen_0, mkFF_i(fg, a), mkFF_i(fg, a6));
    1389         840 :     break;
    1390             :   case t_FF_F2xq:
    1391        3514 :     F2xq_elltwist(gel(e,1), gel(e,2), T, &a, &a6);
    1392        7028 :     V = typ(a)==t_VECSMALL ?
    1393        3514 :           mkvec5(gen_1, mkFF_i(fg, a), gen_0, gen_0, mkFF_i(fg, a6)):
    1394           0 :           mkvec5(gen_0, gen_0, mkFF_i(fg, gel(a,1)), mkFF_i(fg, gel(a,2)), mkFF_i(fg, a6));
    1395        3514 :     break;
    1396             :   default:
    1397       22834 :     Flxq_elltwist(gel(e,1), gel(e,2), T, pp, &a, &a6);
    1398       45668 :     V = typ(a)==t_VECSMALL ?
    1399       30436 :           mkvec5(gen_0, gen_0, gen_0, mkFF_i(fg, a), mkFF_i(fg, a6)):
    1400        7602 :           mkvec5(gen_0, mkFF_i(fg, gel(a,1)), gen_0, gen_0, mkFF_i(fg, a6));
    1401             :   }
    1402       27188 :   return gerepilecopy(av, V);
    1403             : }
    1404             : 
    1405             : static long
    1406           0 : F3x_equalm1(GEN x)
    1407           0 : { return degpol(x)==0 && x[2] == 2; }
    1408             : 
    1409             : GEN
    1410      244034 : FF_ellrandom(GEN E)
    1411             : {
    1412      244034 :   pari_sp av = avma;
    1413      244034 :   GEN fg = ellff_get_field(E), e = ellff_get_a4a6(E), Q;
    1414             :   GEN T,p;
    1415             :   ulong pp;
    1416      244034 :   _getFF(fg,&T,&p,&pp);
    1417      244034 :   switch (fg[1])
    1418             :   {
    1419             :   case t_FF_FpXQ:
    1420         854 :     Q = random_FpXQE(Fq_to_FpXQ(gel(e,1),T,p), Fq_to_FpXQ(gel(e,2),T,p), T, p);
    1421         854 :     Q = FpXQE_changepoint(Q, FqV_to_FpXQV(gel(e,3), T) , T, p);
    1422         854 :     break;
    1423             :   case t_FF_F2xq:
    1424             :     {
    1425      168056 :       long d = F2x_degree(T);
    1426             :       /* if #E(Fq) = 1 return [0] */
    1427      168056 :       if (d<=2 && typ(gel(e,1)) == t_VEC)
    1428             :       { /* over F2 or F4, supersingular */
    1429        1575 :         GEN v = gel(e,1), A6 = gel(e,2), a3 = gel(v,1), A4 = gel(v,2);
    1430        1575 :         if (F2x_equal1(a3) &&
    1431          84 :              ((d==1 && F2x_equal1(A4) && F2x_equal1(A6))
    1432         483 :            || (d==2 && !lgpol(A4)     && F2x_degree(A6)==1))) return ellinf();
    1433             :       }
    1434      168000 :       Q = random_F2xqE(gel(e,1), gel(e,2), T);
    1435      168000 :       Q = F2xqE_changepoint(Q, gel(e,3), T);
    1436      168000 :       break;
    1437             :     }
    1438             :   default:
    1439             :     /* if #E(Fq) = 1 return [0] */
    1440       75124 :     if (pp==3 && degpol(T)==1 && typ(gel(e,1))==t_VECSMALL)
    1441             :     { /* over F3, supersingular */
    1442           0 :       GEN mb4 = gel(e,1), b6 = gel(e,2);
    1443           0 :       if (F3x_equalm1(mb4) && F3x_equalm1(b6)) return ellinf();
    1444             :     }
    1445       75124 :     Q = random_FlxqE(gel(e,1), gel(e,2), T, pp);
    1446       75124 :     Q = FlxqE_changepoint(Q, gel(e,3), T, pp);
    1447             :   }
    1448      243978 :   return gerepilecopy(av, to_FFE(Q, fg));
    1449             : }
    1450             : 
    1451             : GEN
    1452      247450 : FF_ellmul(GEN E, GEN P, GEN n)
    1453             : {
    1454      247450 :   pari_sp av = avma;
    1455      247450 :   GEN fg = ellff_get_field(E), e = ellff_get_a4a6(E), Q;
    1456             :   GEN T,p, Pp, Qp, e3;
    1457             :   ulong pp;
    1458      247450 :   _getFF(fg,&T,&p,&pp);
    1459      247450 :   switch (fg[1])
    1460             :   {
    1461             :   case t_FF_FpXQ:
    1462         854 :     e3 = FqV_to_FpXQV(gel(e,3),T);
    1463         854 :     Pp = FpXQE_changepointinv(RgE_to_FpXQE(P, T, p), e3, T, p);
    1464         854 :     Qp = FpXQE_mul(Pp, n, gel(e,1), T, p);
    1465         854 :     Q = FpXQE_changepoint(Qp, gel(e,3), T, p);
    1466         854 :     break;
    1467             :   case t_FF_F2xq:
    1468      184352 :     Pp = F2xqE_changepointinv(RgE_to_F2xqE(P, T), gel(e,3), T);
    1469      184352 :     Qp = F2xqE_mul(Pp, n, gel(e,1), T);
    1470      184352 :     Q = F2xqE_changepoint(Qp, gel(e,3), T);
    1471      184352 :     break;
    1472             :   default:
    1473       62244 :     Pp = FlxqE_changepointinv(RgE_to_FlxqE(P, T, pp), gel(e,3), T, pp);
    1474       62244 :     Qp = FlxqE_mul(Pp, n, gel(e,1), T, pp);
    1475       62244 :     Q = FlxqE_changepoint(Qp, gel(e,3), T, pp);
    1476             :   }
    1477      247450 :   return gerepilecopy(av, to_FFE(Q, fg));
    1478             : }
    1479             : 
    1480             : GEN
    1481        1750 : FF_ellorder(GEN E, GEN P, GEN o)
    1482             : {
    1483        1750 :   pari_sp av = avma;
    1484        1750 :   GEN fg = ellff_get_field(E), e = ellff_get_a4a6(E);
    1485             :   GEN r,T,p,Pp,e3;
    1486             :   ulong pp;
    1487        1750 :   _getFF(fg,&T,&p,&pp);
    1488        1750 :   switch (fg[1])
    1489             :   {
    1490             :   case t_FF_FpXQ:
    1491          14 :     e3 = FqV_to_FpXQV(gel(e,3),T);
    1492          14 :     Pp = FpXQE_changepointinv(RgE_to_FpXQE(P,T,p), e3, T, p);
    1493          14 :     r = FpXQE_order(Pp, o, gel(e,1), T, p);
    1494          14 :     break;
    1495             :   case t_FF_F2xq:
    1496         280 :     Pp = F2xqE_changepointinv(RgE_to_F2xqE(P,T), gel(e,3), T);
    1497         280 :     r = F2xqE_order(Pp, o, gel(e,1), T);
    1498         280 :     break;
    1499             :   default:
    1500        1456 :     Pp = FlxqE_changepointinv(RgE_to_FlxqE(P,T,pp), gel(e,3), T, pp);
    1501        1456 :     r = FlxqE_order(Pp, o, gel(e,1), T, pp);
    1502             :   }
    1503        1750 :   return gerepileupto(av, r);
    1504             : }
    1505             : 
    1506             : GEN
    1507          91 : FF_elllog(GEN E, GEN P, GEN Q, GEN o)
    1508             : {
    1509          91 :   pari_sp av = avma;
    1510          91 :   GEN fg = ellff_get_field(E), e = ellff_get_a4a6(E);
    1511             :   GEN r,T,p, Pp,Qp, e3;
    1512             :   ulong pp;
    1513          91 :   _getFF(fg,&T,&p,&pp);
    1514          91 :   switch(fg[1])
    1515             :   {
    1516             :   case t_FF_FpXQ:
    1517           0 :     e3 = FqV_to_FpXQV(gel(e,3),T);
    1518           0 :     Pp = FpXQE_changepointinv(RgE_to_FpXQE(P,T,p), e3, T, p);
    1519           0 :     Qp = FpXQE_changepointinv(RgE_to_FpXQE(Q,T,p), e3, T, p);
    1520           0 :     r = FpXQE_log(Pp, Qp, o, gel(e,1), T, p);
    1521           0 :     break;
    1522             :   case t_FF_F2xq:
    1523          42 :     Pp = F2xqE_changepointinv(RgE_to_F2xqE(P,T), gel(e,3), T);
    1524          42 :     Qp = F2xqE_changepointinv(RgE_to_F2xqE(Q,T), gel(e,3), T);
    1525          42 :     r = F2xqE_log(Pp, Qp, o, gel(e,1), T);
    1526          42 :     break;
    1527             :   default:
    1528          49 :     Pp = FlxqE_changepointinv(RgE_to_FlxqE(P,T,pp), gel(e,3), T, pp);
    1529          49 :     Qp = FlxqE_changepointinv(RgE_to_FlxqE(Q,T,pp), gel(e,3), T, pp);
    1530          49 :     r = FlxqE_log(Pp, Qp, o, gel(e,1), T, pp);
    1531             :   }
    1532          91 :   return gerepileupto(av, r);
    1533             : }
    1534             : 
    1535             : GEN
    1536        4984 : FF_ellweilpairing(GEN E, GEN P, GEN Q, GEN m)
    1537             : {
    1538        4984 :   GEN fg = ellff_get_field(E), e = ellff_get_a4a6(E);
    1539             :   GEN r,T,p, Pp,Qp, e3;
    1540             :   ulong pp;
    1541        4984 :   GEN z=_initFF(fg,&T,&p,&pp);
    1542        4984 :   pari_sp av = avma;
    1543        4984 :   switch(fg[1])
    1544             :   {
    1545             :   case t_FF_FpXQ:
    1546           7 :     e3 = FqV_to_FpXQV(gel(e,3),T);
    1547           7 :     Pp = FpXQE_changepointinv(RgE_to_FpXQE(P,T,p), e3, T, p);
    1548           7 :     Qp = FpXQE_changepointinv(RgE_to_FpXQE(Q,T,p), e3, T, p);
    1549           7 :     r = FpXQE_weilpairing(Pp, Qp, m, gel(e,1), T, p);
    1550           7 :     break;
    1551             :   case t_FF_F2xq:
    1552        4963 :     Pp = F2xqE_changepointinv(RgE_to_F2xqE(P,T), gel(e,3), T);
    1553        4963 :     Qp = F2xqE_changepointinv(RgE_to_F2xqE(Q,T), gel(e,3), T);
    1554        4963 :     r = F2xqE_weilpairing(Pp, Qp, m, gel(e,1), T);
    1555        4963 :     break;
    1556             :   default:
    1557          14 :     Pp = FlxqE_changepointinv(RgE_to_FlxqE(P,T,pp), gel(e,3), T, pp);
    1558          14 :     Qp = FlxqE_changepointinv(RgE_to_FlxqE(Q,T,pp), gel(e,3), T, pp);
    1559          14 :     r = FlxqE_weilpairing(Pp, Qp, m, gel(e,1), T, pp);
    1560             :   }
    1561        4984 :   r = gerepileupto(av, r);
    1562        4984 :   return _mkFF(fg,z,r);
    1563             : }
    1564             : 
    1565             : GEN
    1566          91 : FF_elltatepairing(GEN E, GEN P, GEN Q, GEN m)
    1567             : {
    1568          91 :   GEN fg = ellff_get_field(E), e = ellff_get_a4a6(E);
    1569             :   GEN r,T,p, Pp,Qp, e3;
    1570             :   ulong pp;
    1571          91 :   GEN z=_initFF(fg,&T,&p,&pp);
    1572          91 :   pari_sp av = avma;
    1573          91 :   switch(fg[1])
    1574             :   {
    1575             :   case t_FF_FpXQ:
    1576           7 :     e3 = FqV_to_FpXQV(gel(e,3),T);
    1577           7 :     Pp = FpXQE_changepointinv(RgE_to_FpXQE(P,T,p), e3, T, p);
    1578           7 :     Qp = FpXQE_changepointinv(RgE_to_FpXQE(Q,T,p), e3, T, p);
    1579           7 :     r = FpXQE_tatepairing(Pp, Qp, m, gel(e,1), T, p);
    1580           7 :     break;
    1581             :   case t_FF_F2xq:
    1582          28 :     Pp = F2xqE_changepointinv(RgE_to_F2xqE(P,T), gel(e,3), T);
    1583          28 :     Qp = F2xqE_changepointinv(RgE_to_F2xqE(Q,T), gel(e,3), T);
    1584          28 :     r = F2xqE_tatepairing(Pp, Qp, m, gel(e,1), T);
    1585          28 :     break;
    1586             :   default:
    1587          56 :     Pp = FlxqE_changepointinv(RgE_to_FlxqE(P,T,pp), gel(e,3), T, pp);
    1588          56 :     Qp = FlxqE_changepointinv(RgE_to_FlxqE(Q,T,pp), gel(e,3), T, pp);
    1589          56 :     r = FlxqE_tatepairing(Pp, Qp, m, gel(e,1), T, pp);
    1590             :   }
    1591          91 :   r = gerepileupto(av, r);
    1592          91 :   return _mkFF(fg,z,r);
    1593             : }
    1594             : 
    1595             : GEN
    1596      103649 : FFX_roots(GEN Pf, GEN ff)
    1597             : {
    1598      103649 :   pari_sp av = avma;
    1599             :   GEN r,T,p;
    1600             :   ulong pp;
    1601      103649 :   GEN P = FFX_to_raw(Pf, ff);
    1602      103649 :   _getFF(ff,&T,&p,&pp);
    1603      103649 :   switch(ff[1])
    1604             :   {
    1605             :   case t_FF_FpXQ:
    1606          85 :     r = FpXQX_roots(P, T, p);
    1607          85 :     break;
    1608             :   case t_FF_F2xq:
    1609       57512 :     r = F2xqX_roots(P, T);
    1610       57512 :     break;
    1611             :   default:
    1612       46052 :     r = FlxqX_roots(P, T, pp);
    1613             :   }
    1614      103649 :   return gerepilecopy(av, raw_to_FFC(r, ff));
    1615             : }
    1616             : 
    1617             : static GEN
    1618         420 : raw_to_FFX_fact(GEN F, GEN ff)
    1619             : {
    1620             :   GEN y, u, v;
    1621         420 :   GEN P = gel(F,1), E = gel(F,2);
    1622         420 :   long j, l = lg(P);
    1623         420 :   y = cgetg(3,t_MAT);
    1624         420 :   u = cgetg(l,t_COL); gel(y,1) = u;
    1625         420 :   v = cgetg(l,t_COL); gel(y,2) = v;
    1626        2044 :   for (j=1; j<l; j++)
    1627             :   {
    1628        1624 :     gel(u,j) = raw_to_FFX(gel(P,j), ff);
    1629        1624 :     gel(v,j) = utoi(uel(E,j));
    1630             :   }
    1631         420 :   return y;
    1632             : }
    1633             : 
    1634             : static GEN
    1635        2468 : FFX_zero(GEN ff, long v)
    1636             : {
    1637        2468 :   GEN r = cgetg(3,t_POL);
    1638        2468 :   r[1] = evalvarn(v);
    1639        2468 :   gel(r,2) = FF_zero(ff);
    1640        2468 :   return r;
    1641             : }
    1642             : 
    1643             : static GEN
    1644      104315 : FFX_wrap2(GEN Pf, GEN Qf, GEN ff, GEN FpXQX(GEN, GEN, GEN, GEN),
    1645             :           GEN F2xqX(GEN, GEN, GEN), GEN FlxqX(GEN, GEN, GEN, ulong))
    1646             : {
    1647      104315 :   pari_sp av = avma;
    1648             :   GEN r,T,p;
    1649             :   ulong pp;
    1650      104315 :   GEN P = FFX_to_raw(Pf, ff);
    1651      104315 :   GEN Q = FFX_to_raw(Qf, ff);
    1652      104315 :   _getFF(ff,&T,&p,&pp);
    1653      104315 :   switch(ff[1])
    1654             :   {
    1655             :   case t_FF_FpXQ:
    1656        1981 :     r = FpXQX(P, Q, T, p);
    1657        1981 :     break;
    1658             :   case t_FF_F2xq:
    1659       38951 :     r = F2xqX(P, Q, T);
    1660       38951 :     break;
    1661             :   default:
    1662       63383 :     r = FlxqX(P, Q, T, pp);
    1663             :   }
    1664      104315 :   if (!lgpol(r)) { avma = av; return FFX_zero(ff, varn(Pf)); }
    1665      101847 :   return gerepilecopy(av, raw_to_FFX(r, ff));
    1666             : }
    1667             : 
    1668             : GEN
    1669      101886 : FFX_mul(GEN Pf, GEN Qf, GEN ff)
    1670      101886 : { return FFX_wrap2(Pf, Qf, ff, FpXQX_mul, F2xqX_mul, FlxqX_mul); }
    1671             : 
    1672             : GEN
    1673        2261 : FFX_sqr(GEN Pf, GEN ff)
    1674             : {
    1675        2261 :   pari_sp av = avma;
    1676             :   GEN r,T,p;
    1677             :   ulong pp;
    1678        2261 :   GEN P = FFX_to_raw(Pf, ff);
    1679        2261 :   _getFF(ff,&T,&p,&pp);
    1680        2261 :   switch(ff[1])
    1681             :   {
    1682             :   case t_FF_FpXQ:
    1683         196 :     r = FpXQX_sqr(P, T, p);
    1684         196 :     break;
    1685             :   case t_FF_F2xq:
    1686         952 :     r = F2xqX_sqr(P, T);
    1687         952 :     break;
    1688             :   default:
    1689        1113 :     r = FlxqX_sqr(P, T, pp);
    1690             :   }
    1691        2261 :   if (!lgpol(r)) { avma = av; return FFX_zero(ff, varn(Pf)); }
    1692        2261 :   return gerepilecopy(av, raw_to_FFX(r, ff));
    1693             : }
    1694             : 
    1695             : GEN
    1696        2415 : FFX_rem(GEN Pf, GEN Qf, GEN ff)
    1697        2415 : { return FFX_wrap2(Pf, Qf, ff, FpXQX_rem, F2xqX_rem, FlxqX_rem); }
    1698             : 
    1699             : GEN
    1700           0 : FFXQ_sqr(GEN Pf, GEN Qf, GEN ff)
    1701           0 : { return FFX_wrap2(Pf, Qf, ff, FpXQXQ_sqr, F2xqXQ_sqr, FlxqXQ_sqr); }
    1702             : 
    1703             : GEN
    1704          14 : FFXQ_inv(GEN Pf, GEN Qf, GEN ff)
    1705          14 : { return FFX_wrap2(Pf, Qf, ff, FpXQXQ_inv, F2xqXQ_inv, FlxqXQ_inv); }
    1706             : 
    1707             : GEN
    1708          14 : FFXQ_mul(GEN Pf, GEN Qf, GEN Sf, GEN ff)
    1709             : {
    1710          14 :   pari_sp av = avma;
    1711             :   GEN r,T,p;
    1712             :   ulong pp;
    1713          14 :   GEN P = FFX_to_raw(Pf, ff);
    1714          14 :   GEN Q = FFX_to_raw(Qf, ff);
    1715          14 :   GEN S = FFX_to_raw(Sf, ff);
    1716          14 :   _getFF(ff,&T,&p,&pp);
    1717          14 :   switch(ff[1])
    1718             :   {
    1719             :   case t_FF_FpXQ:
    1720           0 :     r = FpXQXQ_mul(P, Q, S, T, p);
    1721           0 :     break;
    1722             :   case t_FF_F2xq:
    1723           7 :     r = F2xqXQ_mul(P, Q, S, T);
    1724           7 :     break;
    1725             :   default:
    1726           7 :     r = FlxqXQ_mul(P, Q, S, T, pp);
    1727             :   }
    1728          14 :   if (!lgpol(r)) { avma = av; return FFX_zero(ff, varn(Pf)); }
    1729          14 :   return gerepilecopy(av, raw_to_FFX(r, ff));
    1730             : }
    1731             : 
    1732             : long
    1733         168 : FFX_ispower(GEN Pf, long k, GEN ff, GEN *pt_r)
    1734             : {
    1735         168 :   pari_sp av = avma;
    1736             :   GEN P,T,p;
    1737             :   ulong pp;
    1738             :   long s;
    1739         168 :   if (degpol(Pf) % k) return 0;
    1740         168 :   P = FFX_to_raw(Pf, ff);
    1741         168 :   _getFF(ff,&T,&p,&pp);
    1742         168 :   switch(ff[1])
    1743             :   {
    1744             :   case t_FF_FpXQ:
    1745          56 :     s = FpXQX_ispower(P, k, T, p, pt_r);
    1746          56 :     break;
    1747             :   case t_FF_F2xq:
    1748          56 :     s = F2xqX_ispower(P, k, T, pt_r);
    1749          56 :     break;
    1750             :   default:
    1751          56 :     s = FlxqX_ispower(P, k, T, pp, pt_r);
    1752             :   }
    1753         168 :   if (s==0) { avma = av; return 0; }
    1754         147 :   if (pt_r)
    1755         147 :     *pt_r = gerepilecopy(av, raw_to_FFX(*pt_r, ff));
    1756           0 :   else avma = av;
    1757         147 :   return 1;
    1758             : }
    1759             : 
    1760             : GEN
    1761         420 : FFX_factor(GEN Pf, GEN ff)
    1762             : {
    1763         420 :   pari_sp av = avma;
    1764             :   GEN r,T,p;
    1765             :   ulong pp;
    1766         420 :   GEN P = FFX_to_raw(Pf, ff);
    1767         420 :   _getFF(ff,&T,&p,&pp);
    1768         420 :   switch(ff[1])
    1769             :   {
    1770             :   case t_FF_FpXQ:
    1771         121 :     r = FpXQX_factor(P, T, p);
    1772         121 :     break;
    1773             :   case t_FF_F2xq:
    1774         133 :     r = F2xqX_factor(P, T);
    1775         133 :     break;
    1776             :   default:
    1777         166 :     r = FlxqX_factor(P, T, pp);
    1778             :   }
    1779         420 :   return gerepilecopy(av, raw_to_FFX_fact(r, ff));
    1780             : }
    1781             : 
    1782             : GEN
    1783           0 : FqX_to_FFX(GEN x, GEN ff)
    1784             : {
    1785             :   long i, lx;
    1786           0 :   GEN y =  cgetg_copy(x,&lx);
    1787           0 :   y[1] = x[1];
    1788           0 :   for (i=2; i<lx; i++) gel(y,i) = Fq_to_FF(gel(x,i), ff);
    1789           0 :   return y;
    1790             : }
    1791             : 
    1792             : GEN
    1793        2271 : ffgen(GEN T, long v)
    1794             : {
    1795        2271 :   GEN A, p = NULL, ff = cgetg(5,t_FFELT);
    1796             :   long d;
    1797        2270 :   switch(typ(T))
    1798             :   {
    1799             :     case t_FFELT:
    1800          28 :       p = FF_p_i(T); T = FF_mod(T); d = degpol(T);
    1801          28 :       break;
    1802             :     case t_POL:
    1803         252 :       d = degpol(T); p = NULL;
    1804         252 :       if (d < 1 || !RgX_is_FpX(T, &p) || !p) pari_err_TYPE("ffgen",T);
    1805         252 :       T = RgX_to_FpX(T, p);
    1806             :       /* testing for irreducibility is too costly */
    1807         252 :       if (!FpX_is_squarefree(T,p)) pari_err_IRREDPOL("ffgen",T);
    1808         245 :       break;
    1809             :     case t_INT:
    1810        1592 :       d = ispseudoprimepower(T,&p);
    1811        1590 :       if (!d) pari_err_PRIME("ffgen",T);
    1812        1590 :       T = init_Fq(p, d, v);
    1813        1594 :       break;
    1814             :     case t_VEC: case t_COL:
    1815         399 :       if (lg(T) == 3) {
    1816         399 :         p = gel(T,1);
    1817         399 :         A = gel(T,2);
    1818         399 :         if (typ(p) == t_INT && typ(A) == t_INT)
    1819             :         {
    1820         399 :           d = itos(A);
    1821         399 :           T = init_Fq(p, d, v);
    1822         399 :           break;
    1823             :         }
    1824             :       }
    1825             :     default:
    1826           0 :       pari_err_TYPE("ffgen",T);
    1827           0 :       return NULL;
    1828             :   }
    1829        2266 :   if (v < 0) v = varn(T);
    1830        2266 :   if (lgefint(p)==3)
    1831             :   {
    1832        1992 :     ulong pp = p[2];
    1833        1992 :     long sv = evalvarn(v);
    1834        1992 :     if (pp==2)
    1835             :     {
    1836         894 :       ff[1] = t_FF_F2xq;
    1837         894 :       T = ZX_to_F2x(T); T[1] = sv;
    1838         893 :       A = polx_F2x(sv); if (d == 1) A = F2x_rem(A, T);
    1839         893 :       p = gen_2;
    1840             :     }
    1841             :     else
    1842             :     {
    1843        1098 :       ff[1] = t_FF_Flxq;
    1844        1098 :       T = ZX_to_Flx(T,pp); T[1] = sv;
    1845        1098 :       A = polx_Flx(sv); if (d == 1) A = Flx_rem(A, T, pp);
    1846        1098 :       p = icopy(p);
    1847             :     }
    1848             :   }
    1849             :   else
    1850             :   {
    1851         274 :     ff[1] = t_FF_FpXQ;
    1852         274 :     setvarn(T,v);
    1853         274 :     A = pol_x(v); if (d == 1) A = FpX_rem(A, T, p);
    1854         274 :     p = icopy(p);
    1855             :   }
    1856        2265 :   gel(ff,2) = A;
    1857        2265 :   gel(ff,3) = T;
    1858        2265 :   gel(ff,4) = p; return ff;
    1859             : }
    1860             : 
    1861             : GEN
    1862        2898 : p_to_FF(GEN p, long v)
    1863             : {
    1864             :   GEN A, T;
    1865        2898 :   GEN ff = cgetg(5,t_FFELT);
    1866        2898 :   if (lgefint(p)==3)
    1867             :   {
    1868        2894 :     ulong pp = p[2];
    1869        2894 :     long sv = evalvarn(v);
    1870        2894 :     if (pp==2)
    1871             :     {
    1872         294 :       ff[1] = t_FF_F2xq;
    1873         294 :       T = polx_F2x(sv);
    1874         294 :       A = pol1_F2x(sv);
    1875         294 :       p = gen_2;
    1876             :     }
    1877             :     else
    1878             :     {
    1879        2600 :       ff[1] = t_FF_Flxq;
    1880        2600 :       T = polx_Flx(sv);
    1881        2600 :       A = pol1_Flx(sv);
    1882        2600 :       p = icopy(p);
    1883             :     }
    1884             :   }
    1885             :   else
    1886             :   {
    1887           4 :     ff[1] = t_FF_FpXQ;
    1888           4 :     T = pol_x(v);
    1889           4 :     A = pol_1(v);
    1890           4 :     p = icopy(p);
    1891             :   }
    1892        2898 :   gel(ff,2) = A;
    1893        2898 :   gel(ff,3) = T;
    1894        2898 :   gel(ff,4) = p; return ff;
    1895             : }
    1896             : GEN
    1897         357 : Tp_to_FF(GEN T, GEN p)
    1898             : {
    1899             :   GEN A, ff;
    1900             :   long v;
    1901         357 :   if (!T) return p_to_FF(p,0);
    1902         245 :   ff = cgetg(5,t_FFELT);
    1903         245 :   v = varn(T);
    1904         245 :   if (lgefint(p)==3)
    1905             :   {
    1906         245 :     ulong pp = p[2];
    1907         245 :     long sv = evalvarn(v);
    1908         245 :     if (pp==2)
    1909             :     {
    1910         119 :       ff[1] = t_FF_F2xq;
    1911         119 :       T = ZX_to_F2x(T);
    1912         119 :       A = pol1_F2x(sv);
    1913         119 :       p = gen_2;
    1914             :     }
    1915             :     else
    1916             :     {
    1917         126 :       ff[1] = t_FF_Flxq;
    1918         126 :       T = ZX_to_Flx(T, pp);
    1919         126 :       A = pol1_Flx(sv);
    1920         126 :       p = icopy(p);
    1921             :     }
    1922             :   }
    1923             :   else
    1924             :   {
    1925           0 :     ff[1] = t_FF_FpXQ;
    1926           0 :     T = ZX_copy(T);
    1927           0 :     A = pol_1(v);
    1928           0 :     p = icopy(p);
    1929             :   }
    1930         245 :   gel(ff,2) = A;
    1931         245 :   gel(ff,3) = T;
    1932         245 :   gel(ff,4) = p; return ff;
    1933             : }
    1934             : 
    1935             : GEN
    1936          28 : fforder(GEN x, GEN o)
    1937             : {
    1938          28 :   if (typ(x)!=t_FFELT) pari_err_TYPE("fforder",x);
    1939          28 :   return FF_order(x,o);
    1940             : }
    1941             : 
    1942             : GEN
    1943         385 : ffprimroot(GEN x, GEN *o)
    1944             : {
    1945         385 :   if (typ(x)!=t_FFELT) pari_err_TYPE("ffprimroot",x);
    1946         385 :   return FF_primroot(x,o);
    1947             : }
    1948             : 
    1949             : GEN
    1950         175 : fflog(GEN x, GEN g, GEN o)
    1951             : {
    1952         175 :   if (typ(x)!=t_FFELT) pari_err_TYPE("fflog",x);
    1953         175 :   if (typ(g)!=t_FFELT) pari_err_TYPE("fflog",g);
    1954         175 :   return FF_log(x,g,o);
    1955             : }
    1956             : 
    1957             : GEN
    1958      144466 : ffrandom(GEN ff)
    1959             : {
    1960             :   ulong pp;
    1961      144466 :   GEN r, T, p, z = _initFF(ff,&T,&p,&pp);
    1962      144466 :   switch(ff[1])
    1963             :   {
    1964             :   case t_FF_FpXQ:
    1965        6940 :     r = random_FpX(degpol(T), varn(T), p);
    1966        6940 :     break;
    1967             :   case t_FF_F2xq:
    1968      117271 :     r = random_F2x(F2x_degree(T), T[1]);
    1969      117271 :     break;
    1970             :   default:
    1971       20255 :     r = random_Flx(degpol(T), T[1], pp);
    1972             :   }
    1973      144466 :   return _mkFF(ff,z,r);
    1974             : }
    1975             : 
    1976             : int
    1977     9258880 : Rg_is_FF(GEN c, GEN *ff)
    1978             : {
    1979     9258880 :   switch(typ(c))
    1980             :   {
    1981             :   case t_FFELT:
    1982       43190 :     if (!*ff) *ff = c;
    1983       42770 :     else if (!FF_samefield(*ff, c)) return 0;
    1984       43190 :     break;
    1985             :   default:
    1986     9215690 :     return 0;
    1987             :   }
    1988       43190 :   return 1;
    1989             : }
    1990             : 
    1991             : int
    1992     9218371 : RgC_is_FFC(GEN x, GEN *ff)
    1993             : {
    1994     9218371 :   long i, lx = lg(x);
    1995     9261561 :   for (i=lx-1; i>0; i--)
    1996     9258880 :     if (!Rg_is_FF(gel(x,i), ff)) return 0;
    1997        2681 :   return (*ff != NULL);
    1998             : }
    1999             : 
    2000             : int
    2001     9216386 : RgM_is_FFM(GEN x, GEN *ff)
    2002             : {
    2003     9216386 :   long j, lx = lg(x);
    2004     9218731 :   for (j=lx-1; j>0; j--)
    2005     9218290 :     if (!RgC_is_FFC(gel(x,j), ff)) return 0;
    2006         441 :   return (*ff != NULL);
    2007             : }
    2008             : 
    2009             : static GEN
    2010        1625 : FqC_to_FpXQC(GEN x, GEN T, GEN p)
    2011             : {
    2012             :   long i, lx;
    2013        1625 :   GEN y = cgetg_copy(x,&lx);
    2014       29103 :   for(i=1; i<lx; i++)
    2015       27478 :     gel(y, i) = Fq_to_FpXQ(gel(x, i), T, p);
    2016        1625 :   return y;
    2017             : }
    2018             : 
    2019             : static GEN
    2020         173 : FqM_to_FpXQM(GEN x, GEN T, GEN p)
    2021             : {
    2022             :   long i, lx;
    2023         173 :   GEN y = cgetg_copy(x,&lx);
    2024        1735 :   for(i=1; i<lx; i++)
    2025        1562 :     gel(y, i) = FqC_to_FpXQC(gel(x, i), T, p);
    2026         173 :   return y;
    2027             : }
    2028             : 
    2029             : /* for functions t_MAT -> t_MAT */
    2030             : static GEN
    2031         294 : FFM_wrap(GEN M, GEN ff, GEN (*Fq)(GEN,GEN,GEN),
    2032             :                        GEN (*Flxq)(GEN,GEN,ulong),
    2033             :                        GEN (*F2xq)(GEN,GEN))
    2034             : {
    2035         294 :   pari_sp av = avma;
    2036             :   ulong pp;
    2037             :   GEN T, p;
    2038         294 :   _getFF(ff,&T,&p,&pp); M = FFM_to_raw(M, ff);
    2039         294 :   switch(ff[1])
    2040             :   {
    2041         112 :   case t_FF_FpXQ: M = Fq(M,T,p); if (M) M = FqM_to_FpXQM(M,T,p);
    2042         112 :                   break;
    2043          70 :   case t_FF_F2xq: M = F2xq(M,T); break;
    2044         112 :   default: M = Flxq(M,T,pp); break;
    2045             :   }
    2046         294 :   if (!M) { avma = av; return NULL; }
    2047         259 :   return gerepilecopy(av, raw_to_FFM(M, ff));
    2048             : }
    2049             : 
    2050             : /* for functions (t_MAT, t_MAT) -> t_MAT */
    2051             : static GEN
    2052         399 : FFM_FFM_wrap(GEN M, GEN N, GEN ff,
    2053             :              GEN (*Fq)(GEN, GEN, GEN, GEN),
    2054             :              GEN (*Flxq)(GEN, GEN, GEN, ulong),
    2055             :              GEN (*F2xq)(GEN, GEN, GEN))
    2056             : {
    2057         399 :   pari_sp av = avma;
    2058             :   ulong pp;
    2059             :   GEN T, p;
    2060         399 :   int is_sqr = M==N;
    2061         399 :   _getFF(ff, &T, &p, &pp);
    2062         399 :   M = FFM_to_raw(M, ff);
    2063         399 :   N = is_sqr? M: FFM_to_raw(N, ff);
    2064         399 :   switch(ff[1])
    2065             :   {
    2066         110 :   case t_FF_FpXQ: M = Fq(M, N, T, p); if (M) M = FqM_to_FpXQM(M, T, p);
    2067         110 :                   break;
    2068          49 :   case t_FF_F2xq: M = F2xq(M, N, T); break;
    2069         240 :   default: M = Flxq(M, N, T, pp); break;
    2070             :   }
    2071         399 :   if (!M) { avma = av; return NULL; }
    2072         308 :   return gerepilecopy(av, raw_to_FFM(M, ff));
    2073             : }
    2074             : 
    2075             : /* for functions (t_MAT, t_COL) -> t_COL */
    2076             : static GEN
    2077         196 : FFM_FFC_wrap(GEN M, GEN C, GEN ff,
    2078             :              GEN (*Fq)(GEN, GEN, GEN, GEN),
    2079             :              GEN (*Flxq)(GEN, GEN, GEN, ulong),
    2080             :              GEN (*F2xq)(GEN, GEN, GEN))
    2081             : {
    2082         196 :   pari_sp av = avma;
    2083             :   ulong pp;
    2084             :   GEN T, p;
    2085         196 :   _getFF(ff, &T, &p, &pp);
    2086         196 :   M = FFM_to_raw(M, ff);
    2087         196 :   C = FFC_to_raw(C, ff);
    2088         196 :   switch(ff[1])
    2089             :   {
    2090          63 :   case t_FF_FpXQ: C = Fq(M, C, T, p); if (C) C = FqC_to_FpXQC(C, T, p);
    2091          63 :                   break;
    2092          70 :   case t_FF_F2xq: C = F2xq(M, C, T); break;
    2093          63 :   default: C = Flxq(M, C, T, pp); break;
    2094             :   }
    2095         196 :   if (!C) { avma = av; return NULL; }
    2096         147 :   return gerepilecopy(av, raw_to_FFC(C, ff));
    2097             : }
    2098             : 
    2099             : GEN
    2100          77 : FFM_ker(GEN M, GEN ff)
    2101          77 : { return FFM_wrap(M,ff, &FqM_ker,&FlxqM_ker,&F2xqM_ker); }
    2102             : GEN
    2103          49 : FFM_image(GEN M, GEN ff)
    2104          49 : { return FFM_wrap(M,ff, &FqM_image,&FlxqM_image,&F2xqM_image); }
    2105             : GEN
    2106         147 : FFM_inv(GEN M, GEN ff)
    2107         147 : { return FFM_wrap(M,ff, &FqM_inv,&FlxqM_inv,&F2xqM_inv); }
    2108             : GEN
    2109          21 : FFM_suppl(GEN M, GEN ff)
    2110          21 : { return FFM_wrap(M,ff, &FqM_suppl,&FlxqM_suppl,&F2xqM_suppl); }
    2111             : 
    2112             : GEN
    2113          70 : FFM_deplin(GEN M, GEN ff)
    2114             : {
    2115          70 :   pari_sp av = avma;
    2116             :   ulong pp;
    2117             :   GEN C, T, p;
    2118          70 :   _getFF(ff, &T, &p, &pp); M = FFM_to_raw(M, ff);
    2119          70 :   switch(ff[1]) {
    2120          28 :   case t_FF_FpXQ: C = FqM_deplin(M, T, p);
    2121          28 :     if (C) C = FqC_to_FpXQC(C, T, p); break;
    2122          14 :   case t_FF_F2xq: C = F2xqM_deplin(M, T); break;
    2123          28 :   default: C = FlxqM_deplin(M, T, pp); break;
    2124             :   }
    2125          70 :   if (!C) { avma = av; return NULL; }
    2126          35 :   return gerepilecopy(av, raw_to_FFC(C, ff));
    2127             : }
    2128             : 
    2129             : GEN
    2130          21 : FFM_indexrank(GEN M, GEN ff)
    2131             : {
    2132          21 :   pari_sp av = avma;
    2133             :   ulong pp;
    2134             :   GEN R, T, p;
    2135          21 :   _getFF(ff,&T,&p,&pp); M = FFM_to_raw(M, ff);
    2136          21 :   switch(ff[1]) {
    2137           7 :   case t_FF_FpXQ: R = FqM_indexrank(M,T,p); break;
    2138           7 :   case t_FF_F2xq: R = F2xqM_indexrank(M,T); break;
    2139           7 :   default: R = FlxqM_indexrank(M,T,pp); break;
    2140             :   }
    2141          21 :   return gerepileupto(av, R);
    2142             : }
    2143             : 
    2144             : long
    2145          63 : FFM_rank(GEN M, GEN ff)
    2146             : {
    2147          63 :   pari_sp av = avma;
    2148             :   long r;
    2149             :   ulong pp;
    2150             :   GEN T, p;
    2151          63 :   _getFF(ff,&T,&p,&pp); M = FFM_to_raw(M, ff);
    2152          63 :   switch(ff[1])
    2153             :   {
    2154          28 :   case t_FF_FpXQ: r = FqM_rank(M,T,p); break;
    2155           7 :   case t_FF_F2xq: r = F2xqM_rank(M,T); break;
    2156          28 :   default: r = FlxqM_rank(M,T,pp); break;
    2157             :   }
    2158          63 :   avma = av; return r;
    2159             : }
    2160             : GEN
    2161          63 : FFM_det(GEN M, GEN ff)
    2162             : {
    2163          63 :   pari_sp av = avma;
    2164             :   ulong pp;
    2165             :   GEN d, T, p;
    2166          63 :   _getFF(ff,&T,&p,&pp); M = FFM_to_raw(M, ff);
    2167          63 :   switch(ff[1])
    2168             :   {
    2169          28 :   case t_FF_FpXQ: d = FqM_det(M,T,p); break;
    2170           7 :   case t_FF_F2xq: d = F2xqM_det(M,T); break;
    2171          28 :   default: d = FlxqM_det(M,T,pp); break;
    2172             :   }
    2173          63 :   return gerepilecopy(av, mkFF_i(ff, d));
    2174             : }
    2175             : 
    2176             : GEN
    2177          56 : FFM_FFC_gauss(GEN M, GEN C, GEN ff)
    2178             : {
    2179          56 :   return FFM_FFC_wrap(M, C, ff, FqM_FqC_gauss,
    2180             :                       FlxqM_FlxqC_gauss, F2xqM_F2xqC_gauss);
    2181             : }
    2182             : 
    2183             : GEN
    2184          63 : FFM_gauss(GEN M, GEN N, GEN ff)
    2185             : {
    2186          63 :   return FFM_FFM_wrap(M, N, ff, FqM_gauss,
    2187             :                       FlxqM_gauss, F2xqM_gauss);
    2188             : }
    2189             : 
    2190             : GEN
    2191          63 : FFM_FFC_invimage(GEN M, GEN C, GEN ff)
    2192             : {
    2193          63 :   return FFM_FFC_wrap(M, C, ff, FqM_FqC_invimage,
    2194             :                       FlxqM_FlxqC_invimage, F2xqM_F2xqC_invimage);
    2195             : }
    2196             : 
    2197             : GEN
    2198         105 : FFM_invimage(GEN M, GEN N, GEN ff)
    2199             : {
    2200         105 :   return FFM_FFM_wrap(M, N, ff, FqM_invimage,
    2201             :                       FlxqM_invimage, F2xqM_invimage);
    2202             : }
    2203             : 
    2204             : GEN
    2205          77 : FFM_FFC_mul(GEN M, GEN C, GEN ff)
    2206             : {
    2207          77 :   return FFM_FFC_wrap(M, C, ff, FqM_FqC_mul,
    2208             :                       FlxqM_FlxqC_mul, F2xqM_F2xqC_mul);
    2209             : }
    2210             : 
    2211             : GEN
    2212         231 : FFM_mul(GEN M, GEN N, GEN ff)
    2213             : {
    2214         231 :   return FFM_FFM_wrap(M, N, ff, FqM_mul, FlxqM_mul, F2xqM_mul);
    2215             : }

Generated by: LCOV version 1.11