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

Generated by: LCOV version 1.13