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 - ellisog.c (source / functions) Hit Total Coverage
Test: PARI/GP v2.10.0 lcov report (development 21059-cbe0d6a) Lines: 925 953 97.1 %
Date: 2017-09-22 06:24:58 Functions: 75 75 100.0 %
Legend: Lines: hit not hit

          Line data    Source code
       1             : /* Copyright (C) 2014 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             : /* Return 1 if the point Q is a Weierstrass (2-torsion) point of the
      18             :  * curve E, return 0 otherwise */
      19             : static long
      20         903 : ellisweierstrasspoint(GEN E, GEN Q)
      21         903 : { return ell_is_inf(Q) || gequal0(ec_dmFdy_evalQ(E, Q)); }
      22             : 
      23             : 
      24             : /* Given an elliptic curve E = [a1, a2, a3, a4, a6] and t,w in the ring of
      25             :  * definition of E, return the curve
      26             :  *  E' = [a1, a2, a3, a4 - 5t, a6 - (E.b2 t + 7w)] */
      27             : static GEN
      28        3255 : make_velu_curve(GEN E, GEN t, GEN w)
      29             : {
      30        3255 :   GEN A4, A6, a1 = ell_get_a1(E), a2 = ell_get_a2(E), a3 = ell_get_a3(E);
      31        3255 :   A4 = gsub(ell_get_a4(E), gmulsg(5L, t));
      32        3255 :   A6 = gsub(ell_get_a6(E), gadd(gmul(ell_get_b2(E), t), gmulsg(7L, w)));
      33        3255 :   return mkvec5(a1,a2,a3,A4,A6);
      34             : }
      35             : 
      36             : /* If phi = (f(x)/h(x)^2, g(x,y)/h(x)^3) is an isogeny, return the
      37             :  * variables x and y in a vecsmall */
      38             : INLINE void
      39        1624 : get_isog_vars(GEN phi, long *vx, long *vy)
      40             : {
      41        1624 :   *vx = varn(gel(phi, 1));
      42        1624 :   *vy = varn(gel(phi, 2));
      43        1624 :   if (*vy == *vx) *vy = gvar2(gel(phi,2));
      44        1624 : }
      45             : 
      46             : static GEN
      47        4592 : RgX_homogenous_evalpow(GEN P, GEN A, GEN B)
      48             : {
      49        4592 :   pari_sp av = avma;
      50             :   long d, i, v;
      51             :   GEN s;
      52        4592 :   if (typ(P)!=t_POL)
      53        1435 :     return mkvec2(P, gen_1);
      54        3157 :   d = degpol(P); v = varn(A);
      55        3157 :   s = scalarpol_shallow(gel(P, d+2), v);
      56       17073 :   for (i = d-1; i >= 0; i--)
      57             :   {
      58       13916 :     s = gadd(gmul(s, A), gmul(gel(B,d+1-i), gel(P,i+2)));
      59       13916 :     if (gc_needed(av,1))
      60             :     {
      61           0 :       if(DEBUGMEM>1) pari_warn(warnmem,"RgX_homogenous_eval(%ld)",i);
      62           0 :       s = gerepileupto(av, s);
      63             :     }
      64             :   }
      65        3157 :   s = gerepileupto(av, s);
      66        3157 :   return mkvec2(s, gel(B,d+1));
      67             : }
      68             : 
      69             : static GEN
      70        1792 : RgXQX_homogenous_evalpow(GEN P, GEN A, GEN B, GEN T)
      71             : {
      72        1792 :   pari_sp av = avma;
      73        1792 :   long i, d = degpol(P), v = varn(A);
      74             :   GEN s;
      75        1792 :   if (signe(P)==0) return mkvec2(pol_0(v), pol_1(v));
      76        1568 :   s = scalarpol_shallow(gel(P, d+2), v);
      77        5229 :   for (i = d-1; i >= 0; i--)
      78             :   {
      79        3661 :     s = RgX_add(RgXQX_mul(s, A, T), RgXQX_RgXQ_mul(gel(B,d+1-i), gel(P,i+2), T));
      80        3661 :     if (gc_needed(av,1))
      81             :     {
      82           0 :       if(DEBUGMEM>1) pari_warn(warnmem,"RgX_homogenous_eval(%ld)",i);
      83           0 :       s = gerepileupto(av, s);
      84             :     }
      85             :   }
      86        1568 :   s = gerepileupto(av, s);
      87        1568 :   return mkvec2(s, gel(B,d+1));
      88             : }
      89             : 
      90             : /* Given isogenies F:E' -> E and G:E'' -> E', return the composite
      91             :  * isogeny F o G:E'' -> E */
      92             : static GEN
      93        1155 : ellcompisog(GEN F, GEN G)
      94             : {
      95        1155 :   pari_sp av = avma;
      96             :   GEN Fv, Gh, Gh2, Gh3, f, g, h, h2, h3, den, num;
      97             :   GEN K, K2, K3, F0, F1, g0, g1, Gp;
      98             :   long v, vx, vy, d;
      99        1155 :   checkellisog(F);
     100        1148 :   checkellisog(G);
     101        1148 :   get_isog_vars(F, &vx, &vy);
     102        1148 :   v = fetch_var_higher();
     103        1148 :   Fv = shallowcopy(gel(F,3)); setvarn(Fv, v);
     104        1148 :   Gh = gel(G,3); Gh2 = gsqr(Gh); Gh3 = gmul(Gh, Gh2);
     105        1148 :   K = gmul(polresultant0(Fv, deg1pol(gneg(Gh2),gel(G,1), v), v, 0), Gh);
     106        1148 :   delete_var();
     107        1148 :   K = RgX_normalize(RgX_div(K, RgX_gcd(K,deriv(K,0))));
     108        1148 :   K2 = gsqr(K); K3 = gmul(K, K2);
     109        1148 :   F0 = polcoeff0(gel(F,2), 0, vy); F1 = polcoeff0(gel(F,2), 1, vy);
     110        1148 :   d = maxss(maxss(degpol(gel(F,1)),degpol(gel(F,3))),maxss(degpol(F0),degpol(F1)));
     111        1148 :   Gp = gpowers(Gh2, d);
     112        1148 :   f  = RgX_homogenous_evalpow(gel(F,1), gel(G,1), Gp);
     113        1148 :   g0 = RgX_homogenous_evalpow(F0, gel(G,1), Gp);
     114        1148 :   g1 = RgX_homogenous_evalpow(F1, gel(G,1), Gp);
     115        1148 :   h =  RgX_homogenous_evalpow(gel(F,3), gel(G,1), Gp);
     116        1148 :   h2 = mkvec2(gsqr(gel(h,1)), gsqr(gel(h,2)));
     117        1148 :   h3 = mkvec2(gmul(gel(h,1),gel(h2,1)), gmul(gel(h,2),gel(h2,2)));
     118        1148 :   f  = gdiv(gmul(gmul(K2, gel(f,1)),gel(h2,2)), gmul(gel(f,2), gel(h2,1)));
     119        1148 :   den = gmul(Gh3, gel(g1,2));
     120        1148 :   num = gadd(gmul(gel(g0,1),den), gmul(gmul(gel(G,2),gel(g1,1)),gel(g0,2)));
     121        1148 :   g = gdiv(gmul(gmul(K3,num),gel(h3,2)),gmul(gmul(gel(g0,2),den), gel(h3,1)));
     122        1148 :   return gerepilecopy(av, mkvec3(f,g,K));
     123             : }
     124             : 
     125             : static GEN
     126        4032 : to_RgX(GEN P, long vx)
     127             : {
     128        4032 :   return typ(P) == t_POL ? lift(P): scalarpol_shallow(lift(P), vx);
     129             : }
     130             : 
     131             : static GEN
     132         448 : divy(GEN P0, GEN P1, GEN Q, GEN T, long vy)
     133             : {
     134         448 :   GEN DP0, P0r = Q_remove_denom(P0, &DP0), P0D;
     135         448 :   GEN DP1, P1r = Q_remove_denom(P1, &DP1), P1D;
     136         448 :   GEN DQ, Qr = Q_remove_denom(Q, &DQ), P2;
     137         448 :   P0D = RgXQX_div(P0r, Qr, T);
     138         448 :   if (DP0) P0D = gdiv(P0D, DP0);
     139         448 :   P1D = RgXQX_div(P1r, Qr, T);
     140         448 :   if (DP1) P1D = gdiv(P1D, DP1);
     141         448 :   P2 = gadd(gmul(P1D, pol_x(vy)), P0D);
     142         448 :   if (DQ) P2 = gmul(P2, DQ);
     143         448 :   return P2;
     144             : }
     145             : 
     146             : static GEN
     147        1400 : ellnfcompisog(GEN nf, GEN F, GEN G)
     148             : {
     149        1400 :   pari_sp av = avma;
     150             :   GEN Fv, Gh, Gh2, Gh3, f, g, gd, h, h21, h22, h31, h32, den;
     151             :   GEN K, K2, K3, F0, F1, G0, G1, g0, g1, Gp;
     152             :   GEN num0, num1, gn0, gn1;
     153             :   GEN g0d, g01, k3h32;
     154             :   GEN T, res;
     155             :   pari_timer ti;
     156             :   long v, vx, vy, d;
     157        1400 :   if (!nf) return ellcompisog(F, G);
     158         448 :   T = nf_get_pol(nf);
     159         448 :   timer_start(&ti);
     160         448 :   checkellisog(F);
     161         448 :   checkellisog(G);
     162         448 :   get_isog_vars(F, &vx, &vy);
     163         448 :   v = fetch_var_higher();
     164         448 :   Fv = shallowcopy(gel(F,3)); setvarn(Fv, v);
     165         448 :   Gh = lift(gel(G,3)); Gh2 = RgXQX_sqr(Gh, T); Gh3 = RgXQX_mul(Gh, Gh2, T);
     166         448 :   res = to_RgX(polresultant0(Fv, deg1pol(gmul(gneg(Gh2),gmodulo(gen_1,T)),gel(G,1), v), v, 0),vx);
     167         448 :   delete_var();
     168         448 :   K = Q_remove_denom(RgXQX_mul(res, Gh, T), NULL);
     169         448 :   if (DEBUGLEVEL) timer_printf(&ti,"ellnfcompisog: resultant");
     170         448 :   K = RgXQX_div(K, nfgcd(K, deriv(K,0), T, NULL), T);
     171         448 :   K = RgX_normalize(K);
     172         448 :   if (DEBUGLEVEL) timer_printf(&ti,"ellnfcompisog: nfgcd");
     173         448 :   K2 = RgXQX_sqr(K, T); K3 = RgXQX_mul(K, K2, T);
     174         448 :   F0 = to_RgX(polcoeff0(gel(F,2), 0, vy), vx);
     175         448 :   F1 = to_RgX(polcoeff0(gel(F,2), 1, vy), vx);
     176         448 :   G0 = to_RgX(polcoeff0(gel(G,2), 0, vy), vx);
     177         448 :   G1 = to_RgX(polcoeff0(gel(G,2), 1, vy), vx);
     178         448 :   d = maxss(maxss(degpol(gel(F,1)),degpol(gel(F,3))),maxss(degpol(F0),degpol(F1)));
     179         448 :   Gp = RgXQX_powers(Gh2, d, T);
     180         448 :   f  = RgXQX_homogenous_evalpow(to_RgX(gel(F,1),vx), gel(G,1), Gp, T);
     181         448 :   g0 = RgXQX_homogenous_evalpow(F0, to_RgX(gel(G,1),vx), Gp, T);
     182         448 :   g1 = RgXQX_homogenous_evalpow(F1, to_RgX(gel(G,1),vx), Gp, T);
     183         448 :   h  = RgXQX_homogenous_evalpow(to_RgX(gel(F,3),vx), gel(G,1), Gp, T);
     184         448 :   if (DEBUGLEVEL) timer_printf(&ti,"ellnfcompisog: evalpow");
     185         448 :   h21 = RgXQX_sqr(gel(h,1),T);
     186         448 :   h22 = RgXQX_sqr(gel(h,2),T);
     187         448 :   h31 = RgXQX_mul(gel(h,1), h21,T);
     188         448 :   h32 = RgXQX_mul(gel(h,2), h22,T);
     189         448 :   if (DEBUGLEVEL) timer_printf(&ti,"h");
     190         448 :   f  = RgXQX_div(RgXQX_mul(RgXQX_mul(K2, gel(f,1), T), h22, T),
     191         448 :                            RgXQX_mul(gel(f,2), h21, T), T);
     192         448 :   if (DEBUGLEVEL) timer_printf(&ti,"f");
     193         448 :   den = RgXQX_mul(Gh3, gel(g1,2), T);
     194         448 :   if (DEBUGLEVEL) timer_printf(&ti,"ellnfcompisog: den");
     195         448 :   g0d = RgXQX_mul(gel(g0,1),den, T);
     196         448 :   g01 = RgXQX_mul(gel(g1,1),gel(g0,2),T);
     197         448 :   num0 = RgX_add(g0d, RgXQX_mul(G0,g01, T));
     198         448 :   num1 = RgXQX_mul(G1,g01, T);
     199         448 :   if (DEBUGLEVEL) timer_printf(&ti,"ellnfcompisog: num");
     200         448 :   k3h32 = RgXQX_mul(K3,h32,T);
     201         448 :   gn0 = RgXQX_mul(num0, k3h32, T);
     202         448 :   gn1 = RgXQX_mul(num1, k3h32, T);
     203         448 :   if (DEBUGLEVEL) timer_printf(&ti,"ellnfcompisog: gn");
     204         448 :   gd = RgXQX_mul(RgXQX_mul(gel(g0,2), den, T), h31, T);
     205         448 :   if (DEBUGLEVEL) timer_printf(&ti,"ellnfcompisog: gd");
     206         448 :   g = divy(gn0, gn1, gd, T, vy);
     207         448 :   if (DEBUGLEVEL) timer_printf(&ti,"ellnfcompisog: divy");
     208         448 :   return gerepilecopy(av, gmul(mkvec3(f,g,K),gmodulo(gen_1,T)));
     209             : }
     210             : 
     211             : /* Given an isogeny phi from ellisogeny() and a point P in the domain of phi,
     212             :  * return phi(P) */
     213             : GEN
     214          70 : ellisogenyapply(GEN phi, GEN P)
     215             : {
     216          70 :   pari_sp ltop = avma;
     217             :   GEN f, g, h, img_f, img_g, img_h, img_h2, img_h3, img, tmp;
     218             :   long vx, vy;
     219          70 :   if (lg(P) == 4) return ellcompisog(phi,P);
     220          49 :   checkellisog(phi);
     221          49 :   checkellpt(P);
     222          42 :   if (ell_is_inf(P)) return ellinf();
     223          28 :   f = gel(phi, 1);
     224          28 :   g = gel(phi, 2);
     225          28 :   h = gel(phi, 3);
     226          28 :   get_isog_vars(phi, &vx, &vy);
     227          28 :   img_h = poleval(h, gel(P, 1));
     228          28 :   if (gequal0(img_h)) { avma = ltop; return ellinf(); }
     229             : 
     230          21 :   img_h2 = gsqr(img_h);
     231          21 :   img_h3 = gmul(img_h, img_h2);
     232          21 :   img_f = poleval(f, gel(P, 1));
     233             :   /* FIXME: This calculation of g is perhaps not as efficient as it could be */
     234          21 :   tmp = gsubst(g, vx, gel(P, 1));
     235          21 :   img_g = gsubst(tmp, vy, gel(P, 2));
     236          21 :   img = cgetg(3, t_VEC);
     237          21 :   gel(img, 1) = gdiv(img_f, img_h2);
     238          21 :   gel(img, 2) = gdiv(img_g, img_h3);
     239          21 :   return gerepileupto(ltop, img);
     240             : }
     241             : 
     242             : /* isog = [f, g, h] = [x, y, 1] = identity */
     243             : static GEN
     244         252 : isog_identity(long vx, long vy)
     245         252 : { return mkvec3(pol_x(vx), pol_x(vy), pol_1(vx)); }
     246             : 
     247             : /* Returns an updated value for isog based (primarily) on tQ and uQ. Used to
     248             :  * iteratively compute the isogeny corresponding to a subgroup generated by a
     249             :  * given point. Ref: Equation 8 in Velu's paper.
     250             :  * isog = NULL codes the identity */
     251             : static GEN
     252         532 : update_isogeny_polys(GEN isog, GEN E, GEN Q, GEN tQ, GEN uQ, long vx, long vy)
     253             : {
     254         532 :   pari_sp ltop = avma, av;
     255         532 :   GEN xQ = gel(Q, 1), yQ = gel(Q, 2);
     256         532 :   GEN rt = deg1pol_shallow(gen_1, gneg(xQ), vx);
     257         532 :   GEN a1 = ell_get_a1(E), a3 = ell_get_a3(E);
     258             : 
     259         532 :   GEN gQx = ec_dFdx_evalQ(E, Q);
     260         532 :   GEN gQy = ec_dFdy_evalQ(E, Q);
     261             :   GEN tmp1, tmp2, tmp3, tmp4, f, g, h, rt_sqr, res;
     262             : 
     263             :   /* g -= uQ * (2 * y + E.a1 * x + E.a3)
     264             :    *   + tQ * rt * (E.a1 * rt + y - yQ)
     265             :    *   + rt * (E.a1 * uQ - gQx * gQy) */
     266         532 :   av = avma;
     267         532 :   tmp1 = gmul(uQ, gadd(deg1pol_shallow(gen_2, gen_0, vy),
     268             :                        deg1pol_shallow(a1, a3, vx)));
     269         532 :   tmp1 = gerepileupto(av, tmp1);
     270         532 :   av = avma;
     271         532 :   tmp2 = gmul(tQ, gadd(gmul(a1, rt),
     272             :                        deg1pol_shallow(gen_1, gneg(yQ), vy)));
     273         532 :   tmp2 = gerepileupto(av, tmp2);
     274         532 :   av = avma;
     275         532 :   tmp3 = gsub(gmul(a1, uQ), gmul(gQx, gQy));
     276         532 :   tmp3 = gerepileupto(av, tmp3);
     277             : 
     278         532 :   if (!isog) isog = isog_identity(vx,vy);
     279         532 :   f = gel(isog, 1);
     280         532 :   g = gel(isog, 2);
     281         532 :   h = gel(isog, 3);
     282         532 :   rt_sqr = gsqr(rt);
     283         532 :   res = cgetg(4, t_VEC);
     284         532 :   av = avma;
     285         532 :   tmp4 = gdiv(gadd(gmul(tQ, rt), uQ), rt_sqr);
     286         532 :   gel(res, 1) = gerepileupto(av, gadd(f, tmp4));
     287         532 :   av = avma;
     288         532 :   tmp4 = gadd(tmp1, gmul(rt, gadd(tmp2, tmp3)));
     289         532 :   gel(res, 2) = gerepileupto(av, gsub(g, gdiv(tmp4, gmul(rt, rt_sqr))));
     290         532 :   av = avma;
     291         532 :   gel(res, 3) = gerepileupto(av, gmul(h, rt));
     292         532 :   return gerepileupto(ltop, res);
     293             : }
     294             : 
     295             : /* Given a point P on E, return the curve E/<P> and, if only_image is zero,
     296             :  * the isogeny pi: E -> E/<P>. The variables vx and vy are used to describe
     297             :  * the isogeny (ignored if only_image is zero) */
     298             : static GEN
     299         427 : isogeny_from_kernel_point(GEN E, GEN P, int only_image, long vx, long vy)
     300             : {
     301         427 :   pari_sp av = avma;
     302             :   GEN isog, EE, f, g, h, h2, h3;
     303         427 :   GEN Q = P, t = gen_0, w = gen_0;
     304             :   long c;
     305         427 :   if (!oncurve(E,P))
     306           7 :     pari_err_DOMAIN("isogeny_from_kernel_point", "point", "not on", E, P);
     307         420 :   if (ell_is_inf(P))
     308             :   {
     309          42 :     if (only_image) return E;
     310          28 :     return mkvec2(E, isog_identity(vx,vy));
     311             :   }
     312             : 
     313         378 :   isog = NULL; c = 1;
     314             :   for (;;)
     315             :   {
     316         903 :     GEN tQ, xQ = gel(Q,1), uQ = ec_2divpol_evalx(E, xQ);
     317         903 :     int stop = 0;
     318         903 :     if (ellisweierstrasspoint(E,Q))
     319             :     { /* ord(P)=2c; take Q=[c]P into consideration and stop */
     320         196 :       tQ = ec_dFdx_evalQ(E, Q);
     321         196 :       stop = 1;
     322             :     }
     323             :     else
     324         707 :       tQ = ec_half_deriv_2divpol_evalx(E, xQ);
     325         903 :     t = gadd(t, tQ);
     326         903 :     w = gadd(w, gadd(uQ, gmul(tQ, xQ)));
     327         903 :     if (!only_image) isog = update_isogeny_polys(isog, E, Q,tQ,uQ, vx,vy);
     328         903 :     if (stop) break;
     329             : 
     330         707 :     Q = elladd(E, P, Q);
     331         707 :     ++c;
     332             :     /* IF x([c]P) = x([c-1]P) THEN [c]P = -[c-1]P and [2c-1]P = 0 */
     333         707 :     if (gequal(gel(Q,1), xQ)) break;
     334         525 :     if (gc_needed(av,1))
     335             :     {
     336           0 :       if(DEBUGMEM>1) pari_warn(warnmem,"isogeny_from_kernel_point");
     337           0 :       gerepileall(av, isog? 4: 3, &Q, &t, &w, &isog);
     338             :     }
     339         525 :   }
     340             : 
     341         378 :   EE = make_velu_curve(E, t, w);
     342         378 :   if (only_image) return EE;
     343             : 
     344         224 :   if (!isog) isog = isog_identity(vx,vy);
     345         224 :   f = gel(isog, 1);
     346         224 :   g = gel(isog, 2);
     347         224 :   if ( ! (typ(f) == t_RFRAC && typ(g) == t_RFRAC))
     348           0 :     pari_err_BUG("isogeny_from_kernel_point (f or g has wrong type)");
     349             : 
     350             :   /* Clean up the isogeny polynomials f, g and h so that the isogeny
     351             :    * is given by (x,y) -> (f(x)/h(x)^2, g(x,y)/h(x)^3) */
     352         224 :   h = gel(isog, 3);
     353         224 :   h2 = gsqr(h);
     354         224 :   h3 = gmul(h, h2);
     355         224 :   f = gmul(f, h2);
     356         224 :   g = gmul(g, h3);
     357         224 :   if (typ(f) != t_POL || typ(g) != t_POL)
     358           0 :     pari_err_BUG("isogeny_from_kernel_point (wrong denominator)");
     359         224 :   return mkvec2(EE, mkvec3(f,g, gel(isog,3)));
     360             : }
     361             : 
     362             : /* Given a t_POL x^n - s1 x^{n-1} + s2 x^{n-2} - s3 x^{n-3} + ...
     363             :  * return the first three power sums (Newton's identities):
     364             :  *   p1 = s1
     365             :  *   p2 = s1^2 - 2 s2
     366             :  *   p3 = (s1^2 - 3 s2) s1 + 3 s3 */
     367             : static void
     368        2891 : first_three_power_sums(GEN pol, GEN *p1, GEN *p2, GEN *p3)
     369             : {
     370        2891 :   long d = degpol(pol);
     371             :   GEN s1, s2, ms3;
     372             : 
     373        2891 :   *p1 = s1 = gneg(RgX_coeff(pol, d-1));
     374             : 
     375        2891 :   s2 = RgX_coeff(pol, d-2);
     376        2891 :   *p2 = gsub(gsqr(s1), gmulsg(2L, s2));
     377             : 
     378        2891 :   ms3 = RgX_coeff(pol, d-3);
     379        2891 :   *p3 = gadd(gmul(s1, gsub(*p2, s2)), gmulsg(-3L, ms3));
     380        2891 : }
     381             : 
     382             : 
     383             : /* Let E and a t_POL h of degree 1 or 3 whose roots are 2-torsion points on E.
     384             :  * - if only_image != 0, return [t, w] used to compute the equation of the
     385             :  *   quotient by the given 2-torsion points
     386             :  * - else return [t,w, f,g,h], along with the contributions f, g and
     387             :  *   h to the isogeny giving the quotient by h. Variables vx and vy are used
     388             :  *   to create f, g and h, or ignored if only_image is zero */
     389             : 
     390             : /* deg h = 1; 2-torsion contribution from Weierstrass point */
     391             : static GEN
     392        1232 : contrib_weierstrass_pt(GEN E, GEN h, long only_image, long vx, long vy)
     393             : {
     394        1232 :   GEN p = ellbasechar(E);
     395        1232 :   GEN a1 = ell_get_a1(E);
     396        1232 :   GEN a3 = ell_get_a3(E);
     397        1232 :   GEN x0 = gneg(constant_coeff(h)); /* h = x - x0 */
     398        1232 :   GEN b = gadd(gmul(a1,x0), a3);
     399             :   GEN y0, Q, t, w, t1, t2, f, g;
     400             : 
     401        1232 :   if (!equalis(p, 2L)) /* char(k) != 2 ==> y0 = -b/2 */
     402        1190 :     y0 = gmul2n(gneg(b), -1);
     403             :   else
     404             :   { /* char(k) = 2 ==> y0 = sqrt(f(x0)) where E is y^2 + h(x) = f(x). */
     405          42 :     if (!gequal0(b)) pari_err_BUG("two_torsion_contrib (a1*x0+a3 != 0)");
     406          42 :     y0 = gsqrt(ec_f_evalx(E, x0), 0);
     407             :   }
     408        1232 :   Q = mkvec2(x0, y0);
     409        1232 :   t = ec_dFdx_evalQ(E, Q);
     410        1232 :   w = gmul(x0, t);
     411        1232 :   if (only_image) return mkvec2(t,w);
     412             : 
     413             :   /* Compute isogeny, f = (x - x0) * t */
     414         462 :   f = deg1pol_shallow(t, gmul(t, gneg(x0)), vx);
     415             : 
     416             :   /* g = (x - x0) * t * (a1 * (x - x0) + (y - y0)) */
     417         462 :   t1 = deg1pol_shallow(a1, gmul(a1, gneg(x0)), vx);
     418         462 :   t2 = deg1pol_shallow(gen_1, gneg(y0), vy);
     419         462 :   g = gmul(f, gadd(t1, t2));
     420         462 :   return mkvec5(t, w, f, g, h);
     421             : }
     422             : /* deg h =3; full 2-torsion contribution. NB: assume h is monic; base field
     423             :  * characteristic is odd or zero (cannot happen in char 2).*/
     424             : static GEN
     425          14 : contrib_full_tors(GEN E, GEN h, long only_image, long vx, long vy)
     426             : {
     427             :   GEN p1, p2, p3, half_b2, half_b4, t, w, f, g;
     428          14 :   first_three_power_sums(h, &p1,&p2,&p3);
     429          14 :   half_b2 = gmul2n(ell_get_b2(E), -1);
     430          14 :   half_b4 = gmul2n(ell_get_b4(E), -1);
     431             : 
     432             :   /* t = 3*(p2 + b4/2) + p1 * b2/2 */
     433          14 :   t = gadd(gmulsg(3L, gadd(p2, half_b4)), gmul(p1, half_b2));
     434             : 
     435             :   /* w = 3 * p3 + p2 * b2/2 + p1 * b4/2 */
     436          14 :   w = gadd(gmulsg(3L, p3), gadd(gmul(p2, half_b2),
     437             :                                 gmul(p1, half_b4)));
     438          14 :   if (only_image) return mkvec2(t,w);
     439             : 
     440             :   /* Compute isogeny */
     441             :   {
     442           7 :     GEN a1 = ell_get_a1(E), a3 = ell_get_a3(E), t1, t2;
     443           7 :     GEN s1 = gneg(RgX_coeff(h, 2));
     444           7 :     GEN dh = RgX_deriv(h);
     445           7 :     GEN psi2xy = gadd(deg1pol_shallow(a1, a3, vx),
     446             :                       deg1pol_shallow(gen_2, gen_0, vy));
     447             : 
     448             :     /* f = -3 (3 x + b2/2 + s1) h + (3 x^2 + (b2/2) x + (b4/2)) h'*/
     449           7 :     t1 = RgX_mul(h, gmulsg(-3, deg1pol(stoi(3), gadd(half_b2, s1), vx)));
     450           7 :     t2 = mkpoln(3, stoi(3), half_b2, half_b4);
     451           7 :     setvarn(t2, vx);
     452           7 :     t2 = RgX_mul(dh, t2);
     453           7 :     f = RgX_add(t1, t2);
     454             : 
     455             :     /* 2g = psi2xy * (f'*h - f*h') - (a1*f + a3*h) * h; */
     456           7 :     t1 = RgX_sub(RgX_mul(RgX_deriv(f), h), RgX_mul(f, dh));
     457           7 :     t2 = RgX_mul(h, RgX_add(RgX_Rg_mul(f, a1), RgX_Rg_mul(h, a3)));
     458           7 :     g = RgX_divs(gsub(gmul(psi2xy, t1), t2), 2L);
     459             : 
     460           7 :     f = RgX_mul(f, h);
     461           7 :     g = RgX_mul(g, h);
     462             :   }
     463           7 :   return mkvec5(t, w, f, g, h);
     464             : }
     465             : 
     466             : /* Given E and a t_POL T whose roots define a subgroup G of E, return the factor
     467             :  * of T that corresponds to the 2-torsion points E[2] \cap G in G */
     468             : INLINE GEN
     469        2884 : two_torsion_part(GEN E, GEN T)
     470        2884 : { return RgX_gcd(T, elldivpol(E, 2, varn(T))); }
     471             : 
     472             : /* Return the jth Hasse derivative of the polynomial f = \sum_{i=0}^n a_i x^i,
     473             :  * i.e. \sum_{i=j}^n a_i \binom{i}{j} x^{i-j}. It is a derivation even when the
     474             :  * coefficient ring has positive characteristic */
     475             : static GEN
     476          98 : derivhasse(GEN f, ulong j)
     477             : {
     478          98 :   ulong i, d = degpol(f);
     479             :   GEN df;
     480          98 :   if (gequal0(f) || d == 0) return pol_0(varn(f));
     481          56 :   if (j == 0) return gcopy(f);
     482          56 :   df = cgetg(2 + (d-j+1), t_POL);
     483          56 :   df[1] = f[1];
     484          56 :   for (i = j; i <= d; ++i) gel(df, i-j+2) = gmul(binomialuu(i,j), gel(f, i+2));
     485          56 :   return normalizepol(df);
     486             : }
     487             : 
     488             : static GEN
     489         784 : non_two_torsion_abscissa(GEN E, GEN h0, GEN x)
     490             : {
     491             :   GEN mp1, dh0, ddh0, t, u, t1, t2, t3;
     492         784 :   long m = degpol(h0);
     493         784 :   mp1 = gel(h0, m + 1); /* negative of first power sum */
     494         784 :   dh0 = RgX_deriv(h0);
     495         784 :   ddh0 = RgX_deriv(dh0);
     496         784 :   t = ec_2divpol_evalx(E, x);
     497         784 :   u = ec_half_deriv_2divpol_evalx(E, x);
     498         784 :   t1 = RgX_sub(RgX_sqr(dh0), RgX_mul(ddh0, h0));
     499         784 :   t2 = RgX_mul(u, RgX_mul(h0, dh0));
     500         784 :   t3 = RgX_mul(RgX_sqr(h0),
     501         784 :                deg1pol_shallow(stoi(2*m), gmulsg(2L, mp1), varn(x)));
     502             :   /* t * (dh0^2 - ddh0*h0) - u*dh0*h0 + (2*m*x - 2*s1) * h0^2); */
     503         784 :   return RgX_add(RgX_sub(RgX_mul(t, t1), t2), t3);
     504             : }
     505             : 
     506             : static GEN
     507        1218 : isog_abscissa(GEN E, GEN kerp, GEN h0, GEN x, GEN two_tors)
     508             : {
     509             :   GEN f0, f2, h2, t1, t2, t3;
     510        1218 :   f0 = (degpol(h0) > 0)? non_two_torsion_abscissa(E, h0, x): pol_0(varn(x));
     511        1218 :   f2 = gel(two_tors, 3);
     512        1218 :   h2 = gel(two_tors, 5);
     513             : 
     514             :   /* Combine f0 and f2 into the final abscissa of the isogeny. */
     515        1218 :   t1 = RgX_mul(x, RgX_sqr(kerp));
     516        1218 :   t2 = RgX_mul(f2, RgX_sqr(h0));
     517        1218 :   t3 = RgX_mul(f0, RgX_sqr(h2));
     518             :   /* x * kerp^2 + f2 * h0^2 + f0 * h2^2 */
     519        1218 :   return RgX_add(t1, RgX_add(t2, t3));
     520             : }
     521             : 
     522             : static GEN
     523        1169 : non_two_torsion_ordinate_char_not2(GEN E, GEN f, GEN h, GEN psi2)
     524             : {
     525        1169 :   GEN a1 = ell_get_a1(E), a3 = ell_get_a3(E);
     526        1169 :   GEN df = RgX_deriv(f), dh = RgX_deriv(h);
     527             :   /* g = df * h * psi2/2 - f * dh * psi2
     528             :    *   - (E.a1 * f + E.a3 * h^2) * h/2 */
     529        1169 :   GEN t1 = RgX_mul(df, RgX_mul(h, RgX_divs(psi2, 2L)));
     530        1169 :   GEN t2 = RgX_mul(f, RgX_mul(dh, psi2));
     531        1169 :   GEN t3 = RgX_mul(RgX_divs(h, 2L),
     532             :                    RgX_add(RgX_Rg_mul(f, a1), RgX_Rg_mul(RgX_sqr(h), a3)));
     533        1169 :   return RgX_sub(RgX_sub(t1, t2), t3);
     534             : }
     535             : 
     536             : /* h = kerq */
     537             : static GEN
     538          49 : non_two_torsion_ordinate_char2(GEN E, GEN h, GEN x, GEN y)
     539             : {
     540          49 :   GEN a1 = ell_get_a1(E), a3 = ell_get_a3(E), a4 = ell_get_a4(E);
     541          49 :   GEN b2 = ell_get_b2(E), b4 = ell_get_b4(E), b6 = ell_get_b6(E);
     542             :   GEN h2, dh, dh2, ddh, D2h, D2dh, H, psi2, u, t, alpha;
     543             :   GEN p1, t1, t2, t3, t4;
     544          49 :   long m, vx = varn(x);
     545             : 
     546          49 :   h2 = RgX_sqr(h);
     547          49 :   dh = RgX_deriv(h);
     548          49 :   dh2 = RgX_sqr(dh);
     549          49 :   ddh = RgX_deriv(dh);
     550          49 :   H = RgX_sub(dh2, RgX_mul(h, ddh));
     551          49 :   D2h = derivhasse(h, 2);
     552          49 :   D2dh = derivhasse(dh, 2);
     553          49 :   psi2 = deg1pol_shallow(a1, a3, vx);
     554          49 :   u = mkpoln(3, b2, gen_0, b6);
     555          49 :   setvarn(u, vx);
     556          49 :   t = deg1pol_shallow(b2, b4, vx);
     557          49 :   alpha = mkpoln(4, a1, a3, gmul(a1, a4), gmul(a3, a4));
     558          49 :   setvarn(alpha, vx);
     559          49 :   m = degpol(h);
     560          49 :   p1 = RgX_coeff(h, m-1); /* first power sum */
     561             : 
     562          49 :   t1 = gmul(gadd(gmul(a1, p1), gmulgs(a3, m)), RgX_mul(h,h2));
     563             : 
     564          49 :   t2 = gmul(a1, gadd(gmul(a1, gadd(y, psi2)), RgX_add(RgX_Rg_add(RgX_sqr(x), a4), t)));
     565          49 :   t2 = gmul(t2, gmul(dh, h2));
     566             : 
     567          49 :   t3 = gadd(gmul(y, t), RgX_add(alpha, RgX_Rg_mul(u, a1)));
     568          49 :   t3 = gmul(t3, RgX_mul(h, H));
     569             : 
     570          49 :   t4 = gmul(u, psi2);
     571          49 :   t4 = gmul(t4, RgX_sub(RgX_sub(RgX_mul(h2, D2dh), RgX_mul(dh, H)),
     572             :                         RgX_mul(h, RgX_mul(dh, D2h))));
     573             : 
     574          49 :   return gadd(t1, gadd(t2, gadd(t3, t4)));
     575             : }
     576             : 
     577             : static GEN
     578        1218 : isog_ordinate(GEN E, GEN kerp, GEN kerq, GEN x, GEN y, GEN two_tors, GEN f)
     579             : {
     580             :   GEN g;
     581        1218 :   if (! equalis(ellbasechar(E), 2L)) {
     582             :     /* FIXME: We don't use (hence don't need to calculate)
     583             :      * g2 = gel(two_tors, 4) when char(k) != 2. */
     584        1169 :     GEN psi2 = ec_dmFdy_evalQ(E, mkvec2(x, y));
     585        1169 :     g = non_two_torsion_ordinate_char_not2(E, f, kerp, psi2);
     586             :   } else {
     587          49 :     GEN h2 = gel(two_tors, 5);
     588          49 :     GEN g2 = gmul(gel(two_tors, 4), RgX_mul(kerq, RgX_sqr(kerq)));
     589          49 :     GEN g0 = non_two_torsion_ordinate_char2(E, kerq, x, y);
     590          49 :     g0 = gmul(g0, RgX_mul(h2, RgX_sqr(h2)));
     591          49 :     g = gsub(gmul(y, RgX_mul(kerp, RgX_sqr(kerp))), gadd(g2, g0));
     592             :   }
     593        1218 :   return g;
     594             : }
     595             : 
     596             : /* Given an elliptic curve E and a polynomial kerp whose roots give the
     597             :  * x-coordinates of a subgroup G of E, return the curve E/G and,
     598             :  * if only_image is zero, the isogeny pi:E -> E/G. Variables vx and vy are
     599             :  * used to describe the isogeny (and are ignored if only_image is zero). */
     600             : static GEN
     601        2884 : isogeny_from_kernel_poly(GEN E, GEN kerp, long only_image, long vx, long vy)
     602             : {
     603             :   long m;
     604        2884 :   GEN b2 = ell_get_b2(E), b4 = ell_get_b4(E), b6 = ell_get_b6(E);
     605             :   GEN p1, p2, p3, x, y, f, g, two_tors, EE, t, w;
     606        2884 :   GEN kerh = two_torsion_part(E, kerp);
     607        2884 :   GEN kerq = RgX_divrem(kerp, kerh, ONLY_DIVIDES);
     608        2884 :   if (!kerq) pari_err_BUG("isogeny_from_kernel_poly");
     609             :   /* isogeny degree: 2*degpol(kerp)+1-degpol(kerh) */
     610        2884 :   m = degpol(kerq);
     611             : 
     612        2884 :   kerp = RgX_Rg_div(kerp, leading_coeff(kerp));
     613        2884 :   kerq = RgX_Rg_div(kerq, leading_coeff(kerq));
     614        2884 :   kerh = RgX_Rg_div(kerh, leading_coeff(kerh));
     615        2884 :   switch(degpol(kerh))
     616             :   {
     617             :   case 0:
     618        2380 :     two_tors = only_image? mkvec2(gen_0, gen_0):
     619         749 :       mkvec5(gen_0, gen_0, pol_0(vx), pol_0(vx), pol_1(vx));
     620        1631 :     break;
     621             :   case 1:
     622        1232 :     two_tors = contrib_weierstrass_pt(E, kerh, only_image,vx,vy);
     623        1232 :     break;
     624             :   case 3:
     625          14 :     two_tors = contrib_full_tors(E, kerh, only_image,vx,vy);
     626          14 :     break;
     627             :   default:
     628           7 :     two_tors = NULL;
     629           7 :     pari_err_DOMAIN("isogeny_from_kernel_poly", "kernel polynomial",
     630             :                     "does not define a subgroup of", E, kerp);
     631             :   }
     632        2877 :   first_three_power_sums(kerq,&p1,&p2,&p3);
     633        2877 :   x = pol_x(vx);
     634        2877 :   y = pol_x(vy);
     635             : 
     636             :   /* t = 6 * p2 + b2 * p1 + m * b4, */
     637        2877 :   t = gadd(gmulsg(6L, p2), gadd(gmul(b2, p1), gmulsg(m, b4)));
     638             : 
     639             :   /* w = 10 * p3 + 2 * b2 * p2 + 3 * b4 * p1 + m * b6, */
     640        2877 :   w = gadd(gmulsg(10L, p3),
     641             :            gadd(gmul(gmulsg(2L, b2), p2),
     642             :                 gadd(gmul(gmulsg(3L, b4), p1), gmulsg(m, b6))));
     643             : 
     644        2877 :   EE = make_velu_curve(E, gadd(t, gel(two_tors, 1)),
     645        2877 :                           gadd(w, gel(two_tors, 2)));
     646        2877 :   if (only_image) return EE;
     647             : 
     648        1218 :   f = isog_abscissa(E, kerp, kerq, x, two_tors);
     649        1218 :   g = isog_ordinate(E, kerp, kerq, x, y, two_tors, f);
     650        1218 :   return mkvec2(EE, mkvec3(f,g,kerp));
     651             : }
     652             : 
     653             : /* Given an elliptic curve E and a subgroup G of E, return the curve
     654             :  * E/G and, if only_image is zero, the isogeny corresponding
     655             :  * to the canonical surjection pi:E -> E/G. The variables vx and
     656             :  * vy are used to describe the isogeny (and are ignored if
     657             :  * only_image is zero). The subgroup G may be given either as
     658             :  * a generating point P on E or as a polynomial kerp whose roots are
     659             :  * the x-coordinates of the points in G */
     660             : GEN
     661        1092 : ellisogeny(GEN E, GEN G, long only_image, long vx, long vy)
     662             : {
     663        1092 :   pari_sp av = avma;
     664             :   GEN j, z;
     665        1092 :   checkell(E);j = ell_get_j(E);
     666        1092 :   if (vx < 0) vx = 0;
     667        1092 :   if (vy < 0) vy = 1;
     668        1092 :   if (varncmp(vx, vy) >= 0)
     669           7 :     pari_err_PRIORITY("ellisogeny", pol_x(vx), "<=", vy);
     670        1085 :   if (!only_image && varncmp(vy, gvar(j)) >= 0)
     671           7 :     pari_err_PRIORITY("ellisogeny", j, ">=", vy);
     672        1078 :   switch(typ(G))
     673             :   {
     674             :   case t_VEC:
     675         441 :     checkellpt(G);
     676         441 :     if (!ell_is_inf(G))
     677             :     {
     678         399 :       GEN x =  gel(G,1), y = gel(G,2);
     679         399 :       if (!only_image)
     680             :       {
     681         245 :         if (varncmp(vy, gvar(x)) >= 0)
     682           7 :           pari_err_PRIORITY("ellisogeny", x, ">=", vy);
     683         238 :         if (varncmp(vy, gvar(y)) >= 0)
     684           7 :           pari_err_PRIORITY("ellisogeny", y, ">=", vy);
     685             :       }
     686             :     }
     687         427 :     z = isogeny_from_kernel_point(E, G, only_image, vx, vy);
     688         420 :     break;
     689             :   case t_POL:
     690         630 :     if (!only_image && varncmp(vy, gvar(constant_coeff(G))) >= 0)
     691           7 :       pari_err_PRIORITY("ellisogeny", constant_coeff(G), ">=", vy);
     692         623 :     z = isogeny_from_kernel_poly(E, G, only_image, vx, vy);
     693         616 :     break;
     694             :   default:
     695           7 :     z = NULL;
     696           7 :     pari_err_TYPE("ellisogeny", G);
     697             :   }
     698        1036 :   return gerepilecopy(av, z);
     699             : }
     700             : 
     701             : static GEN
     702        2338 : trivial_isogeny(void)
     703             : {
     704        2338 :   return mkvec3(pol_x(0), scalarpol(pol_x(1), 0), pol_1(0));
     705             : }
     706             : 
     707             : static GEN
     708         245 : isogeny_a4a6(GEN E)
     709             : {
     710         245 :   GEN a1 = ell_get_a1(E), a3 = ell_get_a3(E), b2 = ell_get_b2(E);
     711         245 :   retmkvec3(deg1pol(gen_1, gdivgs(b2, 12), 0),
     712             :             deg1pol(gdivgs(a1,2), deg1pol(gen_1, gdivgs(a3,2), 1), 0),
     713             :             pol_1(0));
     714             : }
     715             : 
     716             : static GEN
     717         245 : invisogeny_a4a6(GEN E)
     718             : {
     719         245 :   GEN a1 = ell_get_a1(E), a3 = ell_get_a3(E), b2 = ell_get_b2(E);
     720         245 :   retmkvec3(deg1pol(gen_1, gdivgs(b2, -12), 0),
     721             :             deg1pol(gdivgs(a1,-2),
     722             :               deg1pol(gen_1, gadd(gdivgs(a3,-2), gdivgs(gmul(b2,a1), 24)), 1), 0),
     723             :             pol_1(0));
     724             : }
     725             : 
     726             : static GEN
     727        1680 : RgXY_eval(GEN P, GEN x, GEN y)
     728             : {
     729        1680 :   return poleval(poleval(P,x), y);
     730             : }
     731             : 
     732             : static GEN
     733         518 : twistisogeny(GEN iso, GEN d)
     734             : {
     735         518 :   GEN d2 = gsqr(d), d3 = gmul(d, d2);
     736         518 :   return mkvec3(gdiv(gel(iso,1), d2), gdiv(gel(iso,2), d3), gel(iso, 3));
     737             : }
     738             : 
     739             : static GEN
     740        1743 : ellisog_by_Kohel(GEN a4, GEN a6, long n, GEN ker, GEN kert, long flag)
     741             : {
     742        1743 :   GEN E = ellinit(mkvec2(a4, a6), NULL, DEFAULTPREC);
     743        1743 :   GEN F = isogeny_from_kernel_poly(E, ker, flag, 0, 1);
     744        1743 :   GEN Et = ellinit(flag ? F: gel(F, 1), NULL, DEFAULTPREC);
     745        1743 :   GEN c4t = ell_get_c4(Et), c6t = ell_get_c6(Et), jt = ell_get_j(Et);
     746        1743 :   if (!flag)
     747             :   {
     748         518 :     GEN Ft = isogeny_from_kernel_poly(Et, kert, flag, 0, 1);
     749         518 :     GEN isot = twistisogeny(gel(Ft, 2), stoi(n));
     750         518 :     return mkvec5(c4t, c6t, jt, gel(F, 2), isot);
     751             :   }
     752        1225 :   else return mkvec3(c4t, c6t, jt);
     753             : }
     754             : 
     755             : static GEN
     756        1575 : ellisog_by_roots(GEN a4, GEN a6, long n, GEN z, long flag)
     757             : {
     758        1575 :   return ellisog_by_Kohel(a4, a6, n, deg1pol(gen_1, gneg(z), 0),
     759             :                                   deg1pol(gen_1, gmulsg(n, z), 0), flag);
     760             : }
     761             : 
     762             : static GEN
     763        2387 : a4a6_divpol(GEN a4, GEN a6, long n)
     764             : {
     765        2387 :   switch(n)
     766             :   {
     767             :     case 2:
     768         910 :       return mkpoln(4, gen_1, gen_0, a4, a6);
     769             :     case 3:
     770        1477 :       return mkpoln(5, utoi(3), gen_0, gmulgs(a4,6) , gmulgs(a6,12),
     771             :                        gneg(gsqr(a4)));
     772             :   }
     773           0 :   return NULL;
     774             : }
     775             : 
     776             : static GEN
     777        2387 : ellisograph_Kohel_iso(GEN nf, GEN e, long n, GEN z, long flag)
     778             : {
     779             :   long i, r;
     780             :   GEN R, V;
     781        2387 :   GEN c4 = gel(e,1), c6 = gel(e, 2);
     782        2387 :   GEN a4 = gdivgs(c4, -48), a6 = gdivgs(c6, -864);
     783        2387 :   GEN P = a4a6_divpol(a4, a6, n);
     784        2387 :   R = nfroots(nf, z ? RgX_div_by_X_x(P, z, NULL): P);
     785        2387 :   r = lg(R);
     786        2387 :   V = cgetg(r, t_VEC);
     787        3962 :   for (i=1; i < r; i++)
     788        1575 :     gel(V, i) = ellisog_by_roots(a4, a6, n, gel(R, i), flag);
     789        2387 :   return mkvec2(V, R);
     790             : }
     791             : 
     792             : static GEN
     793        2331 : ellisograph_Kohel_r(GEN nf, GEN e, long n, GEN z, long flag)
     794             : {
     795        2331 :   GEN W = ellisograph_Kohel_iso(nf, e, n, z, flag);
     796        2331 :   GEN iso = gel(W, 1), R = gel(W, 2);
     797        2331 :   long i, r = lg(iso);
     798        2331 :   GEN V = cgetg(r, t_VEC);
     799        3850 :   for (i=1; i < r; i++)
     800        1519 :     gel(V, i) = ellisograph_Kohel_r(nf, gel(iso, i), n, gmulgs(gel(R, i), -n), flag);
     801        2331 :   return mkvec2(e, V);
     802             : }
     803             : 
     804             : static GEN
     805         336 : corr(GEN c4, GEN c6)
     806             : {
     807         336 :   GEN c62 = gmul2n(c6, 1);
     808         336 :   return gadd(gdiv(gsqr(c4), c62), gdiv(c62, gmulgs(c4,3)));
     809             : }
     810             : 
     811             : static GEN
     812         336 : elkies98(GEN a4, GEN a6, long l, GEN s, GEN a4t, GEN a6t)
     813             : {
     814             :   GEN C, P, S;
     815             :   long i, n, d;
     816         336 :   d = l == 2 ? 1 : l>>1;
     817         336 :   C = cgetg(d+1, t_VEC);
     818         336 :   gel(C, 1) = gdivgs(gsub(a4, a4t), 5);
     819         336 :   if (d >= 2)
     820         336 :     gel(C, 2) = gdivgs(gsub(a6, a6t), 7);
     821         336 :   if (d >= 3)
     822         224 :     gel(C, 3) = gdivgs(gsub(gsqr(gel(C, 1)), gmul(a4, gel(C, 1))), 3);
     823        2870 :   for (n = 3; n < d; ++n)
     824             :   {
     825        2534 :     GEN s = gen_0;
     826       61390 :     for (i = 1; i < n; i++)
     827       58856 :       s = gadd(s, gmul(gel(C, i), gel(C, n-i)));
     828        2534 :     gel(C, n+1) = gdivgs(gsub(gsub(gmulsg(3, s), gmul(gmulsg((2*n-1)*(n-1), a4), gel(C, n-1))), gmul(gmulsg((2*n-2)*(n-2), a6), gel(C, n-2))), (n-1)*(2*n+5));
     829             :   }
     830         336 :   P = cgetg(d+2, t_VEC);
     831         336 :   gel(P, 1 + 0) = stoi(d);
     832         336 :   gel(P, 1 + 1) = s;
     833         336 :   if (d >= 2)
     834         336 :     gel(P, 1 + 2) = gdivgs(gsub(gel(C, 1), gmulgs(gmulsg(2, a4), d)), 6);
     835        3094 :   for (n = 2; n < d; ++n)
     836        2758 :     gel(P, 1 + n+1) = gdivgs(gsub(gsub(gel(C, n), gmul(gmulsg(4*n-2, a4), gel(P, 1+n-1))), gmul(gmulsg(4*n-4, a6), gel(P, 1+n-2))), 4*n+2);
     837         336 :   S = cgetg(d+3, t_POL);
     838         336 :   S[1] = evalsigne(1) | evalvarn(0);
     839         336 :   gel(S, 2 + d - 0) = gen_1;
     840         336 :   gel(S, 2 + d - 1) = gneg(s);
     841        3430 :   for (n = 2; n <= d; ++n)
     842             :   {
     843        3094 :     GEN s = gen_0;
     844       68362 :     for (i = 1; i <= n; ++i)
     845             :     {
     846       65268 :       GEN p = gmul(gel(P, 1+i), gel(S, 2 + d - (n-i)));
     847       65268 :       s = gadd(s, p);
     848             :     }
     849        3094 :     gel(S, 2 + d - n) = gdivgs(s, -n);
     850             :   }
     851         336 :   return S;
     852             : }
     853             : 
     854             : static GEN
     855         490 : ellisog_by_jt(GEN c4, GEN c6, GEN jt, GEN jtp, GEN s0, long n, long flag)
     856             : {
     857         490 :   GEN jtp2 = gsqr(jtp), den = gmul(jt, gsubgs(jt, 1728));
     858         490 :   GEN c4t = gdiv(jtp2, den);
     859         490 :   GEN c6t = gdiv(gmul(jtp, c4t), jt);
     860         490 :   if (flag)
     861         322 :     return mkvec3(c4t, c6t, jt);
     862             :   else
     863             :   {
     864         168 :     GEN co  = corr(c4, c6);
     865         168 :     GEN cot = corr(c4t, c6t);
     866         168 :     GEN s = gmul2n(gmulgs(gadd(gadd(s0, co), gmulgs(cot,-n)), -n), -2);
     867         168 :     GEN a4  = gdivgs(c4, -48), a6 = gdivgs(c6, -864);
     868         168 :     GEN a4t = gmul(gdivgs(c4t, -48), powuu(n,4)), a6t = gmul(gdivgs(c6t, -864), powuu(n,6));
     869         168 :     GEN ker = elkies98(a4, a6, n, s, a4t, a6t);
     870         168 :     GEN st = gmulgs(s, -n);
     871         168 :     GEN a4tt = gmul(a4,powuu(n,4)), a6tt = gmul(a6,powuu(n,6));
     872         168 :     GEN kert = elkies98(a4t, a6t, n, st, a4tt, a6tt);
     873         168 :     return ellisog_by_Kohel(a4, a6, n, ker, kert, flag);
     874             :   }
     875             : }
     876             : 
     877             : /*
     878             : Based on
     879             : RENE SCHOOF
     880             : Counting points on elliptic curves over finite fields
     881             : Journal de Theorie des Nombres de Bordeaux,
     882             : tome 7, no 1 (1995), p. 219-254.
     883             : <http://www.numdam.org/item?id=JTNB_1995__7_1_219_0>
     884             : */
     885             : 
     886             : static GEN
     887         336 : ellisog_by_j(GEN e, GEN jt, long n, GEN P, long flag)
     888             : {
     889         336 :   pari_sp av = avma;
     890         336 :   GEN c4  = gel(e,1), c6 = gel(e, 2), j = gel(e, 3);
     891         336 :   GEN Px = deriv(P, 0), Py = deriv(P, 1);
     892         336 :   GEN Pxj = RgXY_eval(Px, j, jt), Pyj = RgXY_eval(Py, j, jt);
     893         336 :   GEN Pxx  = deriv(Px, 0), Pxy = deriv(Py, 0), Pyy = deriv(Py, 1);
     894         336 :   GEN Pxxj = RgXY_eval(Pxx,j,jt);
     895         336 :   GEN Pxyj = RgXY_eval(Pxy,j,jt);
     896         336 :   GEN Pyyj = RgXY_eval(Pyy,j,jt);
     897         336 :   GEN c6c4 = gdiv(c6, c4);
     898         336 :   GEN jp = gmul(j, c6c4);
     899         336 :   GEN jtp = gdivgs(gmul(jp, gdiv(Pxj, Pyj)), -n);
     900         336 :   GEN jtpn = gmulgs(jtp, n);
     901         336 :   GEN s0 = gdiv(gadd(gadd(gmul(gsqr(jp),Pxxj),gmul(gmul(jp,jtpn),gmul2n(Pxyj,1))),
     902             :                 gmul(gsqr(jtpn),Pyyj)),gmul(jp,Pxj));
     903         336 :   GEN et = ellisog_by_jt(c4, c6, jt, jtp, s0, n, flag);
     904         336 :   return gerepilecopy(av, et);
     905             : }
     906             : 
     907             : static GEN
     908         623 : ellisograph_iso(GEN nf, GEN e, ulong p, GEN P, GEN oj, long flag)
     909             : {
     910             :   long i, r;
     911             :   GEN Pj, R, V;
     912         623 :   GEN j = gel(e, 3);
     913         623 :   Pj = poleval(P, j);
     914         623 :   R = nfroots(nf,oj ? RgX_div_by_X_x(Pj, oj, NULL):Pj);
     915         623 :   r = lg(R);
     916         623 :   V = cgetg(r, t_VEC);
     917         959 :   for (i=1; i < r; i++)
     918         336 :     gel(V, i) = ellisog_by_j(e, gel(R, i), p, P, flag);
     919         623 :   return V;
     920             : }
     921             : 
     922             : static GEN
     923         581 : ellisograph_r(GEN nf, GEN e, ulong p, GEN P, GEN oj, long flag)
     924             : {
     925         581 :   GEN iso = ellisograph_iso(nf, e, p, P, oj, flag);
     926         581 :   GEN j = gel(e, 3);
     927         581 :   long i, r = lg(iso);
     928         581 :   GEN V = cgetg(r, t_VEC);
     929         875 :   for (i=1; i < r; i++)
     930         294 :     gel(V, i) = ellisograph_r(nf, gel(iso, i), p, P, j, flag);
     931         581 :   return mkvec2(e, V);
     932             : }
     933             : 
     934             : static GEN
     935         581 : ellisograph_a4a6(GEN E, long flag)
     936             : {
     937         581 :   GEN c4 = ell_get_c4(E), c6 = ell_get_c6(E), j = ell_get_j(E);
     938         826 :   return flag ? mkvec3(c4, c6, j):
     939         245 :                 mkvec5(c4, c6, j, isogeny_a4a6(E), invisogeny_a4a6(E));
     940             : }
     941             : 
     942             : static GEN
     943         154 : ellisograph_dummy(GEN E, long n, GEN jt, GEN jtt, GEN s0, long flag)
     944             : {
     945         154 :   GEN c4 = ell_get_c4(E), c6 = ell_get_c6(E);
     946         154 :   GEN c6c4 = gdiv(c6, c4);
     947         154 :   GEN jtp = gmul(c6c4, gdivgs(gmul(jt, jtt), -n));
     948         154 :   GEN iso = ellisog_by_jt(c4, c6, jt, jtp, gmul(s0, c6c4), n, flag);
     949         154 :   GEN v = mkvec2(iso, cgetg(1, t_VEC));
     950         154 :   return mkvec2(ellisograph_a4a6(E, flag), mkvec(v));
     951             : }
     952             : 
     953             : static GEN
     954        1099 : isograph_p(GEN nf, GEN e, ulong p, GEN P, long flag)
     955             : {
     956        1099 :   pari_sp av = avma;
     957             :   GEN iso;
     958        1099 :   if (P)
     959         287 :     iso = ellisograph_r(nf, e, p, P, NULL, flag);
     960             :   else
     961         812 :     iso = ellisograph_Kohel_r(nf, e, p, NULL, flag);
     962        1099 :   return gerepilecopy(av, iso);
     963             : }
     964             : 
     965             : static GEN
     966         273 : ellisograph_p(GEN nf, GEN E, ulong p, long flag)
     967             : {
     968         273 :   pari_sp av = avma;
     969         273 :   GEN e = ellisograph_a4a6(E, flag);
     970         273 :   GEN P = p < 5 ? NULL: polmodular_ZXX(p, 0, 0, 1);
     971         273 :   GEN iso = isograph_p(nf, e, p, P, flag);
     972         273 :   return gerepilecopy(av, iso);
     973             : }
     974             : 
     975             : static long
     976        6825 : etree_nbnodes(GEN T)
     977             : {
     978        6825 :   GEN F = gel(T,2);
     979        6825 :   long n = 1, i, l = lg(F);
     980       11053 :   for (i = 1; i < l; i++)
     981        4228 :     n += etree_nbnodes(gel(F, i));
     982        6825 :   return n;
     983             : }
     984             : 
     985             : static long
     986        3052 : etree_listr(GEN nf, GEN T, GEN V, long n, GEN u, GEN ut)
     987             : {
     988        3052 :   GEN E = gel(T, 1), F = gel(T,2);
     989        3052 :   long i, l = lg(F);
     990        3052 :   GEN iso, isot = NULL;
     991        3052 :   if (lg(E) == 6)
     992             :   {
     993         700 :     iso  = ellnfcompisog(nf,gel(E,4), u);
     994         700 :     isot = ellnfcompisog(nf,ut, gel(E,5));
     995         700 :     gel(V, n) = mkvec5(gel(E,1), gel(E,2), gel(E,3), iso, isot);
     996             :   } else
     997             :   {
     998        2352 :     gel(V, n) = mkvec3(gel(E,1), gel(E,2), gel(E,3));
     999        2352 :     iso = u;
    1000             :   }
    1001        4935 :   for (i = 1; i < l; i++)
    1002        1883 :     n = etree_listr(nf, gel(F, i), V, n + 1, iso, isot);
    1003        3052 :   return n;
    1004             : }
    1005             : 
    1006             : static GEN
    1007        1169 : etree_list(GEN nf, GEN T)
    1008             : {
    1009        1169 :   long n = etree_nbnodes(T);
    1010        1169 :   GEN V = cgetg(n+1, t_VEC);
    1011        1169 :   (void) etree_listr(nf, T, V, 1, trivial_isogeny(), trivial_isogeny());
    1012        1169 :   return V;
    1013             : }
    1014             : 
    1015             : static long
    1016        3052 : etree_distmatr(GEN T, GEN M, long n)
    1017             : {
    1018        3052 :   GEN F = gel(T,2);
    1019        3052 :   long i, j, lF = lg(F), m = n + 1;
    1020        3052 :   GEN V = cgetg(lF, t_VECSMALL);
    1021        3052 :   mael(M, n, n) = 0;
    1022        4935 :   for(i = 1; i < lF; i++)
    1023        1883 :     V[i] = m = etree_distmatr(gel(F,i), M, m);
    1024        4935 :   for(i = 1; i < lF; i++)
    1025             :   {
    1026        1883 :     long mi = i==1 ? n+1: V[i-1];
    1027        4270 :     for(j = mi; j < V[i]; j++)
    1028             :     {
    1029        2387 :       mael(M,n,j) = 1 + mael(M, mi, j);
    1030        2387 :       mael(M,j,n) = 1 + mael(M, j, mi);
    1031             :     }
    1032        4522 :     for(j = 1; j < lF; j++)
    1033        2639 :       if (i != j)
    1034             :       {
    1035         756 :         long i1, j1, mj = j==1 ? n+1: V[j-1];
    1036        1624 :         for (i1 = mi; i1 < V[i]; i1++)
    1037        1960 :           for(j1 = mj; j1 < V[j]; j1++)
    1038        1092 :             mael(M,i1,j1) = 2 + mael(M,mj,j1) + mael(M,i1,mi);
    1039             :       }
    1040             :   }
    1041        3052 :   return m;
    1042             : }
    1043             : 
    1044             : static GEN
    1045        1169 : etree_distmat(GEN T)
    1046             : {
    1047        1169 :   long i, n = etree_nbnodes(T);
    1048        1169 :   GEN M = cgetg(n+1, t_MAT);
    1049        4221 :   for(i = 1; i <= n; i++)
    1050        3052 :     gel(M,i) = cgetg(n+1, t_VECSMALL);
    1051        1169 :   (void)etree_distmatr(T, M, 1);
    1052        1169 :   return M;
    1053             : }
    1054             : 
    1055             : static GEN
    1056         364 : list_to_crv(GEN L)
    1057             : {
    1058             :   long i, l;
    1059         364 :   GEN V = cgetg_copy(L, &l);
    1060        1596 :   for(i=1; i < l; i++)
    1061             :   {
    1062        1232 :     GEN Li = gel(L, i);
    1063        1232 :     GEN e = mkvec2(gdivgs(gel(Li,1), -48), gdivgs(gel(Li,2), -864));
    1064        1232 :     gel(V, i) = lg(Li)==6 ? mkvec3(e, gel(Li,4), gel(Li,5)): e;
    1065             :   }
    1066         364 :   return V;
    1067             : }
    1068             : 
    1069             : static GEN
    1070        1169 : distmat_pow(GEN E, ulong p)
    1071             : {
    1072        1169 :   long i, j, n = lg(E)-1;
    1073        1169 :   GEN M = cgetg(n+1, t_MAT);
    1074        4221 :   for(i = 1; i <= n; i++)
    1075             :   {
    1076        3052 :     gel(M,i) = cgetg(n+1, t_COL);
    1077       11970 :     for(j = 1; j <= n; j++)
    1078        8918 :       gmael(M,i,j) = powuu(p,mael(E,i,j));
    1079             :   }
    1080        1169 :   return M;
    1081             : }
    1082             : 
    1083             : /* Assume there is a single p-isogeny */
    1084             : 
    1085             : static GEN
    1086          84 : isomatdbl(GEN nf, GEN L, GEN M, ulong p, GEN T2, long flag)
    1087             : {
    1088          84 :   long i, j, n = lg(L) -1;
    1089          84 :   GEN P = p > 3 ? polmodular_ZXX(p, 0, 0, 1): NULL;
    1090          84 :   GEN V = cgetg(2*n+1, t_VEC);
    1091          84 :   GEN N = cgetg(2*n+1, t_MAT);
    1092         266 :   for(i=1; i <= n; i++)
    1093         182 :     gel(V, i) = gel(L, i);
    1094         266 :   for(i=1; i <= n; i++)
    1095             :   {
    1096         182 :     GEN e = gel(L, i);
    1097             :     GEN F, E;
    1098         182 :     if (i == 1)
    1099          84 :       F = gmael(T2, 2, 1);
    1100             :     else
    1101             :     {
    1102          98 :       if (p > 3)
    1103          42 :         F = ellisograph_iso(nf, e, p, P, NULL, flag);
    1104             :       else
    1105          56 :         F = gel(ellisograph_Kohel_iso(nf, e, p, NULL, flag), 1);
    1106          98 :       if (lg(F) != 2) pari_err_BUG("isomatdbl");
    1107             :     }
    1108         182 :     E = gel(F, 1);
    1109         182 :     if (flag)
    1110          91 :       gel(V, i+n) = mkvec3(gel(E,1), gel(E,2), gel(E,3));
    1111             :     else
    1112             :     {
    1113          91 :       GEN iso = ellcompisog(gel(E,4), gel(e, 4));
    1114          91 :       GEN isot = ellcompisog(gel(e,5), gel(E, 5));
    1115          91 :       gel(V, i+n) = mkvec5(gel(E,1), gel(E,2), gel(E,3), iso, isot);
    1116             :     }
    1117             :   }
    1118         448 :   for(i=1; i <= 2*n; i++)
    1119         364 :     gel(N, i) = cgetg(2*n+1, t_COL);
    1120         266 :   for(i=1; i <= n; i++)
    1121         588 :     for(j=1; j <= n; j++)
    1122             :     {
    1123         406 :       gcoeff(N,i,j) = gcoeff(N,i+n,j+n) = gcoeff(M,i,j);
    1124         406 :       gcoeff(N,i,j+n) = gcoeff(N,i+n,j) = muliu(gcoeff(M,i,j), p);
    1125             :     }
    1126          84 :   return mkvec2(list_to_crv(V), N);
    1127             : }
    1128             : 
    1129             : static ulong
    1130         329 : ellQ_exceptional_iso(GEN j, GEN *jt, GEN *jtp, GEN *s0)
    1131             : {
    1132         329 :   *jt = j; *jtp = gen_1;
    1133         329 :   if (typ(j)==t_INT)
    1134             :   {
    1135         238 :     long js = itos_or_0(j);
    1136             :     GEN j37;
    1137         238 :     if (js==-32768) { *s0 = mkfracss(-1156,539); return 11; }
    1138         224 :     if (js==-121)
    1139          14 :       { *jt = stoi(-24729001) ; *jtp = mkfracss(4973,5633);
    1140          14 :         *s0 = mkfracss(-1961682050,1204555087); return 11;}
    1141         210 :     if (js==-24729001)
    1142          14 :       { *jt = stoi(-121); *jtp = mkfracss(5633,4973);
    1143          14 :         *s0 = mkfracss(-1961682050,1063421347); return 11;}
    1144         196 :     if (js==-884736)
    1145          14 :       { *s0 = mkfracss(-1100,513); return 19; }
    1146         182 :     j37 = negi(uu32toi(37876312,1780746325));
    1147         182 :     if (js==-9317)
    1148             :     {
    1149          14 :       *jt = j37;
    1150          14 :       *jtp = mkfracss(1984136099,496260169);
    1151          14 :       *s0 = mkfrac(negi(uu32toi(457100760,4180820796UL)),
    1152             :                         uu32toi(89049913, 4077411069UL));
    1153          14 :       return 37;
    1154             :     }
    1155         168 :     if (equalii(j, j37))
    1156             :     {
    1157          14 :       *jt = stoi(-9317);
    1158          14 :       *jtp = mkfrac(utoi(496260169),utoi(1984136099UL));
    1159          14 :       *s0 = mkfrac(negi(uu32toi(41554614,2722784052UL)),
    1160             :                         uu32toi(32367030,2614994557UL));
    1161          14 :       return 37;
    1162             :     }
    1163         154 :     if (js==-884736000)
    1164          14 :     { *s0 = mkfracss(-1073708,512001); return 43; }
    1165         140 :     if (equalii(j, negi(powuu(5280,3))))
    1166          14 :     { *s0 = mkfracss(-176993228,85184001); return 67; }
    1167         126 :     if (equalii(j, negi(powuu(640320,3))))
    1168          14 :     { *s0 = mkfrac(negi(uu32toi(72512,1969695276)), uu32toi(35374,1199927297));
    1169          14 :       return 163; }
    1170             :   } else
    1171             :   {
    1172          91 :     GEN j1 = mkfracss(-297756989,2);
    1173          91 :     GEN j2 = mkfracss(-882216989,131072);
    1174          91 :     if (gequal(j, j1))
    1175             :     {
    1176          14 :       *jt = j2; *jtp = mkfracss(1503991,2878441);
    1177          14 :       *s0 = mkfrac(negi(uu32toi(121934,548114672)),uu32toi(77014,117338383));
    1178          14 :       return 17;
    1179             :     }
    1180          77 :     if (gequal(j, j2))
    1181             :     {
    1182          14 :       *jt = j1; *jtp = mkfracss(2878441,1503991);
    1183          14 :       *s0 = mkfrac(negi(uu32toi(121934,548114672)),uu32toi(40239,4202639633UL));
    1184          14 :       return 17;
    1185             :     }
    1186             :   }
    1187         175 :   return 0;
    1188             : }
    1189             : 
    1190             : static GEN
    1191         259 : mkisomat(ulong p, GEN T)
    1192             : {
    1193         259 :   pari_sp av = avma;
    1194         259 :   GEN L = list_to_crv(etree_list(NULL, T));
    1195         259 :   GEN M = distmat_pow(etree_distmat(T), p);
    1196         259 :   return gerepilecopy(av, mkvec2(L, M));
    1197             : }
    1198             : 
    1199             : static GEN
    1200         826 : mkisomatraw(GEN nf, ulong p, GEN T)
    1201             : {
    1202         826 :   pari_sp av = avma;
    1203         826 :   GEN L = etree_list(nf, T);
    1204         826 :   GEN M = distmat_pow(etree_distmat(T), p);
    1205         826 :   return gerepilecopy(av, mkvec2(L, M));
    1206             : }
    1207             : 
    1208             : static GEN
    1209          84 : mkisomatdbl(ulong p, GEN T, ulong p2, GEN T2, long flag)
    1210             : {
    1211          84 :   GEN L = etree_list(NULL, T);
    1212          84 :   GEN M = distmat_pow(etree_distmat(T), p);
    1213          84 :   return isomatdbl(NULL, L, M, p2, T2, flag);
    1214             : }
    1215             : 
    1216             : /*
    1217             : See
    1218             : M.A Kenku
    1219             : On the number of Q-isomorphism classes of elliptic curves in each Q-isogeny class
    1220             : Journal of Number Theory
    1221             : Volume 15, Issue 2, October 1982, Pages 199-202
    1222             : http://www.sciencedirect.com/science/article/pii/0022314X82900257
    1223             : */
    1224             : 
    1225             : enum { _2 = 1, _3 = 2, _5 = 4, _7 = 8, _13 = 16 };
    1226             : static ulong
    1227         175 : ellQ_goodl(GEN E)
    1228             : {
    1229             :   forprime_t T;
    1230         175 :   long i, CM = ellQ_get_CM(E);
    1231         175 :   ulong mask = 31;
    1232         175 :   GEN disc = ell_get_disc(E);
    1233         175 :   pari_sp av = avma;
    1234         175 :   u_forprime_init(&T, 17UL,ULONG_MAX);
    1235        3675 :   for(i=1; mask && i<=20; i++)
    1236             :   {
    1237        3500 :     ulong p = u_forprime_next(&T);
    1238        3500 :     if (umodiu(disc,p)==0) i--;
    1239             :     else
    1240             :     {
    1241        3500 :       long t = ellap_CM_fast(E, p, CM), D = t*t-4*p;
    1242        3500 :       if (t%2) mask &= ~_2;
    1243        3500 :       if ((mask & _3) && kross(D,3)==-1)  mask &= ~_3;
    1244        3500 :       if ((mask & _5) && kross(D,5)==-1)  mask &= ~_5;
    1245        3500 :       if ((mask & _7) && kross(D,7)==-1)  mask &= ~_7;
    1246        3500 :       if ((mask &_13) && kross(D,13)==-1) mask &= ~_13;
    1247             :     }
    1248             :   }
    1249         175 :   avma = av; return mask;
    1250             : }
    1251             : 
    1252             : static long
    1253          14 : ellQ_goodl_l(GEN E, long l)
    1254             : {
    1255             :   forprime_t T;
    1256             :   long i;
    1257          14 :   GEN disc = ell_get_disc(E);
    1258          14 :   pari_sp av = avma;
    1259          14 :   u_forprime_init(&T, 17UL,ULONG_MAX);
    1260         294 :   for(i=1; i<=20; i++)
    1261             :   {
    1262         280 :     ulong p = u_forprime_next(&T);
    1263         280 :     if (umodiu(disc,p)==0) { i--; continue; }
    1264             :     else
    1265             :     {
    1266         280 :       long t = itos(ellap(E, utoi(p)));
    1267         280 :       if (l==2)
    1268             :       {
    1269         140 :         if (t%2==1) return 0;
    1270             :       }
    1271             :       else
    1272             :       {
    1273         140 :         long D = t*t-4*p;
    1274         140 :         if (kross(D,l)==-1) return 0;
    1275             :       }
    1276         280 :       avma = av;
    1277             :     }
    1278             :   }
    1279          14 :   return 1;
    1280             : }
    1281             : 
    1282             : static ulong
    1283          21 : ellnf_goodl_l(GEN E, GEN v)
    1284             : {
    1285             :   forprime_t T;
    1286             :   long i;
    1287          21 :   GEN nf = ellnf_get_nf(E);
    1288          21 :   GEN disc = ell_get_disc(E);
    1289          21 :   long lv = lg(v);
    1290          21 :   ulong w = 0UL;
    1291          21 :   pari_sp av = avma;
    1292          21 :   u_forprime_init(&T, 17UL,ULONG_MAX);
    1293         448 :   for(i=1; i<=20; i++)
    1294             :   {
    1295         427 :     ulong p = u_forprime_next(&T);
    1296         427 :     GEN pr = idealprimedec(nf, utoi(p));
    1297         427 :     long j, k, lv = lg(v), g = lg(pr)-1;
    1298        1064 :     for (j=1; j<=g; j++)
    1299             :     {
    1300         637 :       GEN prj = gel(pr, j);
    1301         637 :       if (idealval(nf,disc,prj) > 0) {i--; continue;}
    1302             :       else
    1303             :       {
    1304         630 :         long t = itos(ellap(E, prj));
    1305        2730 :         for(k = 1; k < lv; k++)
    1306             :         {
    1307        2100 :           long l = v[k];
    1308        2100 :           if (l==2)
    1309             :           {
    1310         630 :             if (t%2==1) w |= 1<<(k-1);
    1311             :           }
    1312             :           else
    1313             :           {
    1314        1470 :             GEN D = subii(sqrs(t),shifti(pr_norm(prj),2));
    1315        1470 :             if (krois(D,l)==-1) w |= 1<<(k-1);
    1316             :           }
    1317             :         }
    1318             :       }
    1319             :     }
    1320         427 :     avma = av;
    1321             :   }
    1322          21 :   return w^((1UL<<(lv-1))-1);
    1323             : }
    1324             : 
    1325             : static GEN
    1326         602 : ellnf_charpoly(GEN E, GEN pr)
    1327             : {
    1328         602 :   return deg2pol_shallow(gen_1, negi(ellap(E,pr)), pr_norm(pr), 0);
    1329             : }
    1330             : 
    1331             : static GEN
    1332        1204 : RgX_homogenize(GEN P, long v)
    1333             : {
    1334        1204 :   GEN Q = leafcopy(P);
    1335        1204 :   long i, l = lg(P), d = degpol(P);
    1336        1204 :   for (i = 2; i < l; i++) gel(Q,i) = monomial(gel(Q,i), d--, v);
    1337        1204 :   return Q;
    1338             : }
    1339             : 
    1340             : static GEN
    1341        1204 : starlaw(GEN p, GEN q)
    1342             : {
    1343        1204 :   GEN Q = RgX_homogenize(RgX_recip(q), 1);
    1344        1204 :   return ZX_ZXY_resultant(p, Q);
    1345             : }
    1346             : 
    1347             : static GEN
    1348         602 : startor(GEN p, long r)
    1349             : {
    1350         602 :   GEN xr = monomial(gen_1, r, 0);
    1351         602 :   GEN psir = gsub(xr, gen_1);
    1352         602 :   return gsubstpol(starlaw(p, psir),xr,pol_x(0));
    1353             : }
    1354             : 
    1355             : static GEN
    1356         420 : ellnf_get_degree(GEN E, GEN p)
    1357             : {
    1358         420 :   GEN nf = ellnf_get_nf(E);
    1359         420 :   long d = nf_get_degree(nf);
    1360         420 :   GEN dec = idealprimedec(nf, p);
    1361         420 :   long i, l = lg(dec), k;
    1362         420 :   GEN R, starl = deg1pol_shallow(gen_1, gen_m1, 0);
    1363        1022 :   for(i=1; i < l; i++)
    1364             :   {
    1365         602 :     GEN pr = gel(dec,i);
    1366         602 :     GEN q = ellnf_charpoly(E, pr);
    1367         602 :     starl = starlaw(starl, startor(q, 12*pr_get_e(pr)));
    1368             :   }
    1369         420 :   R = p;
    1370        1260 :   for(k=0; 2*k<=d; k++)
    1371         840 :     R = mulii(R, poleval(starl,powiu(p,12*k)));
    1372         420 :   return R;
    1373             : }
    1374             : 
    1375             : /*
    1376             : Based on a GP script by Nicolas Billerey itself
    1377             : based on Th\'eor\`emes 2.4 and 2.8 of the following article:
    1378             : N. Billerey, Crit\`eres d'irr\'eductibilit\'e pour les
    1379             : repr\'esentations des courbes elliptiques,
    1380             : Int. J. Number Theory 7 (2011), no. 4, 1001-1032.
    1381             : */
    1382             : 
    1383             : static GEN
    1384          21 : ellnf_prime_degree(GEN E)
    1385             : {
    1386             :   forprime_t T;
    1387             :   long i;
    1388          21 :   GEN nf = ellnf_get_nf(E);
    1389          21 :   GEN disc = ell_get_disc(E);
    1390          21 :   GEN P, B = gen_0, rB;
    1391          21 :   GEN bad = mulii(nfnorm(nf, disc),nf_get_disc(nf));
    1392          21 :   u_forprime_init(&T, 5UL,ULONG_MAX);
    1393         469 :   for(i=1; i<=20; i++)
    1394             :   {
    1395         448 :     ulong p = u_forprime_next(&T);
    1396         448 :     if (dvdiu(bad, p)) {i--; continue;}
    1397         420 :     B = gcdii(B, ellnf_get_degree(E, utoi(p)));
    1398         420 :     if (Z_issquareall(B,&rB)) B=rB;
    1399             :   }
    1400          21 :   if (signe(B)==0) pari_err_IMPL("ellisomat, CM case");
    1401          21 :   P = vec_to_vecsmall(gel(Z_factor(B),1));
    1402          21 :   return shallowextract(P, utoi(ellnf_goodl_l(E, P)));
    1403             : }
    1404             : 
    1405             : static GEN
    1406         329 : ellQ_isomat(GEN E, long flag)
    1407             : {
    1408         329 :   GEN K = NULL;
    1409         329 :   GEN T2 = NULL, T3 = NULL, T5, T7, T13;
    1410             :   ulong good;
    1411             :   long n2, n3, n5, n7, n13;
    1412             :   GEN jt, jtp, s0;
    1413         329 :   GEN c4 = ell_get_c4(E), c6 = ell_get_c6(E), j = ell_get_j(E);
    1414         329 :   long l = ellQ_exceptional_iso(j, &jt, &jtp, &s0);
    1415         329 :   if (l)
    1416             :   {
    1417             : #if 1
    1418         154 :     return mkisomat(l, ellisograph_dummy(E, l, jt, jtp, s0, flag));
    1419             : #else
    1420             :     return mkisomat(l, ellisograph_p(K, E, l), flag);
    1421             : #endif
    1422             :   }
    1423         175 :   good = ellQ_goodl(ellintegralmodel(E,NULL));
    1424         175 :   if (good & _2)
    1425             :   {
    1426         126 :     T2 = ellisograph_p(K, E, 2, flag);
    1427         126 :     n2 = etree_nbnodes(T2);
    1428         126 :     if (n2>4 || gequalgs(j, 1728) || gequalgs(j, 287496))
    1429          42 :       return mkisomat(2, T2);
    1430          49 :   } else n2 = 1;
    1431         133 :   if (good & _3)
    1432             :   {
    1433          84 :     T3 = ellisograph_p(K, E, 3, flag);
    1434          84 :     n3 = etree_nbnodes(T3);
    1435          84 :     if (n3>1 && n2==2) return mkisomatdbl(3,T3,2,T2, flag);
    1436          42 :     if (n3==2 && n2>1)  return mkisomatdbl(2,T2,3,T3, flag);
    1437          42 :     if (n3>2 || gequal0(j)) return mkisomat(3, T3);
    1438          49 :   } else n3 = 1;
    1439          63 :   if (good & _5)
    1440             :   {
    1441          21 :     T5 = ellisograph_p(K, E, 5, flag);
    1442          21 :     n5 = etree_nbnodes(T5);
    1443          21 :     if (n5>1 && n2>1) return mkisomatdbl(2,T2,5,T5, flag);
    1444          21 :     if (n5>1 && n3>1) return mkisomatdbl(3,T3,5,T5, flag);
    1445           7 :     if (n5>1) return mkisomat(5, T5);
    1446          42 :   } else n5 = 1;
    1447          42 :   if (good & _7)
    1448             :   {
    1449          28 :     T7 = ellisograph_p(K, E, 7, flag);
    1450          28 :     n7 = etree_nbnodes(T7);
    1451          28 :     if (n7>1 && n2>1) return mkisomatdbl(2,T2,7,T7, flag);
    1452           0 :     if (n7>1 && n3>1) return mkisomatdbl(3,T3,7,T7, flag);
    1453           0 :     if (n7>1) return mkisomat(7,T7);
    1454          14 :   } else n7 = 1;
    1455          14 :   if (n2>1) return mkisomat(2,T2);
    1456           0 :   if (n3>1) return mkisomat(3,T3);
    1457           0 :   if (good & _13)
    1458             :   {
    1459           0 :     T13 = ellisograph_p(K, E, 13, flag);
    1460           0 :     n13 = etree_nbnodes(T13);
    1461           0 :     if (n13>1) return mkisomat(13,T13);
    1462           0 :   } else n13 = 1;
    1463           0 :   if (flag)
    1464           0 :     retmkvec2(list_to_crv(mkvec(mkvec3(c4, c6, j))), matid(1));
    1465             :   else
    1466           0 :     retmkvec2(list_to_crv(mkvec(mkvec5(c4, c6, j, isogeny_a4a6(E), invisogeny_a4a6(E)))), matid(1));
    1467             : }
    1468             : 
    1469             : static GEN
    1470         154 : ellnf_isocrv(GEN nf, GEN E, GEN v, GEN PE, long flag)
    1471             : {
    1472         154 :   long i, j, l, lv = lg(v);
    1473             :   GEN L, M;
    1474         154 :   GEN e = ellisograph_a4a6(E, flag);
    1475         154 :   GEN LE = cgetg(lv, t_VEC);
    1476         154 :   long n = 1, k = 2;
    1477         504 :   for (i = 1; i < lv; i++)
    1478             :   {
    1479         350 :     ulong p = uel(v,i);
    1480         350 :     GEN P = gel(PE, i);
    1481         350 :     GEN T = isograph_p(nf, e, p, P, flag);
    1482         350 :     GEN LM = mkisomatraw(nf, p, T);
    1483         350 :     gel(LE, i) = LM;
    1484         350 :     n *= lg(gel(LM,1)) - 1;
    1485             :   }
    1486         154 :   L = cgetg(n+1,t_VEC);
    1487         154 :   M = cgetg(n+1,t_COL);
    1488         154 :   gel(L,1) = e;
    1489         154 :   gel(M,1) = gen_1;
    1490         504 :   for (i = 1; i < lv; i++)
    1491             :   {
    1492         350 :     ulong p = uel(v, i);
    1493         350 :     GEN P = gel(PE, i);
    1494         350 :     GEN Li = gmael(LE, i, 1);
    1495         350 :     GEN Mi = gmael(LE, i, 2);
    1496         350 :     long kk = k;
    1497         350 :     long m = lg(Li);
    1498         952 :     for (j = 2; j < m; j++)
    1499             :     {
    1500         602 :       gel(L, k) = gel(Li, j);
    1501         602 :       gel(M, k) = gmael(Mi,1,j);
    1502         602 :       k++;
    1503             :     }
    1504         826 :     for (l = 2; l < kk; l++)
    1505             :     {
    1506         476 :       GEN e = gel(L, l);
    1507         476 :       GEN T = isograph_p(nf, e, p, P, flag);
    1508         476 :       GEN LMe = mkisomatraw(nf, p, T);
    1509         476 :       GEN Le = gel(LMe, 1);
    1510         476 :       GEN Me = gel(LMe, 2);
    1511         476 :       long m = lg(Le);
    1512        1204 :       for (j = 2; j < m; j++)
    1513             :       {
    1514         728 :         gel(L,k) = gel(Le, j);
    1515         728 :         gel(M,k) = gmul(gmael(Me,1,j),gel(M,l));
    1516         728 :         k++;
    1517             :       }
    1518             :     }
    1519             :   }
    1520         154 :   return mkvec2(L, M);
    1521             : }
    1522             : 
    1523             : static long
    1524        1526 : nfispower(GEN nf, long d, GEN a, GEN b)
    1525             : {
    1526             :   GEN N;
    1527        1526 :   if (gequal(a,b)) return 1;
    1528        1120 :   N = nfroots(nf, gsub(monomial(b, d, 0), monomial(a,0,0)));
    1529        1120 :   return lg(N) > 1;
    1530             : }
    1531             : 
    1532             : static long
    1533        7791 : isomat_eq(GEN nf, GEN e1, GEN e2)
    1534             : {
    1535        7791 :   if (gequal(e1,e2)) return 1;
    1536        7791 :   if (!gequal(gel(e1,3), gel(e2,3))) return 0;
    1537        1526 :   if (gequal0(gel(e1,3)))
    1538           0 :     return nfispower(nf,6,gel(e1,2),gel(e2,2));
    1539        1526 :   if (gequalgs(gel(e1,3),1728))
    1540           0 :     return nfispower(nf,4,gel(e1,1),gel(e2,1));
    1541        1526 :   return nfispower(nf,2,gmul(gel(e1,1),gel(e2,2)),gmul(gel(e1,2),gel(e2,1)));
    1542             : }
    1543             : 
    1544             : static long
    1545        1330 : isomat_find(GEN nf, GEN e, GEN L)
    1546             : {
    1547        1330 :   long i, l = lg(L);
    1548        7791 :   for (i=1; i<l; i++)
    1549        7791 :     if (isomat_eq(nf, e, gel(L,i))) return i;
    1550             :   pari_err_BUG("isomat_find"); return 0; /* LCOV_EXCL_LINE */
    1551             : }
    1552             : 
    1553             : static GEN
    1554         133 : isomat_perm(GEN nf, GEN E, GEN L)
    1555             : {
    1556         133 :   long i, l = lg(E);
    1557         133 :   GEN v = cgetg(l, t_VECSMALL);
    1558        1463 :   for (i=1; i<l; i++)
    1559        1330 :     uel(v, i) = isomat_find(nf, gel(E,i), L);
    1560         133 :   return v;
    1561             : }
    1562             : 
    1563             : static GEN
    1564          21 : ellnf_modpoly(GEN v)
    1565             : {
    1566          21 :   long i, l = lg(v);
    1567          21 :   GEN P = cgetg(l, t_VEC);
    1568          63 :   for(i = 1; i < l; i++)
    1569             :   {
    1570          42 :     ulong p = uel(v, i);
    1571          42 :     gel(P, i) = p<5 ? NULL: polmodular_ZXX(p, 0, 0, 1);
    1572             :   }
    1573          21 :   return P;
    1574             : }
    1575             : 
    1576             : static GEN
    1577          21 : ellnf_isomat(GEN E, long flag)
    1578             : {
    1579          21 :   GEN nf = ellnf_get_nf(E);
    1580          21 :   GEN v = ellnf_prime_degree(E);
    1581          21 :   GEN P = ellnf_modpoly(v);
    1582          21 :   GEN LM = ellnf_isocrv(nf, E, v, P, flag);
    1583          21 :   GEN L = gel(LM,1), M = gel(LM, 2);
    1584          21 :   long i, l = lg(L);
    1585             :   GEN R;
    1586          21 :   R = cgetg(l, t_MAT);
    1587          21 :   gel(R, 1) = M;
    1588         154 :   for(i = 2; i < l; i++)
    1589             :   {
    1590         133 :     GEN Li = gel(L, i);
    1591         133 :     GEN e = mkvec2(gdivgs(gel(Li,1), -48), gdivgs(gel(Li,2), -864));
    1592         133 :     GEN LMi = ellnf_isocrv(nf, ellinit(e, nf, DEFAULTPREC), v, P, 1);
    1593         133 :     GEN LLi = gel(LMi, 1), Mi = gel(LMi, 2);
    1594         133 :     GEN r = isomat_perm(nf, L, LLi);
    1595         133 :     gel(R, i) = vecpermute(Mi, r);
    1596             :   }
    1597          21 :   return mkvec2(list_to_crv(L), R);
    1598             : }
    1599             : 
    1600             : GEN
    1601         364 : ellisomat(GEN E, long p, long flag)
    1602             : {
    1603         364 :   pari_sp av = avma;
    1604             :   GEN r;
    1605         364 :   if (flag < 0 || flag > 1) pari_err_FLAG("ellisomat");
    1606         364 :   if (p < 0) pari_err_PRIME("ellisomat", utoi(p));
    1607         364 :   if (p == 1) { flag = 1; p = 0; } /* for backward compatibility */
    1608         364 :   checkell(E);
    1609         364 :   switch(ell_get_type(E))
    1610             :   {
    1611           0 :     default: pari_err_TYPE("ellisomat",E);
    1612           0 :       return NULL; /* NOT REACHED */
    1613             :     case t_ELL_Q:
    1614         343 :       if (!p) r = ellQ_isomat(E, flag);
    1615             :       else
    1616             :       {
    1617          14 :         if (ellQ_goodl_l(E, p))
    1618          14 :           r = mkisomat(p, ellisograph_p(NULL, E, p, flag));
    1619           0 :         else r = mkvec2(mkvec(ellisograph_a4a6(E, flag)),matid(1));
    1620             :       }
    1621         343 :       break;
    1622             :     case t_ELL_NF:
    1623          21 :       if (!p)
    1624          21 :         r = ellnf_isomat(E, flag);
    1625             :       else
    1626             :       {
    1627           0 :         if (ellnf_goodl_l(E, mkvecsmall(p)))
    1628           0 :           r = mkisomat(p, ellisograph_p(ellnf_get_nf(E), E, p, flag));
    1629             :         else
    1630           0 :           r = mkvec2(mkvec(ellisograph_a4a6(E, flag)),matid(1));
    1631             :       }
    1632             :   }
    1633         364 :   return gerepilecopy(av, r);
    1634             : }

Generated by: LCOV version 1.11