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

Generated by: LCOV version 1.13