Path: blob/master/src/sage/schemes/elliptic_curves/isogeny_small_degree.py
8820 views
r"""1Isogenies of small prime degree.23Functions for the computation of isogenies of small primes4degree. First: `l` = 2, 3, 5, 7, or 13, where the modular curve5`X_0(l)` has genus 0. Second: `l` = 11, 17, 19, 23, 29, 31, 41, 47,659, or 71, where `X_0^+(l)` has genus 0 and `X_0(l)` is elliptic or7hyperelliptic. Also: `l` = 11, 17, 19, 37, 43, 67 or 163 over `\QQ`8(the sporadic cases with only finitely many `j`-invariants each). All9the above only require factorization of a polynomial of degree `l+1`.10Finally, a generic function which works for arbitrary odd primes `l`11(including the characteristic), but requires factorization of the12`l`-division polynomial, of degree `(l^2-1)/2`.131415AUTHORS:1617- John Cremona and Jenny Cooley: 2009-07..11: the genus 0 cases the sporadic cases over `\QQ`.1819- Kimi Tsukazaki and John Cremona: 2013-07: The 10 (hyper)-elliptic20cases and the generic algorithm. See [KT2013]_.2122REFERENCES:2324.. [CW2005] J. E. Cremona and M. Watkins. Computing isogenies of elliptic curves. preprint, 2005.25.. [KT2013] K. Tsukazaki, Explicit Isogenies of Elliptic Curves,26PhD thesis, University of Warwick, 2013.272829"""3031#*****************************************************************************32# Copyright (C) 2012-2013 John Cremona, Jenny Cooley, Kimi Tsukazaki33#34# Distributed under the terms of the GNU General Public License (GPL)35# as published by the Free Software Foundation; either version 2 of36# the License, or (at your option) any later version.37# http://www.gnu.org/licenses/38#*****************************************************************************3940from sage.categories import homset4142from sage.rings.polynomial.polynomial_ring_constructor import PolynomialRing43from sage.rings.polynomial.polynomial_ring import polygen44from sage.rings.all import Integer, ZZ, QQ45from sage.schemes.elliptic_curves.all import EllipticCurve4647from sage.misc.cachefunc import cached_function4849##########################################################################50# The following section is all about computing l-isogenies, where l is51# a prime. The genus 0 cases `l` = 2, 3, 5, 7 and 13 are52# implemented over any field of characteristic not 2, 3 or `l`; over53# `\QQ` the "sporadic" cases `l` = 11, 17, 19, 37, 43, 67 or 163 with54# only finitely many `j`-invariants each. are also implemented.55##########################################################################5657@cached_function58def Fricke_polynomial(l):59r"""60Fricke polynomial for ``l`` =2,3,5,7,13.6162For these primes (and these only) the modular curve `X_0(l)` has63genus zero, and its field is generated by a single modular64function called the Fricke module (or Hauptmodul), `t`. There is65a classical choice of such a generator `t` in each case, and the66`j`-function is a rational function of `t` of degree `l+1` of the67form `P(t)/t` where `P` is a polynomial of degree `l+1`. Up to68scaling, `t` is determined by the condition that the ramification69points above `j=\infty` are `t=0` (with ramification degree `1`)70and `t=\infty` (with degree `l`). The ramification above `j=0`71and `j=1728` may be seen in the factorizations of `j(t)` and72`k(t)` where `k=j-1728`.7374OUTPUT:7576The polynomial `P(t)` as an element of `\ZZ[t]`.7778TESTS::7980sage: from sage.schemes.elliptic_curves.isogeny_small_degree import Fricke_polynomial81sage: Fricke_polynomial(2)82t^3 + 48*t^2 + 768*t + 409683sage: Fricke_polynomial(3)84t^4 + 36*t^3 + 270*t^2 + 756*t + 72985sage: Fricke_polynomial(5)86t^6 + 30*t^5 + 315*t^4 + 1300*t^3 + 1575*t^2 + 750*t + 12587sage: Fricke_polynomial(7)88t^8 + 28*t^7 + 322*t^6 + 1904*t^5 + 5915*t^4 + 8624*t^3 + 4018*t^2 + 748*t + 4989sage: Fricke_polynomial(13)90t^14 + 26*t^13 + 325*t^12 + 2548*t^11 + 13832*t^10 + 54340*t^9 + 157118*t^8 + 333580*t^7 + 509366*t^6 + 534820*t^5 + 354536*t^4 + 124852*t^3 + 15145*t^2 + 746*t + 1391"""92Zt = PolynomialRing(ZZ,'t')93t = Zt.gen()94if l==2: return (t+16)**395elif l==3: return (t+3)**3*(t+27)96elif l==5: return (t**2+10*t+5)**397elif l==7: return (t**2+5*t+1)**3 * (t**2+13*t+49)98elif l==13: return (t**2+5*t+13)*(t**4+7*t**3+20*t**2+19*t+1)**399else:100raise ValueError("The only genus zero primes are 2, 3, 5, 7 or 13.")101102@cached_function103def Fricke_module(l):104r"""105Fricke module for ``l`` =2,3,5,7,13.106107For these primes (and these only) the modular curve `X_0(l)` has108genus zero, and its field is generated by a single modular109function called the Fricke module (or Hauptmodul), `t`. There is110a classical choice of such a generator `t` in each case, and the111`j`-function is a rational function of `t` of degree `l+1` of the112form `P(t)/t` where `P` is a polynomial of degree `l+1`. Up to113scaling, `t` is determined by the condition that the ramification114points above `j=\infty` are `t=0` (with ramification degree `1`)115and `t=\infty` (with degree `l`). The ramification above `j=0`116and `j=1728` may be seen in the factorizations of `j(t)` and117`k(t)` where `k=j-1728`.118119OUTPUT:120121The rational function `P(t)/t`.122123TESTS::124125sage: from sage.schemes.elliptic_curves.isogeny_small_degree import Fricke_module126sage: Fricke_module(2)127(t^3 + 48*t^2 + 768*t + 4096)/t128sage: Fricke_module(3)129(t^4 + 36*t^3 + 270*t^2 + 756*t + 729)/t130sage: Fricke_module(5)131(t^6 + 30*t^5 + 315*t^4 + 1300*t^3 + 1575*t^2 + 750*t + 125)/t132sage: Fricke_module(7)133(t^8 + 28*t^7 + 322*t^6 + 1904*t^5 + 5915*t^4 + 8624*t^3 + 4018*t^2 + 748*t + 49)/t134sage: Fricke_module(13)135(t^14 + 26*t^13 + 325*t^12 + 2548*t^11 + 13832*t^10 + 54340*t^9 + 157118*t^8 + 333580*t^7 + 509366*t^6 + 534820*t^5 + 354536*t^4 + 124852*t^3 + 15145*t^2 + 746*t + 13)/t136"""137try:138t = PolynomialRing(QQ,'t').gen()139return Fricke_polynomial(l) / t140except ValueError:141raise ValueError("The only genus zero primes are 2, 3, 5, 7 or 13.")142143@cached_function144def Psi(l, use_stored=True):145r"""146Generic kernel polynomial for genus zero primes.147148For each of the primes `l` for which `X_0(l)` has genus zero149(namely `l=2,3,5,7,13`), we may define an elliptic curve `E_t`150over `\QQ(t)`, with coefficients in `\ZZ[t]`, which has good151reduction except at `t=0` and `t=\infty` (which lie above152`j=\infty`) and at certain other values of `t` above `j=0` when153`l=3` (one value) or `l\equiv1\pmod{3}` (two values) and above154`j=1728` when `l=2` (one value) or `l\equiv1 \pmod{4}` (two155values). (These exceptional values correspond to endomorphisms of156`E_t` of degree `l`.) The `l`-division polynomial of `E_t` has a157unique factor of degree `(l-1)/2` (or 1 when `l=2`), with158coefficients in `\ZZ[t]`, which we call the Generic Kernel159Polynomial for `l`. These are used, by specialising `t`, in the160function :meth:`isogenies_prime_degree_genus_0`, which also has to161take into account the twisting factor between `E_t` for a specific162value of `t` and the short Weierstrass form of an elliptic curve163with `j`-invariant `j(t)`. This enables the computation of the164kernel polynomials of isogenies without having to compute and165factor division polynomials.166167All of this data is quickly computed from the Fricke modules, except168that for `l=13` the factorization of the Generic Division Polynomial169takes a long time, so the value have been precomputed and cached; by170default the cached values are used, but the code here will recompute171them when ``use_stored`` is ``False``, as in the doctests.172173INPUT:174175- ``l`` -- either 2, 3, 5, 7, or 13.176177- ``use_stored`` (boolean, default True) -- If True, use178precomputed values, otherwise compute them on the fly.179180.. note:181182This computation takes a negligible time for `l=2,3,5,7`183but more than 100s for `l=13`. The reason184for allowing dynamic computation here instead of just using185precomputed values is for testing.186187TESTS::188189sage: from sage.schemes.elliptic_curves.isogeny_small_degree import Fricke_module, Psi190sage: assert Psi(2, use_stored=True) == Psi(2, use_stored=False)191sage: assert Psi(3, use_stored=True) == Psi(3, use_stored=False)192sage: assert Psi(5, use_stored=True) == Psi(5, use_stored=False)193sage: assert Psi(7, use_stored=True) == Psi(7, use_stored=False)194sage: assert Psi(13, use_stored=True) == Psi(13, use_stored=False) # not tested (very long time)195"""196if not l in [2,3,5,7,13]:197raise ValueError("Genus zero primes are 2, 3, 5, 7 or 13.")198199R = PolynomialRing(ZZ,2,'Xt')200X,t = R.gens()201202if use_stored:203if l==2:204return X + t + 64205if l==3:206return X + t + 27207if l==5:208return X**2 + 2*X*(t**2 + 22*t + 125)+ (t**2 + 22*t + 89) * (t**2 + 22*t + 125)209if l==7:210return (X**3 + 3*(t**2 + 13*t + 49)*X**2211+ 3*(t**2 + 13*t + 33)*(t**2 + 13*t + 49)*X212+ (t**2 + 13*t + 49)*(t**4 + 26*t**3 + 219*t**2 + 778*t + 881))213if l==13:214return (t**24 + 66*t**23 + 2091*t**22 + 6*X*t**20 + 42582*t**21 + 330*X*t**19 + 627603*t**20 + 8700*X*t**18 + 7134744*t**19 + 15*X**2*t**16 + 146886*X*t**17 + 65042724*t**18 + 660*X**2*t**15 + 1784532*X*t**16 + 487778988*t**17 + 13890*X**2*t**14 + 16594230*X*t**15 + 3061861065*t**16 + 20*X**3*t**12 + 186024*X**2*t**13 + 122552328*X*t**14 + 16280123754*t**15 + 660*X**3*t**11 + 1774887*X**2*t**12 + 735836862*X*t**13 + 73911331425*t**14 + 10380*X**3*t**10 + 12787272*X**2*t**11 + 3646188342*X*t**12 + 287938949178*t**13 + 15*X**4*t**8 + 102576*X**3*t**9 + 71909658*X**2*t**10 + 15047141292*X*t**11 + 964903805434*t**12 + 330*X**4*t**7 + 707604*X**3*t**8 + 321704316*X**2*t**9 + 51955096824*X*t**10 + 2781843718722*t**11 + 3435*X**4*t**6 + 3582876*X**3*t**7 + 1155971196*X**2*t**8 + 150205315932*X*t**9 + 6885805359741*t**10 + 6*X**5*t**4 + 21714*X**4*t**5 + 13632168*X**3*t**6 + 3343499244*X**2*t**7 + 362526695094*X*t**8 + 14569390179114*t**9 + 66*X**5*t**3 + 90660*X**4*t**4 + 39215388*X**3*t**5 + 7747596090*X**2*t**6 + 725403501318*X*t**7 + 26165223178293*t**8 + 336*X**5*t**2 + 255090*X**4*t**3 + 84525732*X**3*t**4 + 14206132008*X**2*t**5 + 1189398495432*X*t**6 + 39474479008356*t**7 + X**6 + 858*X**5*t + 472143*X**4*t**2 + 132886992*X**3*t**3 + 20157510639*X**2*t**4 + 1569568001646*X*t**5 + 49303015587132*t**6 + 1014*X**5 + 525954*X**4*t + 144222780*X**3*t**2 + 21320908440*X**2*t**3 + 1622460290100*X*t**4 + 49941619724976*t**5 + 272259*X**4 + 96482100*X**3*t + 15765293778*X**2*t**2 + 1260038295438*X*t**3 + 39836631701295*t**4 + 29641924*X**3 + 7210949460*X**2*t + 686651250012*X*t**2 + 23947528862166*t**3 + 1506392823*X**2 + 231462513906*X*t + 10114876838391*t**2 + 35655266790*X + 2644809206442*t + 317295487717)215# The coefficients for l=13 are:216# X**6: 1217# X**5: (6) * (t**2 + 5*t + 13) * (t**2 + 6*t + 13)218# X**4: (3) * (t**2 + 5*t + 13) * (t**2 + 6*t + 13) * (5*t**4 + 55*t**3 + 260*t**2 + 583*t + 537)219# X**3: (4) * (t**2 + 5*t + 13) * (t**2 + 6*t + 13)**2 * (5*t**6 + 80*t**5 + 560*t**4 + 2214*t**3 + 5128*t**2 + 6568*t + 3373)220# X**2: (3) * (t**2 + 5*t + 13)**2 * (t**2 + 6*t + 13)**2 * (5*t**8 + 110*t**7 + 1045*t**6 + 5798*t**5 + 20508*t**4 + 47134*t**3 + 67685*t**2 + 54406*t + 17581)221# X**1: (6) * (t**2 + 5*t + 13)**2 * (t**2 + 6*t + 13)**3 * (t**10 + 27*t**9 + 316*t**8 + 2225*t**7 + 10463*t**6 + 34232*t**5 + 78299*t**4 + 122305*t**3 + 122892*t**2 + 69427*t + 16005)222# X**0: (t**2 + 5*t + 13)**2 * (t**2 + 6*t + 13)**3 * (t**14 + 38*t**13 + 649*t**12 + 6844*t**11 + 50216*t**10 + 271612*t**9 + 1115174*t**8 + 3520132*t**7 + 8549270*t**6 + 15812476*t**5 + 21764840*t**4 + 21384124*t**3 + 13952929*t**2 + 5282630*t + 854569)223#224225# Here the generic kernel polynomials are actually calculated:226j = Fricke_module(l)227k = j-1728228from sage.misc.all import prod229f = prod( [p for p,e in j.factor() if e==3]230+[p for p,e in k.factor() if e==2])231A4 = -3*t**2*j*k // f**2232A6 = -2*t**3*j*k**2 // f**3233E = EllipticCurve([0,0,0,A4,A6])234assert E.j_invariant() == j235return E.division_polynomial(l,X).factor()[0][0]236237238def isogenies_prime_degree_genus_0(E, l=None):239"""240Returns list of ``l`` -isogenies with domain ``E``.241242INPUT:243244- ``E`` -- an elliptic curve.245246- ``l`` -- either None or 2, 3, 5, 7, or 13.247248OUTPUT:249250(list) When ``l`` is None a list of all isogenies of degree 2, 3,2515, 7 and 13, otherwise a list of isogenies of the given degree.252253.. note::254255This function would normally be invoked indirectly via256``E.isogenies_prime_degree(l)``, which automatically calls the257appropriate function.258259ALGORITHM:260261Cremona and Watkins [CW2005]_. See also [KT2013]_, Chapter 4.262263EXAMPLES::264265sage: from sage.schemes.elliptic_curves.isogeny_small_degree import isogenies_prime_degree_genus_0266sage: E = EllipticCurve([0,12])267sage: isogenies_prime_degree_genus_0(E, 5)268[]269270sage: E = EllipticCurve('1450c1')271sage: isogenies_prime_degree_genus_0(E)272[Isogeny of degree 3 from Elliptic Curve defined by y^2 + x*y = x^3 + x^2 + 300*x - 1000 over Rational Field to Elliptic Curve defined by y^2 + x*y = x^3 + x^2 - 5950*x - 182250 over Rational Field]273274sage: E = EllipticCurve('50a1')275sage: isogenies_prime_degree_genus_0(E)276[Isogeny of degree 3 from Elliptic Curve defined by y^2 + x*y + y = x^3 - x - 2 over Rational Field to Elliptic Curve defined by y^2 + x*y + y = x^3 - 126*x - 552 over Rational Field,277Isogeny of degree 5 from Elliptic Curve defined by y^2 + x*y + y = x^3 - x - 2 over Rational Field to Elliptic Curve defined by y^2 + x*y + y = x^3 - 76*x + 298 over Rational Field]278"""279if not l in [2, 3, 5, 7, 13, None]:280raise ValueError("%s is not a genus 0 prime."%l)281F = E.base_ring()282j = E.j_invariant()283if F.characteristic() in [2, 3, l]:284raise NotImplementedError("2, 3, 5, 7 and 13-isogenies are not yet implemented in characteristic 2 and 3, and when the characteristic is the same as the degree of the isogeny.")285if l==2:286return isogenies_2(E)287if l==3:288return isogenies_3(E)289if j==F(0):290if l==5:291return isogenies_5_0(E)292if l==7:293return isogenies_7_0(E)294if l==13:295return isogenies_13_0(E)296if j==F(1728):297if l==5:298return isogenies_5_1728(E)299if l==7:300return isogenies_7_1728(E)301if l==13:302return isogenies_13_1728(E)303304if l != None:305R = PolynomialRing(F,'t')306t = R.gen()307f = R(Fricke_polynomial(l))308t_list = (f-j*t).roots(multiplicities=False)309t_list.sort()310# The generic kernel polynomial applies to a standard curve311# E_t with the correct j-invariant; we must compute the312# appropriate twising factor to scale X by:313c4, c6 = E.c_invariants()314T = c4/(3*c6)315jt = Fricke_module(l)316kt = jt-1728317from sage.misc.all import prod318psi = Psi(l)319X = t320f = R(prod( [p for p,e in jt.factor() if e==3]321+[p for p,e in kt.factor() if e==2]))322kernels = [R(psi(X*T*(j-1728)*t0/f(t0),t0)) for t0 in t_list]323kernels = [ker.monic() for ker in kernels]324E1 = EllipticCurve([-27*c4,-54*c6])325w = E.isomorphism_to(E1)326model = "minimal" if F is QQ else None327isogs = [E1.isogeny(kernel=ker, model=model) for ker in kernels]328[isog.set_pre_isomorphism(w) for isog in isogs]329return isogs330331if l == None:332return sum([isogenies_prime_degree_genus_0(E, l) for l in [2,3,5,7,13]],[])333334335# The following code computes data to be used in336# isogenies_sporadic_Q. Over Q there are only finitely many337# j-invariants of curves with l-isogenies where l is not equal to 2,338# 3, 5, 7 or 13. In these cases l is equal to 11, 17, 19, 37, 43, 67339# or 163. We refer to these l as "sporadic".340341# sporadic_j is a dictionary holding for each possible sporadic342# j-invariant, the unique l such that an l-isogeny exists.343sporadic_j = {344QQ(-121) : 11,345QQ(-32768) : 11,346QQ(-24729001) : 11,347QQ(-297756989)/2 : 17,348QQ(-882216989)/131072 : 17,349QQ(-884736) : 19,350QQ(-9317) : 37,351QQ(-162677523113838677) : 37,352QQ(-884736000) : 43,353QQ(-147197952000) : 67,354QQ(-262537412640768000) : 163355}356357@cached_function358def _sporadic_Q_data(j):359"""360Returns technical data used in computing sporadic isogenies over `\QQ`.361362INPUT:363364- ``j`` -- The `j`-invariant of a sporadic curve, i.e. one of the365keys of ``sporadic_j``.366367OUTPUT:368369``([a4,a6],coeffs)`` where ``[a4,a6]`` are the coefficients of a370short Weierstrass equation of an elliptic curve E with j(E)=``j``,371and ``coeffs`` is a list of coefficients of a polynomial defining372the kernel of an l-isogeny from E.373374Whenever we have a curve of j-invariant ``j``, we can compute the375corresponding l-isogeny by just scaling ``coeffs`` by the right376twisting factor and using the result as a kernel-polynomial.377378ALGORITHM:379380For small l it works fine to factor the l-division polynomial, but381this takes a long time for the larger l and is a very bad idea for382l=163; hence we use floating point arithmetic with a precision383which is known to work. This idea was suggested by Samir Siksek.384385TESTS::386387sage: from sage.schemes.elliptic_curves.isogeny_small_degree import sporadic_j, _sporadic_Q_data388sage: [_sporadic_Q_data(j) for j in sorted(sporadic_j.keys()) if j != -262537412640768000]389[([-269675595, -1704553285050],390[-855506888466179262477032094260950275409164148942611063430052125977143159,391-1469048260972089939455942042937882262144594798448952781325533511718750,392-1171741935131505774747142644126089902595908234671576131857702734375,393-574934780393177024547076427530739751753985644656221274606250000,394-193516922725803688001809624711400287605136013195315374687500,395-47085563820928456130325308223963045033502182349693125000,396-8472233937388712980597845725196873697064639957437500,397-1124815211213953261752081095348112305023653750000,398-105684015609077608033913080859605951322531250,399-5911406027236569746089675554748135312500,40022343907270397352965399097794968750,40143602171843758666292581116410000,4025054350766002463251474186500,403350135768194635636171000,40416633063574896677300,405549939627039600,40612182993865,407163170,4081]),409([-117920, 15585808],410[426552448394636714720553816389274308035895411389805883034985546818882031845376,411-55876556222880738651382959148329502876096075327084935039031884373558741172224,4123393295715290183821010552313572221545212247684503012173117764703828786020352,413-125729166452196578653551230178028570067747190427221869867485520072257044480,4143121342502030777257351089270834971957072933779704445667351054593298530304,415-52544031605544530265465344472543470442324636919759253720520768014516224,416532110915869155495738137756847596184665209453108323879594125221167104,417-399031158106622651277981701966309467713625045637309782055519780864,418-101914346170769215732007802723651742508893380955930030421292613632,4192296526155500449624398016447877283594461904009374321659789443072,420-31950871094301541469458501953701002806003991982768349794795520,421329792235011603804948028315065667439678526339671142107709440,422-2655636715955021784085217734679612378726442691190553837568,42316825164648840434987220620681420687654501026066872664064,424-81705027839007003131400500185224450729843244954288128,425273656504606483403474090105104132405333665144373248,426-320807702482945680116212224172370503903312084992,427-3166683390779345463318656135338172047199043584,42827871349428383710305216046431806697565585408,429-132774697798318602604125735604528772808704,430436096215568182871014215818309741314048,431-964687143341252402362763535357837312,432942144169187362941776488535425024,4332794850106281773765892648206336,434-17236916236678037389276086272,43550979778712911923486851072,436-105035658611718440992768,437161833913559276412928,438-188675698610077696,439163929317513984,440-102098677888,44142387952,442-10184,4431]),444([-13760, 621264],445[-1961864562041980324821547425314935668736,446784270445793223959453256359333693751296,447-120528107728500223255333768387027271680,44810335626145581464192664472924270362624,449-568426570575654606865505142156820480,45021261993723422650574629752537088000,451-544630471727787626557612832587776,4528870521306520473088172555763712,453-54993059067301585878494740480,454-1434261324709904840432549888,45550978938193065926383894528,456-845761855773797582372864,4578627493611216601088000,458-48299605284169187328,459-32782260293713920,4603415534989828096,461-34580115625984,462199359712512,463-730488128,4641658080,465-2064,4661]),467([-3940515, 3010787550],468[-6458213126940667330314375,46934699336325466068070000,470-72461450055340471500,47168342601718080000,472-15140380554450,473-25802960400,47423981220,475-8160,4761]),477([-38907, -2953962], [-20349931239, -424530315, -134838, 53658, 429, 1]),478([-608, 5776],479[-34162868224,480-8540717056,4816405537792,482-1123778560,48384283392,484-2033152,485-92416,4866992,487-152,4881]),489([-9504, 365904], [1294672896, -92835072, 1463616, 7920, -264, 1]),490([-10395, 444150],491[-38324677699334121599624973029296875,492-17868327793500376961572310472656250,4932569568362004197901139023084765625,494-95128267987528547588017818750000,495-822168183291347061312510937500,496134395594560592096297190625000,497-2881389756919344324888937500,498-2503855007083401977250000,499922779077075655997443750,500-11503912310262102937500,501-18237870962450291250,5021457822151548910000,503-10087015556047500,504-13677678063000,505490243338900,506-2461460400,5075198445,508-4410,5091]),510([-856035, -341748450],511[103687510635057329105625,512961598491955315190000,5131054634146768300500,514-6553122389064000,515-14554350284850,516-2046589200,51713185540,5188160,5191]),520([-3267, -280962], [1480352841, -56169531, -2829222, 10890, 429, 1])]521"""522from sage.rings.all import RealField523from sage.misc.all import prod524ell = sporadic_j[j]525E = EllipticCurve(j=j).short_weierstrass_model()526a4a6 = list(E.ainvs())[3:]527L = E.period_lattice()528pr = 100529if ell==163:530pr=1000531elif ell>30:532pr=300533w1, w2 = L.basis(prec=pr)534X = polygen(RealField(pr),'X')535w = w1 # real period536if j in [-121, -24729001, -162677523113838677, QQ(-882216989)/131072]:537w = 2*w2-w1 # imaginary period538kerpol = prod(([X-L.elliptic_exponential(n*w/ell)[0] for n in range(1,(ell+1)//2)]))539kerpolcoeffs = [c.real().round() for c in list(kerpol)]540return (a4a6,kerpolcoeffs)541542def isogenies_sporadic_Q(E, l=None):543"""544Returns list of ``l`` -isogenies with domain ``E`` (defined over `\QQ`).545546Returns a list of sporadic l-isogenies from E (l = 11, 17, 19, 37,54743, 67 or 163). Only for elliptic curves over `\QQ`.548549INPUT:550551- ``E`` -- an elliptic curve defined over `\QQ`.552553- ``l`` -- either None or a prime number.554555OUTPUT:556557(list) If ``l`` is None, a list of all isogenies with domain ``E``558and of degree 11, 17, 19, 37, 43, 67 or 163; otherwise a list of559isogenies of the given degree.560561.. note::562563This function would normally be invoked indirectly via564``E.isogenies_prime_degree(l)``, which automatically calls the appropriate565function.566567EXAMPLES::568569sage: from sage.schemes.elliptic_curves.isogeny_small_degree import isogenies_sporadic_Q570sage: E = EllipticCurve('121a1')571sage: isogenies_sporadic_Q(E, 11)572[Isogeny of degree 11 from Elliptic Curve defined by y^2 + x*y + y = x^3 + x^2 - 30*x - 76 over Rational Field to Elliptic Curve defined by y^2 + x*y + y = x^3 + x^2 - 305*x + 7888 over Rational Field]573sage: isogenies_sporadic_Q(E, 13)574[]575sage: isogenies_sporadic_Q(E, 17)576[]577sage: isogenies_sporadic_Q(E)578[Isogeny of degree 11 from Elliptic Curve defined by y^2 + x*y + y = x^3 + x^2 - 30*x - 76 over Rational Field to Elliptic Curve defined by y^2 + x*y + y = x^3 + x^2 - 305*x + 7888 over Rational Field]579580sage: E = EllipticCurve([1, 1, 0, -660, -7600])581sage: isogenies_sporadic_Q(E, 17)582[Isogeny of degree 17 from Elliptic Curve defined by y^2 + x*y = x^3 + x^2 - 660*x - 7600 over Rational Field to Elliptic Curve defined by y^2 + x*y = x^3 + x^2 - 878710*x + 316677750 over Rational Field]583sage: isogenies_sporadic_Q(E)584[Isogeny of degree 17 from Elliptic Curve defined by y^2 + x*y = x^3 + x^2 - 660*x - 7600 over Rational Field to Elliptic Curve defined by y^2 + x*y = x^3 + x^2 - 878710*x + 316677750 over Rational Field]585sage: isogenies_sporadic_Q(E, 11)586[]587588sage: E = EllipticCurve([0, 0, 1, -1862, -30956])589sage: isogenies_sporadic_Q(E, 11)590[]591sage: isogenies_sporadic_Q(E, 19)592[Isogeny of degree 19 from Elliptic Curve defined by y^2 + y = x^3 - 1862*x - 30956 over Rational Field to Elliptic Curve defined by y^2 + y = x^3 - 672182*x + 212325489 over Rational Field]593sage: isogenies_sporadic_Q(E)594[Isogeny of degree 19 from Elliptic Curve defined by y^2 + y = x^3 - 1862*x - 30956 over Rational Field to Elliptic Curve defined by y^2 + y = x^3 - 672182*x + 212325489 over Rational Field]595596sage: E = EllipticCurve([0, -1, 0, -6288, 211072])597sage: E.conductor()59819600599sage: isogenies_sporadic_Q(E,37)600[Isogeny of degree 37 from Elliptic Curve defined by y^2 = x^3 - x^2 - 6288*x + 211072 over Rational Field to Elliptic Curve defined by y^2 = x^3 - x^2 - 163137088*x - 801950801728 over Rational Field]601602sage: E = EllipticCurve([1, 1, 0, -25178045, 48616918750])603sage: E.conductor()604148225605sage: isogenies_sporadic_Q(E,37)606[Isogeny of degree 37 from Elliptic Curve defined by y^2 + x*y = x^3 + x^2 - 25178045*x + 48616918750 over Rational Field to Elliptic Curve defined by y^2 + x*y = x^3 + x^2 - 970*x - 13075 over Rational Field]607608sage: E = EllipticCurve([-3440, 77658])609sage: E.conductor()610118336611sage: isogenies_sporadic_Q(E,43)612[Isogeny of degree 43 from Elliptic Curve defined by y^2 = x^3 - 3440*x + 77658 over Rational Field to Elliptic Curve defined by y^2 = x^3 - 6360560*x - 6174354606 over Rational Field]613614sage: E = EllipticCurve([-29480, -1948226])615sage: E.conductor()616287296617sage: isogenies_sporadic_Q(E,67)618[Isogeny of degree 67 from Elliptic Curve defined by y^2 = x^3 - 29480*x - 1948226 over Rational Field to Elliptic Curve defined by y^2 = x^3 - 132335720*x + 585954296438 over Rational Field]619620sage: E = EllipticCurve([-34790720, -78984748304])621sage: E.conductor()622425104623sage: isogenies_sporadic_Q(E,163)624[Isogeny of degree 163 from Elliptic Curve defined by y^2 = x^3 - 34790720*x - 78984748304 over Rational Field to Elliptic Curve defined by y^2 = x^3 - 924354639680*x + 342062961763303088 over Rational Field]625"""626if E.base_ring() != QQ:627raise ValueError("The elliptic curve must be defined over QQ.")628j = E.j_invariant()629j = QQ(j)630if (j not in sporadic_j631or (l is not None and sporadic_j[j] != l)):632return []633634data = _sporadic_Q_data(j)635Ew = E.short_weierstrass_model()636E_to_Ew = E.isomorphism_to(Ew)637c4, c6 = Ew.c_invariants()638(a4,a6), f = data639d = (c6*a4)/(18*c4*a6) # twisting factor640R = PolynomialRing(E.base_field(),'X')641n = len(f)642ker = R([d**(n-i-1) * f[i] for i in range(n)])643isog = Ew.isogeny(kernel=ker, degree=l, model="minimal", check=False)644isog.set_pre_isomorphism(E_to_Ew)645return [isog]646647648def isogenies_2(E):649"""650Returns a list of all 2-isogenies with domain ``E``.651652INPUT:653654- ``E`` -- an elliptic curve.655656OUTPUT:657658(list) 2-isogenies with domain ``E``. In general these are659normalised, but over `\QQ` the codomain is a minimal model.660661EXAMPLES::662663sage: from sage.schemes.elliptic_curves.isogeny_small_degree import isogenies_2664sage: E = EllipticCurve('14a1'); E665Elliptic Curve defined by y^2 + x*y + y = x^3 + 4*x - 6 over Rational Field666sage: [phi.codomain().ainvs() for phi in isogenies_2(E)]667[(1, 0, 1, -36, -70)]668669sage: E = EllipticCurve([1,2,3,4,5]); E670Elliptic Curve defined by y^2 + x*y + 3*y = x^3 + 2*x^2 + 4*x + 5 over Rational Field671sage: [phi.codomain().ainvs() for phi in isogenies_2(E)]672[]673sage: E = EllipticCurve(QQbar, [9,8]); E674Elliptic Curve defined by y^2 = x^3 + 9*x + 8 over Algebraic Field675sage: isogenies_2(E) # not implemented676"""677f2 = E.division_polynomial(2)678x2 = f2.roots(multiplicities=False)679x2.sort()680x = f2.parent().gen()681ff = [x-x2i for x2i in x2]682model = "minimal" if E.base_field() is QQ else None683isogs = [E.isogeny(f, model=model) for f in ff]684return isogs685686def isogenies_3(E):687"""688Returns a list of all 3-isogenies with domain ``E``.689690INPUT:691692- ``E`` -- an elliptic curve.693694OUTPUT:695696(list) 3-isogenies with domain ``E``. In general these are697normalised, but over `\QQ` the codomain is a minimal model.698699EXAMPLES::700701sage: from sage.schemes.elliptic_curves.isogeny_small_degree import isogenies_3702sage: E = EllipticCurve(GF(17), [1,1])703sage: [phi.codomain().ainvs() for phi in isogenies_3(E)]704[(0, 0, 0, 9, 7), (0, 0, 0, 0, 1)]705706sage: E = EllipticCurve(GF(17^2,'a'), [1,1])707sage: [phi.codomain().ainvs() for phi in isogenies_3(E)]708[(0, 0, 0, 9, 7), (0, 0, 0, 0, 1), (0, 0, 0, 5*a + 1, a + 13), (0, 0, 0, 12*a + 6, 16*a + 14)]709710sage: E = EllipticCurve('19a1')711sage: [phi.codomain().ainvs() for phi in isogenies_3(E)]712[(0, 1, 1, 1, 0), (0, 1, 1, -769, -8470)]713714sage: E = EllipticCurve([1,1])715sage: [phi.codomain().ainvs() for phi in isogenies_3(E)]716[]717"""718f3 = E.division_polynomial(3)719x3 = f3.roots(multiplicities=False)720x3.sort()721x = f3.parent().gen()722ff = [x-x3i for x3i in x3]723model = "minimal" if E.base_field() is QQ else None724isogs = [E.isogeny(f, model=model) for f in ff]725return isogs726727# 6 special cases: `l` = 5, 7, 13 and `j` = 0, 1728.728729def isogenies_5_0(E):730"""731Returns a list of all the 5-isogenies with domain ``E`` when the732j-invariant is 0.733734OUTPUT:735736(list) 5-isogenies with codomain E. In general these are737normalised, but over `\QQ` the codomain is a minimal model.738739.. note::740741This implementation requires that the characteristic is not 2,7423 or 5.743744.. note::745746This function would normally be invoked indirectly via ``E.isogenies_prime_degree(5)``.747748EXAMPLES::749750sage: from sage.schemes.elliptic_curves.isogeny_small_degree import isogenies_5_0751sage: E = EllipticCurve([0,12])752sage: isogenies_5_0(E)753[]754755sage: E = EllipticCurve(GF(13^2,'a'),[0,-3])756sage: isogenies_5_0(E)757[Isogeny of degree 5 from Elliptic Curve defined by y^2 = x^3 + 10 over Finite Field in a of size 13^2 to Elliptic Curve defined by y^2 = x^3 + (4*a+6)*x + (2*a+10) over Finite Field in a of size 13^2, Isogeny of degree 5 from Elliptic Curve defined by y^2 = x^3 + 10 over Finite Field in a of size 13^2 to Elliptic Curve defined by y^2 = x^3 + (12*a+5)*x + (2*a+10) over Finite Field in a of size 13^2, Isogeny of degree 5 from Elliptic Curve defined by y^2 = x^3 + 10 over Finite Field in a of size 13^2 to Elliptic Curve defined by y^2 = x^3 + (10*a+2)*x + (2*a+10) over Finite Field in a of size 13^2, Isogeny of degree 5 from Elliptic Curve defined by y^2 = x^3 + 10 over Finite Field in a of size 13^2 to Elliptic Curve defined by y^2 = x^3 + (3*a+12)*x + (11*a+12) over Finite Field in a of size 13^2, Isogeny of degree 5 from Elliptic Curve defined by y^2 = x^3 + 10 over Finite Field in a of size 13^2 to Elliptic Curve defined by y^2 = x^3 + (a+4)*x + (11*a+12) over Finite Field in a of size 13^2, Isogeny of degree 5 from Elliptic Curve defined by y^2 = x^3 + 10 over Finite Field in a of size 13^2 to Elliptic Curve defined by y^2 = x^3 + (9*a+10)*x + (11*a+12) over Finite Field in a of size 13^2]758759sage: K.<a> = NumberField(x**6-320*x**3-320)760sage: E = EllipticCurve(K,[0,0,1,0,0])761sage: isogenies_5_0(E)762[Isogeny of degree 5 from Elliptic Curve defined by y^2 + y = x^3 over Number Field in a with defining polynomial x^6 - 320*x^3 - 320 to Elliptic Curve defined by y^2 = x^3 + (a^5-400*a^2)*x + (280*a^3-3120) over Number Field in a with defining polynomial x^6 - 320*x^3 - 320,763Isogeny of degree 5 from Elliptic Curve defined by y^2 + y = x^3 over Number Field in a with defining polynomial x^6 - 320*x^3 - 320 to Elliptic Curve defined by y^2 = x^3 + (23/2*a^5-3700*a^2)*x + (-280*a^3+86480) over Number Field in a with defining polynomial x^6 - 320*x^3 - 320]764765"""766F = E.base_field()767if E.j_invariant() != 0:768raise ValueError("j-invariant must be 0.")769if F.characteristic() in [2,3,5]:770raise NotImplementedError("Not implemented in characteristic 2, 3 or 5.")771if not F(5).is_square():772return []773Ew = E.short_weierstrass_model()774a = Ew.a6()775x = polygen(F)776betas = (x**6-160*a*x**3-80*a**2).roots(multiplicities=False)777betas.sort()778if len(betas)==0:779return []780gammas = [(beta**2 *(beta**3-140*a))/(120*a) for beta in betas]781model = "minimal" if F is QQ else None782isogs = [Ew.isogeny(x**2+beta*x+gamma, model=model) for beta,gamma in zip(betas,gammas)]783iso = E.isomorphism_to(Ew)784[isog.set_pre_isomorphism(iso) for isog in isogs]785return isogs786787def isogenies_5_1728(E):788"""789Returns a list of 5-isogenies with domain ``E`` when the j-invariant is7901728.791792OUTPUT:793794(list) 5-isogenies with codomain E. In general these are795normalised; but if `-1` is a square then there are two796endomorphisms of degree `5`, for which the codomain is the same as797the domain curve; and over `\QQ`, the codomain is a minimal model.798799.. note::800801This implementation requires that the characteristic is not 2,8023 or 5.803804.. note::805806This function would normally be invoked indirectly via ``E.isogenies_prime_degree(5)``.807808EXAMPLES::809810sage: from sage.schemes.elliptic_curves.isogeny_small_degree import isogenies_5_1728811sage: E = EllipticCurve([7,0])812sage: isogenies_5_1728(E)813[]814815sage: E = EllipticCurve(GF(13),[11,0])816sage: isogenies_5_1728(E)817[Isogeny of degree 5 from Elliptic Curve defined by y^2 = x^3 + 11*x over Finite Field of size 13 to Elliptic Curve defined by y^2 = x^3 + 11*x over Finite Field of size 13,818Isogeny of degree 5 from Elliptic Curve defined by y^2 = x^3 + 11*x over Finite Field of size 13 to Elliptic Curve defined by y^2 = x^3 + 11*x over Finite Field of size 13]819820An example of endomorphisms of degree 5::821822sage: K.<i> = QuadraticField(-1)823sage: E = EllipticCurve(K,[0,0,0,1,0])824sage: isogenies_5_1728(E)825[Isogeny of degree 5 from Elliptic Curve defined by y^2 = x^3 + x over Number Field in i with defining polynomial x^2 + 1 to Elliptic Curve defined by y^2 = x^3 + x over Number Field in i with defining polynomial x^2 + 1,826Isogeny of degree 5 from Elliptic Curve defined by y^2 = x^3 + x over Number Field in i with defining polynomial x^2 + 1 to Elliptic Curve defined by y^2 = x^3 + x over Number Field in i with defining polynomial x^2 + 1]827sage: _[0].rational_maps()828(((4/25*i + 3/25)*x^5 + (4/5*i - 2/5)*x^3 - x)/(x^4 + (-4/5*i + 2/5)*x^2 + (-4/25*i - 3/25)),829((11/125*i + 2/125)*x^6*y + (-23/125*i + 64/125)*x^4*y + (141/125*i + 162/125)*x^2*y + (3/25*i - 4/25)*y)/(x^6 + (-6/5*i + 3/5)*x^4 + (-12/25*i - 9/25)*x^2 + (2/125*i - 11/125)))830831An example of 5-isogenies over a number field::832833sage: K.<a> = NumberField(x**4+20*x**2-80)834sage: K(5).is_square() #necessary but not sufficient!835True836sage: E = EllipticCurve(K,[0,0,0,1,0])837sage: isogenies_5_1728(E)838[Isogeny of degree 5 from Elliptic Curve defined by y^2 = x^3 + x over Number Field in a with defining polynomial x^4 + 20*x^2 - 80 to Elliptic Curve defined by y^2 = x^3 + (-20*a^2-39)*x + (35*a^3+112*a) over Number Field in a with defining polynomial x^4 + 20*x^2 - 80,839Isogeny of degree 5 from Elliptic Curve defined by y^2 = x^3 + x over Number Field in a with defining polynomial x^4 + 20*x^2 - 80 to Elliptic Curve defined by y^2 = x^3 + (-20*a^2-39)*x + (-35*a^3-112*a) over Number Field in a with defining polynomial x^4 + 20*x^2 - 80]840"""841F = E.base_field()842if E.j_invariant() != 1728:843raise ValueError("j-invariant must be 1728.")844if F.characteristic() in [2,3,5]:845raise NotImplementedError("Not implemented in characteristic 2, 3 or 5.")846model = "minimal" if F is QQ else None847# quick test for a negative answer (from Fricke module)848square5 = F(5).is_square()849square1 = F(-1).is_square()850if not square5 and not square1:851return []852Ew = E.short_weierstrass_model()853iso = E.isomorphism_to(Ew)854a = Ew.a4()855x = polygen(F)856isogs = []857# 2 cases858# Type 1: if -1 is a square we have 2 endomorphisms859if square1:860i = F(-1).sqrt()861isogs = [Ew.isogeny(f) for f in [x**2+a/(1+2*i), x**2+a/(1-2*i)]]862[isog.set_post_isomorphism(isog.codomain().isomorphism_to(E)) for isog in isogs]863# Type 2: if 5 is a square we have up to 4 (non-endomorphism) isogenies864if square5:865betas = (x**4+20*a*x**2-80*a**2).roots(multiplicities=False)866betas.sort()867gammas = [a*(beta**2-2)/6 for beta in betas]868isogs += [Ew.isogeny(x**2+beta*x+gamma, model=model) for beta,gamma in zip(betas,gammas)]869[isog.set_pre_isomorphism(iso) for isog in isogs]870return isogs871872def isogenies_7_0(E):873"""874Returns list of all 7-isogenies from E when the j-invariant is 0.875876OUTPUT:877878(list) 7-isogenies with codomain E. In general these are879normalised; but if `-3` is a square then there are two880endomorphisms of degree `7`, for which the codomain is the same as881the domain; and over `\QQ`, the codomain is a minimal model.882883.. note::884885This implementation requires that the characteristic is not 2,8863 or 7.887888.. note::889890This function would normally be invoked indirectly via ``E.isogenies_prime_degree(7)``.891892EXAMPLES:893894First some examples of endomorphisms::895896sage: from sage.schemes.elliptic_curves.isogeny_small_degree import isogenies_7_0897sage: K.<r> = QuadraticField(-3)898sage: E = EllipticCurve(K, [0,1])899sage: isogenies_7_0(E)900[Isogeny of degree 7 from Elliptic Curve defined by y^2 = x^3 + 1 over Number Field in r with defining polynomial x^2 + 3 to Elliptic Curve defined by y^2 = x^3 + 1 over Number Field in r with defining polynomial x^2 + 3,901Isogeny of degree 7 from Elliptic Curve defined by y^2 = x^3 + 1 over Number Field in r with defining polynomial x^2 + 3 to Elliptic Curve defined by y^2 = x^3 + 1 over Number Field in r with defining polynomial x^2 + 3]902903sage: E = EllipticCurve(GF(13^2,'a'),[0,-3])904sage: isogenies_7_0(E)905[Isogeny of degree 7 from Elliptic Curve defined by y^2 = x^3 + 10 over Finite Field in a of size 13^2 to Elliptic Curve defined by y^2 = x^3 + 10 over Finite Field in a of size 13^2, Isogeny of degree 7 from Elliptic Curve defined by y^2 = x^3 + 10 over Finite Field in a of size 13^2 to Elliptic Curve defined by y^2 = x^3 + 10 over Finite Field in a of size 13^2]906907Now some examples of 7-isogenies which are not endomorphisms::908909sage: K = GF(101)910sage: E = EllipticCurve(K, [0,1])911sage: isogenies_7_0(E)912[Isogeny of degree 7 from Elliptic Curve defined by y^2 = x^3 + 1 over Finite Field of size 101 to Elliptic Curve defined by y^2 = x^3 + 55*x + 100 over Finite Field of size 101, Isogeny of degree 7 from Elliptic Curve defined by y^2 = x^3 + 1 over Finite Field of size 101 to Elliptic Curve defined by y^2 = x^3 + 83*x + 26 over Finite Field of size 101]913914Examples over a number field::915916sage: from sage.schemes.elliptic_curves.isogeny_small_degree import isogenies_7_0917sage: E = EllipticCurve('27a1').change_ring(QuadraticField(-3,'r'))918sage: isogenies_7_0(E)919[Isogeny of degree 7 from Elliptic Curve defined by y^2 + y = x^3 + (-7) over Number Field in r with defining polynomial x^2 + 3 to Elliptic Curve defined by y^2 + y = x^3 + (-7) over Number Field in r with defining polynomial x^2 + 3,920Isogeny of degree 7 from Elliptic Curve defined by y^2 + y = x^3 + (-7) over Number Field in r with defining polynomial x^2 + 3 to Elliptic Curve defined by y^2 + y = x^3 + (-7) over Number Field in r with defining polynomial x^2 + 3]921922sage: K.<a> = NumberField(x^6 + 1512*x^3 - 21168)923sage: E = EllipticCurve(K, [0,1])924sage: isogs = isogenies_7_0(E)925sage: [phi.codomain().a_invariants() for phi in isogs]926[(0, 0, 0, -5/294*a^5 - 300/7*a^2, -55/2*a^3 - 1133),927(0, 0, 0, -295/1176*a^5 - 5385/14*a^2, 55/2*a^3 + 40447)]928sage: [phi.codomain().j_invariant() for phi in isogs]929[158428486656000/7*a^3 - 313976217600000,930-158428486656000/7*a^3 - 34534529335296000]931"""932if E.j_invariant()!=0:933raise ValueError("j-invariant must be 0.")934F = E.base_field()935if F.characteristic() in [2,3,7]:936raise NotImplementedError("Not implemented when the characteristic of the base field is 2, 3 or 7.")937x = polygen(F)938Ew = E.short_weierstrass_model()939iso = E.isomorphism_to(Ew)940a = Ew.a6()941model = "minimal" if F is QQ else None942943# there will be 2 endomorphisms if -3 is a square:944945ts = (x**2+3).roots(multiplicities=False)946ts.sort()947kers = [7*x-(2+6*t) for t in ts]948kers = [k(x**3/a).monic() for k in kers]949isogs = [Ew.isogeny(k,model=model) for k in kers]950if len(isogs)>0:951[endo.set_post_isomorphism(endo.codomain().isomorphism_to(E)) for endo in isogs]952953# we may have up to 6 other isogenies:954ts = (x**2-21).roots(multiplicities=False)955for t0 in ts:956s3 = a/(28+6*t0)957ss = (x**3-s3).roots(multiplicities=False)958ss.sort()959ker = x**3 - 2*t0*x**2 - 4*t0*x + 4*t0 + 28960kers = [ker(x/s).monic() for s in ss]961isogs += [Ew.isogeny(k, model=model) for k in kers]962963[isog.set_pre_isomorphism(iso) for isog in isogs]964return isogs965966def isogenies_7_1728(E):967"""968Returns list of all 7-isogenies from E when the j-invariant is 1728.969970OUTPUT:971972(list) 7-isogenies with codomain E. In general these are973normalised; but over `\QQ` the codomain is a minimal model.974975.. note::976977This implementation requires that the characteristic is not 2,9783, or 7.979980.. note::981982This function would normally be invoked indirectly via ``E.isogenies_prime_degree(7)``.983984EXAMPLES::985986sage: from sage.schemes.elliptic_curves.isogeny_small_degree import isogenies_7_1728987sage: E = EllipticCurve(GF(47), [1, 0])988sage: isogenies_7_1728(E)989[Isogeny of degree 7 from Elliptic Curve defined by y^2 = x^3 + x over Finite Field of size 47 to Elliptic Curve defined by y^2 = x^3 + 26 over Finite Field of size 47,990Isogeny of degree 7 from Elliptic Curve defined by y^2 = x^3 + x over Finite Field of size 47 to Elliptic Curve defined by y^2 = x^3 + 21 over Finite Field of size 47]991992An example in characteristic 53 (for which an earlier implementation did not work)::993994sage: from sage.schemes.elliptic_curves.isogeny_small_degree import isogenies_7_1728995sage: E = EllipticCurve(GF(53), [1, 0])996sage: isogenies_7_1728(E)997[]998sage: E = EllipticCurve(GF(53^2,'a'), [1, 0])999sage: [iso.codomain().ainvs() for iso in isogenies_7_1728(E)]1000[(0, 0, 0, 36, 19*a + 15), (0, 0, 0, 36, 34*a + 38), (0, 0, 0, 33, 39*a + 28), (0, 0, 0, 33, 14*a + 25), (0, 0, 0, 19, 45*a + 16), (0, 0, 0, 19, 8*a + 37), (0, 0, 0, 3, 45*a + 16), (0, 0, 0, 3, 8*a + 37)]10011002::10031004sage: K.<a> = NumberField(x^8 + 84*x^6 - 1890*x^4 + 644*x^2 - 567)1005sage: E = EllipticCurve(K, [1, 0])1006sage: isogs = isogenies_7_1728(E)1007sage: [phi.codomain().a_invariants() for phi in isogs]1008[(0,10090,10100,101135/636*a^6 + 55/12*a^4 - 79135/636*a^2 + 1127/212,1012155/636*a^7 + 245/12*a^5 - 313355/636*a^3 - 3577/636*a),1013(0,10140,10150,101635/636*a^6 + 55/12*a^4 - 79135/636*a^2 + 1127/212,1017-155/636*a^7 - 245/12*a^5 + 313355/636*a^3 + 3577/636*a)]1018sage: [phi.codomain().j_invariant() for phi in isogs]1019[-526110256146528/53*a^6 + 183649373229024*a^4 - 3333881559996576/53*a^2 + 2910267397643616/53,1020-526110256146528/53*a^6 + 183649373229024*a^4 - 3333881559996576/53*a^2 + 2910267397643616/53]1021sage: E1 = isogs[0].codomain()1022sage: E2 = isogs[1].codomain()1023sage: E1.is_isomorphic(E2)1024False1025sage: E1.is_quadratic_twist(E2)1026-11027"""1028if E.j_invariant()!=1728:1029raise ValueError("j_invariant must be 1728 (in base field).")1030F = E.base_field()1031if F.characteristic() in [2,3,7]:1032raise NotImplementedError("Not implemented when the characteristic of the base field is 2, 3 or 7.")1033Ew = E.short_weierstrass_model()1034iso = E.isomorphism_to(Ew)1035a = Ew.a4()10361037ts = (Fricke_module(7)-1728).numerator().roots(F,multiplicities=False)1038if len(ts)==0:1039return []1040ts.sort()1041isogs = []1042model = "minimal" if F is QQ else None1043x = polygen(F)1044for t0 in ts:1045s2 = a/t01046ss = (x**2-s2).roots(multiplicities=False)1047ss.sort()1048ker = 9*x**3 + (-3*t0**3 - 36*t0**2 - 123*t0)*x**2 + (-8*t0**3 - 101*t0**2 - 346*t0 + 35)*x - 7*t0**3 - 88*t0**2 - 296*t0 + 2810491050kers = [ker(x/s) for s in ss]1051isogs += [Ew.isogeny(k.monic(), model=model) for k in kers]1052[isog.set_pre_isomorphism(iso) for isog in isogs]1053return isogs10541055def isogenies_13_0(E):1056"""1057Returns list of all 13-isogenies from E when the j-invariant is 0.10581059OUTPUT:10601061(list) 13-isogenies with codomain E. In general these are1062normalised; but if `-3` is a square then there are two1063endomorphisms of degree `13`, for which the codomain is the same1064as the domain.10651066.. note::10671068This implementation requires that the characteristic is not 2,10693 or 13.10701071.. note::10721073This function would normally be invoked indirectly via ``E.isogenies_prime_degree(13)``.10741075EXAMPLES::10761077sage: from sage.schemes.elliptic_curves.isogeny_small_degree import isogenies_13_010781079Endomorphisms of degree 13 will exist when -3 is a square::10801081sage: K.<r> = QuadraticField(-3)1082sage: E = EllipticCurve(K, [0, r]); E1083Elliptic Curve defined by y^2 = x^3 + r over Number Field in r with defining polynomial x^2 + 31084sage: isogenies_13_0(E)1085[Isogeny of degree 13 from Elliptic Curve defined by y^2 = x^3 + r over Number Field in r with defining polynomial x^2 + 3 to Elliptic Curve defined by y^2 = x^3 + r over Number Field in r with defining polynomial x^2 + 3,1086Isogeny of degree 13 from Elliptic Curve defined by y^2 = x^3 + r over Number Field in r with defining polynomial x^2 + 3 to Elliptic Curve defined by y^2 = x^3 + r over Number Field in r with defining polynomial x^2 + 3]1087sage: isogenies_13_0(E)[0].rational_maps()1088(((7/338*r + 23/338)*x^13 + (-164/13*r - 420/13)*x^10 + (720/13*r + 3168/13)*x^7 + (3840/13*r - 576/13)*x^4 + (4608/13*r + 2304/13)*x)/(x^12 + (4*r + 36)*x^9 + (1080/13*r + 3816/13)*x^6 + (2112/13*r - 5184/13)*x^3 + (-17280/169*r - 1152/169)), ((18/2197*r + 35/2197)*x^18*y + (23142/2197*r + 35478/2197)*x^15*y + (-1127520/2197*r - 1559664/2197)*x^12*y + (-87744/2197*r + 5992704/2197)*x^9*y + (-6625152/2197*r - 9085824/2197)*x^6*y + (-28919808/2197*r - 2239488/2197)*x^3*y + (-1990656/2197*r - 3870720/2197)*y)/(x^18 + (6*r + 54)*x^15 + (3024/13*r + 11808/13)*x^12 + (31296/13*r + 51840/13)*x^9 + (487296/169*r - 2070144/169)*x^6 + (-940032/169*r + 248832/169)*x^3 + (1990656/2197*r + 3870720/2197)))10891090An example of endomorphisms over a finite field::10911092sage: K = GF(19^2,'a')1093sage: E = EllipticCurve(j=K(0)); E1094Elliptic Curve defined by y^2 = x^3 + 1 over Finite Field in a of size 19^21095sage: isogenies_13_0(E)1096[Isogeny of degree 13 from Elliptic Curve defined by y^2 = x^3 + 1 over Finite Field in a of size 19^2 to Elliptic Curve defined by y^2 = x^3 + 1 over Finite Field in a of size 19^2,1097Isogeny of degree 13 from Elliptic Curve defined by y^2 = x^3 + 1 over Finite Field in a of size 19^2 to Elliptic Curve defined by y^2 = x^3 + 1 over Finite Field in a of size 19^2]1098sage: isogenies_13_0(E)[0].rational_maps()1099((6*x^13 - 6*x^10 - 3*x^7 + 6*x^4 + x)/(x^12 - 5*x^9 - 9*x^6 - 7*x^3 + 5), (-8*x^18*y - 9*x^15*y + 9*x^12*y - 5*x^9*y + 5*x^6*y - 7*x^3*y + 7*y)/(x^18 + 2*x^15 + 3*x^12 - x^9 + 8*x^6 - 9*x^3 + 7))11001101A previous implementation did not work in some characteristics::11021103sage: K = GF(29)1104sage: E = EllipticCurve(j=K(0))1105sage: isogenies_13_0(E)1106[Isogeny of degree 13 from Elliptic Curve defined by y^2 = x^3 + 1 over Finite Field of size 29 to Elliptic Curve defined by y^2 = x^3 + 26*x + 12 over Finite Field of size 29, Isogeny of degree 13 from Elliptic Curve defined by y^2 = x^3 + 1 over Finite Field of size 29 to Elliptic Curve defined by y^2 = x^3 + 16*x + 28 over Finite Field of size 29]11071108::11091110sage: K = GF(101)1111sage: E = EllipticCurve(j=K(0)); E.ainvs()1112(0, 0, 0, 0, 1)1113sage: [phi.codomain().ainvs() for phi in isogenies_13_0(E)]1114[(0, 0, 0, 64, 36), (0, 0, 0, 42, 66)]11151116::11171118sage: x = polygen(QQ)1119sage: f = x^12 + 78624*x^9 - 130308048*x^6 + 2270840832*x^3 - 545001799681120sage: K.<a> = NumberField(f)1121sage: E = EllipticCurve(j=K(0)); E.ainvs()1122(0, 0, 0, 0, 1)1123sage: [phi.codomain().ainvs() for phi in isogenies_13_0(E)]1124[(0, 0, 0, -739946459/23857162861049856*a^11 - 2591641747/1062017577504*a^8 + 16583647773233/4248070310016*a^5 - 14310911337/378211388*a^2, 26146225/4248070310016*a^9 + 7327668845/14750244132*a^6 + 174618431365/756422776*a^3 - 378332499709/94552847), (0, 0, 0, 3501275/5964290715262464*a^11 + 24721025/531008788752*a^8 - 47974903745/1062017577504*a^5 - 6773483100/94552847*a^2, 6699581/4248070310016*a^9 + 1826193509/14750244132*a^6 - 182763866047/756422776*a^3 - 321460597/94552847)]1125"""1126if E.j_invariant()!=0:1127raise ValueError("j-invariant must be 0.")1128F = E.base_field()1129if F.characteristic() in [2,3,13]:1130raise NotImplementedError("Not implemented when the characteristic of the base field is 2, 3 or 13.")1131Ew = E.short_weierstrass_model()1132iso = E.isomorphism_to(Ew)1133a = Ew.a6()1134model = "minimal" if F is QQ else None1135x = polygen(F)11361137# there will be 2 endomorphisms if -3 is a square:1138ts = (x**2+3).roots(multiplicities=False)1139ts.sort()1140kers = [13*x**2 + (78*t + 26)*x + 24*t + 40 for t in ts]1141kers = [k(x**3/a).monic() for k in kers]1142isogs = [Ew.isogeny(k,model=model) for k in kers]1143if len(isogs)>0:1144[endo.set_post_isomorphism(endo.codomain().isomorphism_to(E)) for endo in isogs]11451146# we may have up to 12 other isogenies:1147ts = (x**4 + 7*x**3 + 20*x**2 + 19*x + 1).roots(multiplicities=False)1148ts.sort()1149for t0 in ts:1150s3 = a / (6*t0**3 + 32*t0**2 + 68*t0 + 4)1151ss = (x**3-s3).roots(multiplicities=False)1152ss.sort()1153ker = (x**6 + (20*t0**3 + 106*t0**2 + 218*t0 + 4)*x**51154+ (-826*t0**3 - 4424*t0**2 - 9244*t0 - 494)*x**41155+ (13514*t0**3 + 72416*t0**2 + 151416*t0 + 8238)*x**31156+ (-101948*t0**3 - 546304*t0**2 - 1142288*t0 - 62116)*x**21157+ (354472*t0**3 + 1899488*t0**2 + 3971680*t0 + 215960)*x1158- 459424*t0**3 - 2461888*t0**2 - 5147648*t0 - 279904)1159kers = [ker(x/s).monic() for s in ss]1160isogs += [Ew.isogeny(k, model=model) for k in kers]11611162[isog.set_pre_isomorphism(iso) for isog in isogs]11631164return isogs11651166def isogenies_13_1728(E):1167"""1168Returns list of all 13-isogenies from E when the j-invariant is 1728.11691170OUTPUT:11711172(list) 13-isogenies with codomain E. In general these are1173normalised; but if `-1` is a square then there are two1174endomorphisms of degree `13`, for which the codomain is the same1175as the domain; and over `\QQ`, the codomain is a minimal model.11761177.. note::11781179This implementation requires that the characteristic is not11802, 3 or 13.11811182.. note::11831184This function would normally be invoked indirectly via ``E.isogenies_prime_degree(13)``.11851186EXAMPLES::11871188sage: from sage.schemes.elliptic_curves.isogeny_small_degree import isogenies_13_172811891190sage: K.<i> = QuadraticField(-1)1191sage: E = EllipticCurve([0,0,0,i,0]); E.ainvs()1192(0, 0, 0, i, 0)1193sage: isogenies_13_1728(E)1194[Isogeny of degree 13 from Elliptic Curve defined by y^2 = x^3 + i*x over Number Field in i with defining polynomial x^2 + 1 to Elliptic Curve defined by y^2 = x^3 + i*x over Number Field in i with defining polynomial x^2 + 1,1195Isogeny of degree 13 from Elliptic Curve defined by y^2 = x^3 + i*x over Number Field in i with defining polynomial x^2 + 1 to Elliptic Curve defined by y^2 = x^3 + i*x over Number Field in i with defining polynomial x^2 + 1]11961197::11981199sage: K = GF(83)1200sage: E = EllipticCurve(K, [0,0,0,5,0]); E.ainvs()1201(0, 0, 0, 5, 0)1202sage: isogenies_13_1728(E)1203[]1204sage: K = GF(89)1205sage: E = EllipticCurve(K, [0,0,0,5,0]); E.ainvs()1206(0, 0, 0, 5, 0)1207sage: isogenies_13_1728(E)1208[Isogeny of degree 13 from Elliptic Curve defined by y^2 = x^3 + 5*x over Finite Field of size 89 to Elliptic Curve defined by y^2 = x^3 + 5*x over Finite Field of size 89,1209Isogeny of degree 13 from Elliptic Curve defined by y^2 = x^3 + 5*x over Finite Field of size 89 to Elliptic Curve defined by y^2 = x^3 + 5*x over Finite Field of size 89]12101211::12121213sage: K = GF(23)1214sage: E = EllipticCurve(K, [1,0])1215sage: isogenies_13_1728(E)1216[Isogeny of degree 13 from Elliptic Curve defined by y^2 = x^3 + x over Finite Field of size 23 to Elliptic Curve defined by y^2 = x^3 + 16 over Finite Field of size 23, Isogeny of degree 13 from Elliptic Curve defined by y^2 = x^3 + x over Finite Field of size 23 to Elliptic Curve defined by y^2 = x^3 + 7 over Finite Field of size 23]12171218::12191220sage: x = polygen(QQ)1221sage: f = x^12 + 1092*x^10 - 432432*x^8 + 6641024*x^6 - 282896640*x^4 - 149879808*x^2 - 3493601281222sage: K.<a> = NumberField(f)1223sage: E = EllipticCurve(K, [1,0])1224sage: [phi.codomain().ainvs() for phi in isogenies_13_1728(E)]1225[(0,12260,12270,122811090413835/20943727039698624*a^10 + 32280103535965/55849938772529664*a^8 - 355655987835845/1551387188125824*a^6 + 19216954517530195/5235931759924656*a^4 - 1079766118721735/5936430566808*a^2 + 156413528482727/8080141604822,1229214217013065/82065216155553792*a^11 + 1217882637605/427423000810176*a^9 - 214645003230565/189965778137856*a^7 + 22973355421236025/1282269002430528*a^5 - 2059145797340695/2544184528632*a^3 - 23198483147321/989405094468*a),1230(0,12310,12320,123311090413835/20943727039698624*a^10 + 32280103535965/55849938772529664*a^8 - 355655987835845/1551387188125824*a^6 + 19216954517530195/5235931759924656*a^4 - 1079766118721735/5936430566808*a^2 + 156413528482727/8080141604822,1234-214217013065/82065216155553792*a^11 - 1217882637605/427423000810176*a^9 + 214645003230565/189965778137856*a^7 - 22973355421236025/1282269002430528*a^5 + 2059145797340695/2544184528632*a^3 + 23198483147321/989405094468*a)]12351236"""1237if E.j_invariant()!=1728:1238raise ValueError("j-invariant must be 1728.")1239F = E.base_field()1240if F.characteristic() in [2, 3, 13]:1241raise NotImplementedError("Not implemented when the characteristic of the base field is 2, 3 or 13.")1242Ew = E.short_weierstrass_model()1243iso = E.isomorphism_to(Ew)1244a = Ew.a4()1245model = "minimal" if F is QQ else None1246x = polygen(F)12471248# we will have two endomorphisms if -1 is a square:1249ts = (x**2+1).roots(multiplicities=False)1250ts.sort()1251kers = [13*x**3 + (-26*i - 13)*x**2 + (-52*i - 13)*x - 2*i - 3 for i in ts]1252kers = [k(x**2/a).monic() for k in kers]1253isogs = [Ew.isogeny(k,model=model) for k in kers]1254if len(isogs)>0:1255[endo.set_post_isomorphism(endo.codomain().isomorphism_to(E)) for endo in isogs]12561257# we may have up to 12 other isogenies:12581259ts = (x**6 + 10*x**5 + 46*x**4 + 108*x**3 + 122*x**2 + 38*x - 1).roots(multiplicities=False)1260ts.sort()1261for t0 in ts:1262s2 = a/(66*t0**5 + 630*t0**4 + 2750*t0**3 + 5882*t0**2 + 5414*t0 + 162)1263ss = (x**2-s2).roots(multiplicities=False)1264ss.sort()1265ker = (x**6 + (-66*t0**5 - 630*t0**4 - 2750*t0**3 - 5882*t0**21266- 5414*t0 - 162)*x**5 + (-21722*t0**5 - 205718*t0**4 -1267890146*t0**3 - 1873338*t0**2 - 1652478*t0 + 61610)*x**41268+ (-3391376*t0**5 - 32162416*t0**4 - 139397232*t0**3 -1269294310576*t0**2 - 261885968*t0 + 6105552)*x**3 +1270(-241695080*t0**5 - 2291695976*t0**4 - 9930313256*t0**31271- 20956609720*t0**2 - 18625380856*t0 + 469971320)*x**2 +1272(-8085170432*t0**5 - 76663232384*t0**4 -1273332202985024*t0**3 - 701103233152*t0**2 -1274623190845440*t0 + 15598973056)*x - 101980510208*t0**5 -1275966973468160*t0**4 - 4190156868352*t0**3 -12768843158270336*t0**2 - 7860368751232*t0 + 196854655936)12771278kers = [ker(x/s).monic() for s in ss]1279isogs += [Ew.isogeny(k, model=model) for k in kers]12801281[isog.set_pre_isomorphism(iso) for isog in isogs]12821283return isogs12841285# List of primes l for which X_0(l) is (hyper)elliptic and X_0^+(l) has genus 012861287hyperelliptic_primes = [11, 17, 19, 23, 29, 31, 41, 47, 59, 71]12881289@cached_function1290def _hyperelliptic_isogeny_data(l):1291r"""1292Helper function for elliptic curve isogenies.12931294INPUT:12951296- ``l`` -- a prime in [11, 17, 19, 23, 29, 31, 41, 47, 59, 71]12971298OUTPUT:12991300- A dict holding a collection of precomputed data needed for computing `l`-isogenies.13011302EXAMPLES::13031304sage: from sage.schemes.elliptic_curves.isogeny_small_degree import _hyperelliptic_isogeny_data1305sage: HID = _hyperelliptic_isogeny_data(11)1306sage: HID['A2']130755*u - 331308sage: HID['A4']1309-183*u^2 + 738*u - 180*v - 1351310sage: HID['A6']13111330*u^3 - 11466*u^2 + 1332*u*v + 2646*u - 1836*v + 18901312sage: HID['alpha']1313u^11 - 55*u^10 + 1188*u^9 - 12716*u^8 + 69630*u^7 - 177408*u^6 + 133056*u^5 + 132066*u^4 - 187407*u^3 + 40095*u^2 + 24300*u - 67501314sage: HID['beta']1315u^9 - 47*u^8 + 843*u^7 - 7187*u^6 + 29313*u^5 - 48573*u^4 + 10665*u^3 + 27135*u^2 - 12150*u1316sage: HID['hyper_poly']1317u^4 - 16*u^3 + 2*u^2 + 12*u - 713181319sage: _hyperelliptic_isogeny_data(37)1320Traceback (most recent call last):1321...1322ValueError: 37 must be one of [11, 17, 19, 23, 29, 31, 41, 47, 59, 71].13231324"""1325if not l in hyperelliptic_primes:1326raise ValueError("%s must be one of %s."%(l,hyperelliptic_primes))1327data = {}1328Zu = PolynomialRing(ZZ,'u')1329Zuv = PolynomialRing(ZZ,['u','v'])1330Zxuv = PolynomialRing(ZZ,['x','u','v'])1331x,u,v = Zxuv.gens()1332if l == 11:1333data['hyper_poly'] = Zu([-7, 12, 2, -16, 1])1334data['A2'] = Zu([-33, 55])1335data['A4'] = Zuv(Zu([-135, 738, -183])+v*Zu([-180]))1336data['A6'] = Zuv(Zu([1890, 2646, -11466, 1330]) + v*Zu([-1836, 1332]))1337data['alpha'] = Zu([-6750, 24300, 40095, -187407, 132066, 133056, -177408, 69630, -12716, 1188, -55, 1])1338data['beta'] = Zu([0, -12150, 27135, 10665, -48573, 29313, -7187, 843, -47, 1])1339#beta factors as (u - 15) * (u - 6) * (u - 3) * (u - 1) * u * (u**2 - 12*u - 9) * (u**2 - 10*u + 5)1340return data1341if l == 17:1342data['hyper_poly'] = Zu([-8, 4, -3, -10, 1])1343data['A2'] = Zu([68, -204, 136])1344data['A4'] = Zuv(Zu([60, 720, -2595, 2250, -435]) + v*Zu([-360, 792, -432]))1345data['A6'] = Zuv(Zu([-8512, 22608, -5064, -57528, 87288, -43704, 4912] ) + v*Zu( [2520, -15372, 28098, -20160, 4914]))1346data['alpha'] = Zu([16000, -67200, 2720, 557600, -1392232, 1073992, 1104830, -3131026, 2450210, 73746, -1454945, 1110355, -424065, 95659, -13243, 1105, -51, 1])1347data['beta'] = Zu([0, 22400, -105920, 146208, 111616, -593800, 680948, -102282, -457950, 468035, -219274, 58549, -9374, 889, -46, 1])1348#beta factors as (u - 10) * (u - 5) * (u - 2) * (u - 1) * u * (u + 1) * (u**2 - 10*u + 7) * (u**2 - 6*u - 4) * (u**2 - 4*u + 2) * (u**3 - 9*u**2 + 8*u - 4)1349data['endo'] = 17*x**8 + 17*(-4*u + 4)*v*x**6 + 17*(4*u + 6)*v**2*x**4 + 17*(4*u + 4)*v**3*x**2 + (-4*u + 1)*v**41350data['endo_u'] = 11351return data1352if l == 19:1353data['hyper_poly'] = Zu([-8, 20, -8, -8, 1])1354data['A2'] = Zu([-114, 57, 171])1355data['A4'] = Zuv(Zu([-1020, 444, 2733, 726, -543]) + v*Zu([-180, -720, -540]))1356data['A6'] = Zuv(Zu([-10080, 21816, 54324, -37386, -86742, -20070, 6858]) + v*Zu([-2968, -13748, -11284, 6356, 6860]))1357data['alpha'] = Zu([16000, -22400, -337440, 475456, 1562104, -1988616, -3025294, 3245960, 2833014, -2420087, -1140950, 932406, 129580, -180443, 21090, 11153, -4066, 570, -38, 1])1358data['beta'] = Zu([0, 33600, -8160, -292400, 23472, 791244, 39282, -847909, -47024, 392654, -24046, -82469, 19162, 4833, -2652, 446, -34, 1])1359#beta factors as (u - 7) * (u - 2) * (u - 1) * u * (u + 1) * (u**2 - 8*u - 4) * (u**2 - 6*u - 15) * (u**2 - 5*u - 5) * (u**2 - 5*u + 2) * (u**2 - 2*u - 4) * (u**2 + u - 1)1360data['endo'] = 19*x**9 + 19*(-12*u - 24)*v*x**6 + 19*(-24*u - 24)*v**2*x**3 + (96*u - 224)*v**31361data['endo_u'] = -11362return data1363if l == 23:1364data['hyper_poly'] = Zu([-7, 10, -11, 2, 2, -8, 1])1365data['A2'] = Zu([69, -230, 253])1366data['A4'] = Zuv(Zu([405, 180, -930, 2820, -795]) + v*Zu([360, -792]))1367data['A6'] = Zuv(Zu([-15498, 34020, -36918, -8120, 51114, -72492, 12166]) + v*Zu([-1080, 7704, -24840, 12168]))1368data['alpha'] = Zu([-6750, 48600, -83835, -170775, 1115109, -2492280, 2732814, -116403, -4877702, 8362616, -6612454, 302266, 5423124, -6447728, 3209696, 336674, -1470068, 953856, -336927, 74221, -10465, 920, -46, 1])1369data['beta'] = Zu( [0, 12150, -72495, 168588, -144045, -254034, 930982, -1256170, 604358, 693650, -1563176, 1271974, -225188, -444070, 421050, -184350, 47754, -7696, 759, -42, 1])1370#beta factors as (u - 5) * (u - 3) * (u - 2) * (u - 1) * u * (u + 1) * (u**2 - 8*u + 3) * (u**2 - 6*u - 9) * (u**3 - 7*u**2 + 3*u - 5) * (u**3 - 7*u**2 + 7*u - 3) * (u**4 - 4*u**3 - 1)1371return data1372if l == 29:1373data['hyper_poly'] = Zu([-7, 8, 8, 2, -12, -4, 1])1374data['A2'] = Zu([-174, -232, 348, 406])1375data['A4'] = Zuv(Zu([-1215, -3096, 132, 7614, 6504, -360, -1263] ) + v*Zu( [180, -720, -2160, -1260]))1376data['A6'] = Zuv(Zu([-18900, -63504, 24696, 285068, 285264, -185136, -506268, -275520, 504, 24388] ) + v*Zu( [4482, -2448, -59868, -94968, -18144, 48276, 24390]))1377data['alpha'] = Zu([-6750, -12150, 281880, 570024, -1754181, -5229135, 2357613, 19103721, 9708910, -31795426, -38397537, 19207947, 54103270, 9216142, -37142939, -18871083, 14041394, 10954634, -3592085, -3427365, 853818, 622398, -189399, -53679, 26680, -580, -1421, 319, -29, 1])1378data['beta'] = Zu([0, -24300, -57510, 257850, 839187, -373185, -3602119, -2371192, 5865017, 8434433, -2363779, -10263744, -2746015, 5976011, 3151075, -2093854, -1356433, 569525, 299477, -129484, -28279, 19043, -895, -1076, 273, -27, 1])1379#beta factors as (u - 3) * (u - 1) * u * (u + 1) * (u + 2) * (u**2 - 6*u + 2) * (u**2 - 5*u - 5) * (u**2 - 5*u + 3) * (u**2 - 3*u - 9) * (u**2 - u - 3) * (u**2 - u - 1) * (u**2 + u - 1) * (u**3 - 4*u**2 - 6*u - 5) * (u**4 - 2*u**3 - 5*u**2 - 4*u - 1)1380data['endo'] = 29*x**14 + 29*(-14*u + 3)*v*x**12 + 29*(-20*u + 73)*v**2*x**10 + 29*(-58*u + 115)*v**3*x**8 + 29*(-56*u + 59)*v**4*x**6 + 29*(30*u + 1)*v**5*x**4 + 29*(12*u - 5)*v**6*x**2 + (2*u + 5)*v**71381data['endo_u'] = -11382return data1383if l == 31:1384data['hyper_poly'] = Zu([-3, -14, -11, 18, 6, -8, 1])1385data['A2'] = Zu([558, 837, -1488, 465])1386data['A4'] = Zuv(Zu([-4140, -12468, 15189, 16956, -27054, 11184, -1443]) + v*Zu([2160, -7560, 6120, -1440]))1387data['A6'] = Zuv(Zu([71280, 592056, -108324, -2609730, 2373048, 1282266, -2793204, 1530882, -356976, 29790]) + v*Zu([-81312, 181664, 294728, -868392, 701400, -238840, 29792]))1388data['alpha'] = Zu([108000, 475200, -7053120, -27353408, 90884374, 303670296, -665806437, -1361301729, 3259359840, 2249261823, -9368721606, 2279583264, 13054272515, -12759480061, -4169029296, 14390047139, -7803693550, -2988803682, 6239473912, -3296588360, 134066754, 908915598, -685615437, 294482733, -87483178, 18983315, -3052818, 361336, -30659, 1767, -62, 1])1389data['beta'] = Zu([0, 712800, 1216080, -18430560, -15262464, 168899202, -12931221, -720077416, 624871714, 1239052988, -2259335558, 68648452, 2679085427, -2318039014, -229246628, 1710545918, -1243026758, 211524870, 296674626, -291810274, 145889932, -48916468, 11793961, -2085662, 269348, -24778, 1540, -58, 1])1390#beta factors as (u - 3) * (u - 2) * (u - 1) * u * (u + 1) * (u**2 - 8*u + 11) * (u**2 - 7*u + 2) * (u**2 - 5*u - 2) * (u**2 - 5*u + 5) * (u**2 - 4*u - 4) * (u**2 - 4*u - 1) * (u**2 - 2*u - 1) * (u**2 - u - 1) * (u**3 - 9*u**2 + 21*u - 15) * (u**4 - 8*u**3 + 8*u**2 + 12*u - 9)1391data['endo'] = 31*x**15 + 31*(-66*u + 86)*v*x**12 + 31*(168*u + 280)*v**2*x**9 + 31*(576*u + 1792)*v**3*x**6 + 31*(384*u + 896)*v**4*x**3 + (-3072*u - 2048)*v**51392data['endo_u'] = 21393return data1394if l == 41:1395data['hyper_poly'] = Zu([-8, -20, -15, 8, 20, 10, -8, -4, 1])1396data['A2'] = Zu([328, 656, -656, -1148, 820])1397data['A4'] = Zuv(Zu([-1380, -4008, 1701, 10872, 6144, -18378, -2160, 9732, -2523]) + v*Zu([720, -1440, -2160, 5400, -2520]))1398data['A6'] = Zuv(Zu([4480, 155616, 16080, -550720, -343968, 832680, 938632, -621648, -1468608, 953920, 427632, -413016, 68920]) + v*Zu([-14616, 6804, 96390, -2016, -324324, 184464, 260568, -276192, 68922]))1399data['alpha'] = Zu([16000, 67200, -465760, -2966432, -1742664, 20985112, 46140990, -31732934, -217030548, -147139488, 436080674, 745775322, -271341362, -1542677562, -605560447, 1832223375, 1772593672, -1270633050, -2400692229, 343522723, 2179745361, 282422801, -1503727029, -421357697, 879637411, 261059095, -462271351, -61715127, 193718727, -24135265, -49355103, 20512341, 3613289, -4706595, 1099661, 163057, -162483, 46617, -7544, 738, -41, 1])1400data['beta'] = Zu([0, 44800, 167040, -447040, -2734272, -1104272, 13488360, 21067652, -24681704, -83929974, -8986886, 169059382, 127641266, -196479899, -283039783, 124573790, 366614063, -12946368, -332987597, -58867672, 241909907, 60568430, -155045647, -17919564, 79114945, -12025938, -24060781, 11190142, 1979597, -2931764, 750233, 110144, -122263, 37484, -6439, 666, -39, 1])1401#beta factors as (u - 5) * (u - 2) * (u - 1) * u * (u + 1) * (u**2 - 5*u + 5) * (u**2 - 3*u - 7) * (u**2 - 2*u - 4) * (u**2 - 2*u - 1) * (u**2 - u - 1) * (u**2 - 2) * (u**2 + u - 1) * (u**3 - 3*u**2 - 5*u - 2) * (u**3 - 2*u**2 - 2*u - 1) * (u**4 - 6*u**3 + 5*u**2 + 2*u - 1) * (u**4 - 5*u**3 + u**2 + 4) * (u**4 - 4*u**3 + 2)1402data['endo'] = 41*x**20 + 41*(-12*u - 22)*v*x**18 + 41*(-252*u - 247)*v**2*x**16 + 41*(-176*u - 424)*v**3*x**14 + 41*(464*u - 254)*v**4*x**12 + 41*(1688*u - 868)*v**5*x**10 + 41*(1720*u - 1190)*v**6*x**8 + 41*(528*u - 232)*v**7*x**6 + 41*(16*u + 29)*v**8*x**4 + 41*(20*u + 10)*v**9*x**2 + (4*u + 5)*v**101403data['endo_u'] = 11404return data1405if l == 47:1406data['hyper_poly'] = Zu([-11, 28, -38, 30, -13, -16, 19, -24, 11, -6, 1])1407data['A2'] = Zu([376, -1504, 2209, -1598, 1081])1408data['A4'] = Zuv(Zu([2400, -4080, -1440, 18000, -26355, 34740, -22050, 12900, -3315]) + v*Zu([1152, -3384, 3672, -3312]))1409data['A6'] = Zuv(Zu([-119504, 606336, -1505280, 2109392, -1509360, -515808, 2920702, -4614012, 4334322, -3260312, 1571442, -622428, 103822]) + v*Zu([2016, 48384, -235872, 438984, -627480, 503496, -311976, 103824]))1410data['alpha'] = Zu([-65536, 688128, -2502656, -96256, 38598656, -187217920, 508021120, -845669120, 552981696, 1469334304, -5945275904, 11705275552, -14673798654, 9100068184, 8421580132, -34288012648, 56657584158, -60426283952, 36612252089, 9942017442, -60791892299, 93046207239, -92028642340, 59196883097, -10454018992, -33364599371, 57280402355, -57873890484, 41879296232, -20241250112, 2065827049, 8435506655, -11611941072, 10182603298, -7040645261, 4071881378, -2013138357, 856757031, -313468474, 97893151, -25770006, 5617769, -990431, 136864, -14194, 1034, -47, 1])1411data['beta'] = Zu([0, 114688, -1114112, 4854784, -11205632, 7426048, 42663936, -182555136, 394092544, -508851472, 213245648, 743315936, -2203729384, 3409478688, -3280008936, 1139839970, 2576264698, -6272528962, 8005203155, -6671665088, 2744569094, 1996771588, -5520074039, 6637395180, -5455622885, 3028415830, -601645255, -1012737914, 1632999370, -1525982346, 1093778952, -644352392, 319489974, -134176208, 47566499, -14083902, 3424200, -667810, 101271, -11438, 901, -44, 1])1412#beta factors as (u - 4) * (u - 2) * (u - 1) * u * (u + 1) * (u**2 - 5*u + 2) * (u**2 - 2*u - 1) * (u**3 - 5*u**2 + 5*u - 7) * (u**3 - 4*u**2 + 3*u - 4) * (u**3 - 4*u**2 + 3*u - 1) * (u**3 - 3*u**2 + 2*u - 4) * (u**3 - 2*u**2 + 2*u - 2) * (u**3 + u + 1) * (u**4 - 4*u**3 - 2*u**2 - 4) * (u**5 - 5*u**4 + 5*u**3 - 11*u**2 + 6*u - 4) * (u**6 - 4*u**5 + 2*u**4 - 4*u**3 - u**2 + 4*u - 2)1413return data1414if l == 59:1415data['hyper_poly'] = Zu([-8, -4, 20, -24, -3, 40, -62, 40, 3, -28, 22, -8, 1])1416data['A2'] = Zu([590, -1475, -295, 4130, -4425, 1711])1417data['A4'] = Zuv(Zu([-2460, 8844, -3843, -20718, 57153, -50418, -12600, 72762, -69339, 30978, -5223]) + v*Zu([900, 360, -7560, 10800, -5220]))1418data['A6'] = Zuv(Zu([25760, -373560, 568020, 1147870, -4634370, 5318070, 1631996, -14270202, 21535998, -14119408, -2820102, 14275410, -13535292, 6790074, -1847898, 205378]) + v*Zu([-23688, 27972, 183708, -696024, 721980, 453600, -1925028, 2039184, -1027404, 205380]))1419data['alpha'] = Zu([16000, -67200, -783520, 5573376, -5127336, -60792184, 241324042, -170978932, -1262437160, 4310971231, -3953349811, -10887235780, 41679530185, -51342089572, -33068562195, 230682514316, -372641172307, 121615007703, 682044179678, -1549365239197, 1373184591667, 614906882627, -3566756201696, 4920423266916, -2342393877496, -3589340274442, 8772457933356, -8488557160148, 1742977715620, 7131088674129, -11643540780203, 8512399456274, -315658868113, -6917286294515, 8713332734648, -5190227733987, -54249978263, 3397583328372, -3658171840037, 1987950394792, -179519591637, -748989116551, 800595050760, -459184355769, 134398080099, 28871590941, -64236756338, 46651654354, -23352309386, 9059054346, -2830320860, 721829600, -150487052, 25475079, -3452149, 365800, -29205, 1652, -59, 1])1420data['beta'] = Zu([0, -56000, 320800, 391440, -7693120, 21125500, 11515130, -204780145, 486681785, -102547033, -2147060784, 5552726794, -4419031758, -9431888681, 33728080307, -42367773552, -2994127157, 105330637610, -188172973931, 127559513693, 123083802224, -421097252069, 490425751691, -161944881372, -408669953969, 799965143719, -668167261718, 69589638764, 563644022562, -787681290965, 505670881115, 2900924856, -364669742737, 407962360532, -223582547975, 9985786664, 102435489491, -105519055992, 58212400117, -14331637533, -6742538722, 10205452686, -6853903214, 3244679736, -1188153136, 347102566, -81626216, 15409226, -2307408, 268126, -23322, 1429, -55, 1])1421#beta factors as (u - 2) * (u - 1) * u * (u + 1) * (u**2 - 4*u - 1) * (u**2 - 3*u - 5) * (u**2 - 3*u - 2) * (u**2 - 3*u + 1) * (u**2 - u - 1) * (u**3 - 6*u**2 + 10*u - 7) * (u**3 - 5*u**2 + 7*u - 5) * (u**3 - 3*u**2 + 2*u - 1) * (u**3 - u**2 + 1) * (u**4 - 5*u**3 + 4*u**2 - 1) * (u**4 - 4*u**3 + 3*u**2 + 2*u - 4) * (u**4 - 3*u**3 - u - 1) * (u**4 - u**3 + 2*u - 1) * (u**5 - 6*u**4 + 10*u**3 - 11*u**2 + 8*u - 4) * (u**6 - 5*u**5 + 5*u**4 - 5*u**2 + 5*u - 5)1422return data1423if l == 71:1424data['hyper_poly'] = Zu([-7, 6, -27, 40, -58, 66, -66, 40, 15, -48, 66, -66, 37, -10, 1])1425data['A2'] = Zu([213, -1420, 4260, -4970, 9940, -9088, 2485])1426data['A4'] = Zuv(Zu([2565, -10008, 18024, -26532, 23208, 7584, -104418, 189432, -251736, 275148, -182232, 60144, -7563]) + v*Zu([720, -4320, 7560, -20160, 23040, -7560]))1427data['A6'] = Zuv(Zu([-69930, 382536, -1898568, 5206124, -11813256, 23115792, -35705670, 44318064, -41531952, 20674360, 23881872, -77986944, 114989770, -124612152, 103122936, -59431204, 21485688, -4294416, 357910]) + v*Zu([18576, -53856, 57672, 161856, -961920, 3199176, -5706288, 8032896, -9352584, 6786720, -2505888, 357912]))1428data['alpha'] = Zu([-6750, 97200, -603855, 2263977, -4854483, -2486349, 75190491, -399596520, 1441975423, -4089818964, 9450153463, -17516526653, 23635982289, -11859874932, -53385529273, 230566737711, -585283867605, 1136695427037, -1753961304140, 2020891913264, -1147488305875, -1930304898882, 8102336330029, -17218530732347, 27006964902986, -32365758791872, 25902000374138, -468390635342, -46332664858222, 107139839089502, -162234735929274, 182582147217312, -140033523896938, 22513210292184, 152367877270246, -334009986053250, 451855980915164, -443144048889720, 284518400252142, -11142427766850, -289840331821002, 512373447321402, -576967281819172, 466024421705696, -230395084854230, -36287337331916, 241209603962570, -330646545417814, 304702155703516, -205131886553392, 87504290135653, 5131997859077, -54867900326127, 66216047255551, -54817285755105, 36239054778472, -20052219750661, 9464634765852, -3841191816845, 1343947848527, -405138280373, 104923131180, -23228729413, 4364552115, -689157169, 90223321, -9613968, 812240, -52327, 2414, -71, 1])1429data['beta'] = Zu([0, 12150, -163215, 1115640, -5311143, 18820224, -50700172, 99823812, -102454041, -183909134, 1354660714, -4462311942, 10695310224, -20015395554, 28262441676, -23240987282, -17879387475, 124501604946, -315187724212, 564766450688, -765154573538, 705985549104, -115433273216, -1206098873334, 3175185881748, -5228317292044, 6292310032120, -5077451367560, 719644756530, 6451571564682, -14460150103020, 19999710623352, -19681838601268, 11819712227412, 2180981559572, -17790742756618, 29025463386612, -31179247603548, 23207078145510, -8345354986332, -7468523752270, 18486966963350, -21719818051100, 17831212433536, -10100011266030, 2336962513536, 2906983627184, -4989755986066, 4711466210012, -3361479243242, 1952316811463, -948555371584, 389878900245, -136099552242, 40341734984, -10121407164, 2136756509, -376218102, 54551634, -6399080, 591884, -41538, 2078, -66, 1])1430#beta factors as (u - 3) * (u - 2) * (u - 1) * u * (u + 1) * (u**2 - 5*u + 5) * (u**2 - 3*u + 1) * (u**2 - 2*u - 1) * (u**2 - u - 1) * (u**3 - 5*u**2 + 5*u - 3) * (u**3 - 4*u**2 - 1) * (u**3 - 2*u**2 - 1) * (u**4 - 6*u**3 + 7*u**2 + 6*u - 9) * (u**4 - 5*u**3 + 4*u**2 + u + 3) * (u**4 - 5*u**3 + 6*u**2 - 3*u + 5) * (u**4 - 4*u**3 + u**2 - 4*u + 1) * (u**4 - 4*u**3 + 2*u**2 - u + 1) * (u**4 - 2*u**3 - 3*u**2 - 2*u - 1) * (u**4 - 2*u**3 + u - 1) * (u**6 - 5*u**5 + 8*u**4 - 7*u**3 + 6*u**2 - 3*u + 1) * (u**8 - 6*u**7 + 9*u**6 - 2*u**5 + 2*u**3 - 9*u**2 + 2*u - 1)1431return data14321433@cached_function1434def Psi2(l):1435"""1436Returns the generic kernel polynomial for hyperelliptic `l`-isogenies.14371438INPUT:14391440- ``l`` -- either 11, 17, 19, 23, 29, 31, 41, 47, 59, or 71.14411442OUTPUT:14431444The generic `l`-kernel polynomial.14451446TESTS::14471448sage: from sage.schemes.elliptic_curves.isogeny_small_degree import Psi21449sage: Psi2(11)1450x^5 - 55*x^4*u + 994*x^3*u^2 - 8774*x^2*u^3 + 41453*x*u^4 - 928945/11*u^5 + 33*x^4 + 276*x^3*u - 7794*x^2*u^2 + 4452*x*u^3 + 1319331/11*u^4 + 216*x^3*v - 4536*x^2*u*v + 31752*x*u^2*v - 842616/11*u^3*v + 162*x^3 + 38718*x^2*u - 610578*x*u^2 + 33434694/11*u^3 - 4536*x^2*v + 73872*x*u*v - 2745576/11*u^2*v - 16470*x^2 + 580068*x*u - 67821354/11*u^2 - 185976*x*v + 14143896/11*u*v + 7533*x - 20437029/11*u - 12389112/11*v + 19964151/1114511452"""1453if not l in hyperelliptic_primes:1454raise ValueError("%s must be one of %s."%(l,hyperelliptic_primes))14551456data = _hyperelliptic_isogeny_data(l)1457R = PolynomialRing(QQ,['x','u'])1458x, u = R.gens()1459L = PolynomialRing(R,'y')1460y = L.gen()1461K = R.extension(y**2-R(data['hyper_poly']),name = 'v')1462v = K.gen()1463from sage.categories.homset import Hom1464h = Hom(K,K)(-v)14651466A = K(data['A4'])1467B = K(data['A6'])1468Abar = h(A)*l**21469Bbar = -h(B)*l**31470s1 = K(data['A2'])14711472d = (l-1)//21473s = [1]1474t = [d,s1,((1-10*d)*A - Abar)*(1/QQ(30))]1475t += [((1-28*d)*B - 42*t[1]*A - Bbar)*(1/QQ(70))]1476c = [0,6*t[2] + 2*A*t[0],10*t[3] + 6*A*t[1] + 4*B*t[0]]1477for n in range(2,d):1478k = sum(c[i]*c[n-i] for i in range(1,n))1479c += [(3*k-(2*n-1)*(n-1)*A*c[n-1]-(2*n-2)*(n-2)*B*c[n-2])*(1/QQ((n-1)*(2*n+5)))]1480for n in range(3,d):1481t += [(c[n]-(4*n-2)*A*t[n-1]-(4*n-4)*B*t[n-2])*(1/QQ(4*n+2))]1482for n in range(1,d+1):1483s += [(-1/QQ(n))*sum((-1)**i*t[i]*s[n-i] for i in range(1,n+1))]1484psi = sum((-1)**i*s[i]*x**(d-i) for i in range(0,d+1))1485R = PolynomialRing(QQ,['x','u','v'])1486return R(psi)148714881489def isogenies_prime_degree_genus_plus_0(E, l=None):1490"""1491Returns list of ``l`` -isogenies with domain ``E``.14921493INPUT:14941495- ``E`` -- an elliptic curve.14961497- ``l`` -- either None or 11, 17, 19, 23, 29, 31, 41, 47, 59, or 71.14981499OUTPUT:15001501(list) When ``l`` is None a list of all isogenies of degree 11, 17, 19, 23,150229, 31, 41, 47, 59, or 71, otherwise a list of isogenies of the given degree.15031504.. note::15051506This function would normally be invoked indirectly via1507``E.isogenies_prime_degree(l)``, which automatically calls the appropriate function.15081509ALGORITHM:15101511See [KT2013]_, Chapter 5.15121513EXAMPLES::15141515sage: from sage.schemes.elliptic_curves.isogeny_small_degree import isogenies_prime_degree_genus_plus_015161517sage: E = EllipticCurve('121a1')1518sage: isogenies_prime_degree_genus_plus_0(E, 11)1519[Isogeny of degree 11 from Elliptic Curve defined by y^2 + x*y + y = x^3 + x^2 - 30*x - 76 over Rational Field to Elliptic Curve defined by y^2 + x*y + y = x^3 + x^2 - 305*x + 7888 over Rational Field]15201521sage: E = EllipticCurve([1, 1, 0, -660, -7600])1522sage: isogenies_prime_degree_genus_plus_0(E, 17)1523[Isogeny of degree 17 from Elliptic Curve defined by y^2 + x*y = x^3 + x^2 - 660*x - 7600 over Rational Field to Elliptic Curve defined by y^2 + x*y = x^3 + x^2 - 878710*x + 316677750 over Rational Field]15241525sage: E = EllipticCurve([0, 0, 1, -1862, -30956])1526sage: isogenies_prime_degree_genus_plus_0(E, 19)1527[Isogeny of degree 19 from Elliptic Curve defined by y^2 + y = x^3 - 1862*x - 30956 over Rational Field to Elliptic Curve defined by y^2 + y = x^3 - 672182*x + 212325489 over Rational Field]15281529sage: K = QuadraticField(-295,'a')1530sage: a = K.gen()1531sage: E = EllipticCurve_from_j(-484650135/16777216*a + 4549855725/16777216)1532sage: isogenies_prime_degree_genus_plus_0(E, 23)1533[Isogeny of degree 23 from Elliptic Curve defined by y^2 = x^3 + (-14460494784192904095/140737488355328*a+270742665778826768325/140737488355328)*x + (37035998788154488846811217135/590295810358705651712*a-1447451882571839266752561148725/590295810358705651712) over Number Field in a with defining polynomial x^2 + 295 to Elliptic Curve defined by y^2 = x^3 + (-5130542435555445498495/140737488355328*a+173233955029127361005925/140737488355328)*x + (-1104699335561165691575396879260545/590295810358705651712*a+3169785826904210171629535101419675/590295810358705651712) over Number Field in a with defining polynomial x^2 + 295]15341535sage: K = QuadraticField(-199,'a')1536sage: a = K.gen()1537sage: E = EllipticCurve_from_j(94743000*a + 269989875)1538sage: isogenies_prime_degree_genus_plus_0(E, 29)1539[Isogeny of degree 29 from Elliptic Curve defined by y^2 = x^3 + (-153477413215038000*a+5140130723072965125)*x + (297036215130547008455526000*a+2854277047164317800973582250) over Number Field in a with defining polynomial x^2 + 199 to Elliptic Curve defined by y^2 = x^3 + (251336161378040805000*a-3071093219933084341875)*x + (-8411064283162168580187643221000*a+34804337770798389546017184785250) over Number Field in a with defining polynomial x^2 + 199]15401541sage: K = QuadraticField(253,'a')1542sage: a = K.gen()1543sage: E = EllipticCurve_from_j(208438034112000*a - 3315409892960000)1544sage: isogenies_prime_degree_genus_plus_0(E, 31)1545[Isogeny of degree 31 from Elliptic Curve defined by y^2 = x^3 + (4146345122185433034677956608000*a-65951656549965037259634800640000)*x + (-18329111516954473474583425393698245080252416000*a+291542366110383928366510368064204147260129280000) over Number Field in a with defining polynomial x^2 - 253 to Elliptic Curve defined by y^2 = x^3 + (200339763852548615776123686912000*a-3186599019027216904280948275200000)*x + (7443671791411479629112717260182286294850207744000*a-118398847898864757209685951728838895495168655360000) over Number Field in a with defining polynomial x^2 - 253]15461547sage: E = EllipticCurve_from_j(GF(5)(1))1548sage: isogenies_prime_degree_genus_plus_0(E, 41)1549[Isogeny of degree 41 from Elliptic Curve defined by y^2 = x^3 + x + 2 over Finite Field of size 5 to Elliptic Curve defined by y^2 = x^3 + x + 3 over Finite Field of size 5, Isogeny of degree 41 from Elliptic Curve defined by y^2 = x^3 + x + 2 over Finite Field of size 5 to Elliptic Curve defined by y^2 = x^3 + x + 3 over Finite Field of size 5]15501551sage: K = QuadraticField(5,'a')1552sage: a = K.gen()1553sage: E = EllipticCurve_from_j(184068066743177379840*a - 411588709724712960000)1554sage: isogenies_prime_degree_genus_plus_0(E, 47) # long time (4.3s)1555[Isogeny of degree 47 from Elliptic Curve defined by y^2 = x^3 + (454562028554080355857852049849975895490560*a-1016431595837124114668689286176511361024000)*x + (-249456798429896080881440540950393713303830363999480904280965120*a+557802358738710443451273320227578156598454035482869042774016000) over Number Field in a with defining polynomial x^2 - 5 to Elliptic Curve defined by y^2 = x^3 + (39533118442361013730577638493616965245992960*a-88398740199669828340617478832005245173760000)*x + (214030321479466610282320528611562368963830105830555363061803253760*a-478586348074220699687616322532666163722004497458452316582576128000) over Number Field in a with defining polynomial x^2 - 5]15561557sage: K = QuadraticField(-66827,'a')1558sage: a = K.gen()1559sage: E = EllipticCurve_from_j(-98669236224000*a + 4401720074240000)1560sage: isogenies_prime_degree_genus_plus_0(E, 59) # long time (25s, 2012)1561[Isogeny of degree 59 from Elliptic Curve defined by y^2 = x^3 + (2605886146782144762297974784000*a+1893681048912773634944634716160000)*x + (-116918454256410782232296183198067568744071168000*a+17012043538294664027185882358514011304812871680000) over Number Field in a with defining polynomial x^2 + 66827 to Elliptic Curve defined by y^2 = x^3 + (-19387084027159786821400775098368000*a-4882059104868154225052787156713472000)*x + (-25659862010101415428713331477227179429538847260672000*a-2596038148441293485938798119003462972840818381946880000) over Number Field in a with defining polynomial x^2 + 66827]15621563sage: E = EllipticCurve_from_j(GF(13)(5))1564sage: isogenies_prime_degree_genus_plus_0(E, 71) # long time1565[Isogeny of degree 71 from Elliptic Curve defined by y^2 = x^3 + x + 4 over Finite Field of size 13 to Elliptic Curve defined by y^2 = x^3 + 10*x + 7 over Finite Field of size 13, Isogeny of degree 71 from Elliptic Curve defined by y^2 = x^3 + x + 4 over Finite Field of size 13 to Elliptic Curve defined by y^2 = x^3 + 10*x + 7 over Finite Field of size 13]15661567sage: E = EllipticCurve(GF(13),[0,1,1,1,0])1568sage: isogenies_prime_degree_genus_plus_0(E)1569[Isogeny of degree 17 from Elliptic Curve defined by y^2 + y = x^3 + x^2 + x over Finite Field of size 13 to Elliptic Curve defined by y^2 + y = x^3 + x^2 + 10*x + 1 over Finite Field of size 13,1570Isogeny of degree 17 from Elliptic Curve defined by y^2 + y = x^3 + x^2 + x over Finite Field of size 13 to Elliptic Curve defined by y^2 + y = x^3 + x^2 + 12*x + 4 over Finite Field of size 13,1571Isogeny of degree 29 from Elliptic Curve defined by y^2 + y = x^3 + x^2 + x over Finite Field of size 13 to Elliptic Curve defined by y^2 + y = x^3 + x^2 + 12*x + 6 over Finite Field of size 13,1572Isogeny of degree 29 from Elliptic Curve defined by y^2 + y = x^3 + x^2 + x over Finite Field of size 13 to Elliptic Curve defined by y^2 + y = x^3 + x^2 + 5*x + 6 over Finite Field of size 13,1573Isogeny of degree 41 from Elliptic Curve defined by y^2 + y = x^3 + x^2 + x over Finite Field of size 13 to Elliptic Curve defined by y^2 + y = x^3 + x^2 + 12*x + 4 over Finite Field of size 13,1574Isogeny of degree 41 from Elliptic Curve defined by y^2 + y = x^3 + x^2 + x over Finite Field of size 13 to Elliptic Curve defined by y^2 + y = x^3 + x^2 + 5*x + 6 over Finite Field of size 13]15751576"""1577if l is None:1578return sum([isogenies_prime_degree_genus_plus_0(E, l) for l in hyperelliptic_primes],[])15791580if not l in hyperelliptic_primes:1581raise ValueError("%s must be one of %s."%(l,hyperelliptic_primes))15821583F = E.base_ring()1584j = E.j_invariant()1585if F.characteristic() in [2, 3, l]:1586raise NotImplementedError("11, 17, 19, 23, 29, 31, 41, 47, 59, and 71-isogenies are not yet implemented in characteristic 2 and 3, and when the characteristic is the same as the degree of the isogeny.")15871588if j == F(0):1589return isogenies_prime_degree_genus_plus_0_j0(E, l)1590if j == F(1728):1591return isogenies_prime_degree_genus_plus_0_j1728(E, l)15921593Fu = PolynomialRing(F,'u')1594u = Fu.gen()1595Fuv = PolynomialRing(F,['u','v'])1596Fxuv = PolynomialRing(F,['x','u','v'])1597X = u1598data = _hyperelliptic_isogeny_data(l)1599a = Fu(data['alpha'])1600b = Fu(data['beta'])1601f = Fu(data['hyper_poly'])1602P = a1603Q = Fu((a**2 - f*b**2)/Fu(4))1604u_list = (j**2-P*j+Q).roots(multiplicities=False)16051606S = []1607for u0 in u_list:1608if b(u0) == 0:1609S += [[u0,v0] for v0 in (X**2-f(u0)).roots(multiplicities=False)]1610else:1611S += [[u0,(2*j-a(u0))/b(u0)]]1612if len(S) == 0:1613return []1614S.sort()16151616c4, c6 = E.c_invariants()1617b2 = E.b2()1618kernels = []16191620psi = Fxuv(Psi2(l))1621for u0, v0 in S:1622A4 = Fuv(data['A4'])(u0,v0) #non-zero since j!=01623A6 = Fuv(data['A6'])(u0,v0) #non-zero since j!=17281624T = (c4*A6)/(2*c6*A4)1625kernels += [psi((36*X+3*b2)*T,u0,v0).monic()]1626return [E.isogeny(ker) for ker in kernels]1627162816291630def isogenies_prime_degree_genus_plus_0_j0(E, l):1631"""1632Returns a list of hyperelliptic ``l`` -isogenies with domain ``E`` when `j(E)=0`.16331634INPUT:16351636- ``E`` -- an elliptic curve with j-invariant 0.16371638- ``l`` -- 11, 17, 19, 23, 29, 31, 41, 47, 59, or 71.16391640OUTPUT:16411642(list) a list of all isogenies of degree 11, 17, 19, 23, 29, 31, 41, 47, 59, or 71.16431644.. note::16451646This implementation requires that the characteristic is not 2, 3 or ``l``.16471648.. note::16491650This function would normally be invoked indirectly via ``E.isogenies_prime_degree(l)``.16511652EXAMPLES::16531654sage: from sage.schemes.elliptic_curves.isogeny_small_degree import isogenies_prime_degree_genus_plus_0_j016551656sage: u = polygen(QQ)1657sage: K.<a> = NumberField(u^4+228*u^3+486*u^2-540*u+225)1658sage: E = EllipticCurve(K,[0,-121/5*a^3-20691/5*a^2-29403/5*a+3267])1659sage: isogenies_prime_degree_genus_plus_0_j0(E,11)1660[Isogeny of degree 11 from Elliptic Curve defined by y^2 = x^3 + (-121/5*a^3-20691/5*a^2-29403/5*a+3267) over Number Field in a with defining polynomial x^4 + 228*x^3 + 486*x^2 - 540*x + 225 to Elliptic Curve defined by y^2 = x^3 + (-44286*a^2+178596*a-32670)*x + (-17863351/5*a^3+125072739/5*a^2-74353653/5*a-682803) over Number Field in a with defining polynomial x^4 + 228*x^3 + 486*x^2 - 540*x + 225, Isogeny of degree 11 from Elliptic Curve defined by y^2 = x^3 + (-121/5*a^3-20691/5*a^2-29403/5*a+3267) over Number Field in a with defining polynomial x^4 + 228*x^3 + 486*x^2 - 540*x + 225 to Elliptic Curve defined by y^2 = x^3 + (-3267*a^3-740157*a^2+600039*a-277695)*x + (-17863351/5*a^3-4171554981/5*a^2+3769467867/5*a-272366523) over Number Field in a with defining polynomial x^4 + 228*x^3 + 486*x^2 - 540*x + 225]16611662sage: E = EllipticCurve(GF(5^6,'a'),[0,1])1663sage: isogenies_prime_degree_genus_plus_0_j0(E,17)1664[Isogeny of degree 17 from Elliptic Curve defined by y^2 = x^3 + 1 over Finite Field in a of size 5^6 to Elliptic Curve defined by y^2 = x^3 + 2 over Finite Field in a of size 5^6, Isogeny of degree 17 from Elliptic Curve defined by y^2 = x^3 + 1 over Finite Field in a of size 5^6 to Elliptic Curve defined by y^2 = x^3 + 2 over Finite Field in a of size 5^6, Isogeny of degree 17 from Elliptic Curve defined by y^2 = x^3 + 1 over Finite Field in a of size 5^6 to Elliptic Curve defined by y^2 = x^3 + 2 over Finite Field in a of size 5^6, Isogeny of degree 17 from Elliptic Curve defined by y^2 = x^3 + 1 over Finite Field in a of size 5^6 to Elliptic Curve defined by y^2 = x^3 + 2 over Finite Field in a of size 5^6, Isogeny of degree 17 from Elliptic Curve defined by y^2 = x^3 + 1 over Finite Field in a of size 5^6 to Elliptic Curve defined by y^2 = x^3 + 2 over Finite Field in a of size 5^6, Isogeny of degree 17 from Elliptic Curve defined by y^2 = x^3 + 1 over Finite Field in a of size 5^6 to Elliptic Curve defined by y^2 = x^3 + 2 over Finite Field in a of size 5^6, Isogeny of degree 17 from Elliptic Curve defined by y^2 = x^3 + 1 over Finite Field in a of size 5^6 to Elliptic Curve defined by y^2 = x^3 + 2 over Finite Field in a of size 5^6, Isogeny of degree 17 from Elliptic Curve defined by y^2 = x^3 + 1 over Finite Field in a of size 5^6 to Elliptic Curve defined by y^2 = x^3 + 2 over Finite Field in a of size 5^6, Isogeny of degree 17 from Elliptic Curve defined by y^2 = x^3 + 1 over Finite Field in a of size 5^6 to Elliptic Curve defined by y^2 = x^3 + 2 over Finite Field in a of size 5^6, Isogeny of degree 17 from Elliptic Curve defined by y^2 = x^3 + 1 over Finite Field in a of size 5^6 to Elliptic Curve defined by y^2 = x^3 + 2 over Finite Field in a of size 5^6, Isogeny of degree 17 from Elliptic Curve defined by y^2 = x^3 + 1 over Finite Field in a of size 5^6 to Elliptic Curve defined by y^2 = x^3 + 2 over Finite Field in a of size 5^6, Isogeny of degree 17 from Elliptic Curve defined by y^2 = x^3 + 1 over Finite Field in a of size 5^6 to Elliptic Curve defined by y^2 = x^3 + 2 over Finite Field in a of size 5^6, Isogeny of degree 17 from Elliptic Curve defined by y^2 = x^3 + 1 over Finite Field in a of size 5^6 to Elliptic Curve defined by y^2 = x^3 + 2 over Finite Field in a of size 5^6, Isogeny of degree 17 from Elliptic Curve defined by y^2 = x^3 + 1 over Finite Field in a of size 5^6 to Elliptic Curve defined by y^2 = x^3 + 2 over Finite Field in a of size 5^6, Isogeny of degree 17 from Elliptic Curve defined by y^2 = x^3 + 1 over Finite Field in a of size 5^6 to Elliptic Curve defined by y^2 = x^3 + 2 over Finite Field in a of size 5^6, Isogeny of degree 17 from Elliptic Curve defined by y^2 = x^3 + 1 over Finite Field in a of size 5^6 to Elliptic Curve defined by y^2 = x^3 + 2 over Finite Field in a of size 5^6, Isogeny of degree 17 from Elliptic Curve defined by y^2 = x^3 + 1 over Finite Field in a of size 5^6 to Elliptic Curve defined by y^2 = x^3 + 2 over Finite Field in a of size 5^6, Isogeny of degree 17 from Elliptic Curve defined by y^2 = x^3 + 1 over Finite Field in a of size 5^6 to Elliptic Curve defined by y^2 = x^3 + 2 over Finite Field in a of size 5^6]1665"""1666if not l in hyperelliptic_primes:1667raise ValueError("%s must be one of %s."%(l,hyperelliptic_primes))1668F = E.base_field()1669if E.j_invariant() != 0:1670raise ValueError,("j-invariant must be 0.")1671if F.characteristic() in [2,3,l]:1672raise NotImplementedError("Not implemented in characteristic 2, 3 or l.")16731674Fu = PolynomialRing(F,'u')1675u = Fu.gen()1676Fuv = PolynomialRing(F,['u','v'])1677Fxuv = PolynomialRing(F,['x','u','v'])1678X = u1679data = _hyperelliptic_isogeny_data(l)1680a = Fu(data['alpha'])1681b = Fu(data['beta'])1682f = Fu(data['hyper_poly'])1683Q = Fu((a**2 - f*b**2)/Fu(4))1684u_list = Q.roots(multiplicities=False)1685c6, b2 = E.c6(), E.b2()1686kernels = []16871688if l % 3 == 1 and F(-3).is_square():1689p = F(-3).sqrt()1690endo = Fxuv(data['endo'])1691kernels += [endo(36*X+3*b2,p,-54*c6).monic(), endo(36*X+3*b2,-p,-54*c6).monic()]16921693S = []1694for u0 in u_list:1695if l % 3 == 1 and u0 == F(data['endo_u']):1696continue1697if b(u0) == 0:1698S += [[u0,v0] for v0 in (X**2-f(u0)).roots(multiplicities=False)]1699else:1700S += [[u0,-a(u0)/b(u0)]]1701if len(S)==0 and len(kernels) == 0:1702return []1703S.sort()17041705psi = Fxuv(Psi2(l))1706for u0,v0 in S:1707A6 = Fuv(data['A6'])(u0,v0) # non-zero since j!=17281708kernels += [psi((36*X+3*b2)*T,u0,v0).monic() for T in (X**3-A6/(-54*c6)).roots(multiplicities=False)]1709return [E.isogeny(ker) for ker in kernels]17101711def isogenies_prime_degree_genus_plus_0_j1728(E, l):1712"""1713Returns a list of ``l`` -isogenies with domain ``E`` when `j(E)=1728`.17141715INPUT:17161717- ``E`` -- an elliptic curve with j-invariant 1728.17181719- ``l`` -- 11, 17, 19, 23, 29, 31, 41, 47, 59, or 71.17201721OUTPUT:17221723(list) a list of all isogenies of degree 11, 17, 19, 23, 29, 31, 41, 47, 59, or 71.17241725.. note::17261727This implementation requires that the characteristic is not 2, 3 or ``l``.17281729.. note::17301731This function would normally be invoked indirectly via ``E.isogenies_prime_degree(l)``.17321733EXAMPLES::17341735sage: from sage.schemes.elliptic_curves.isogeny_small_degree import isogenies_prime_degree_genus_plus_0_j172817361737sage: u = polygen(QQ)1738sage: K.<a> = NumberField(u^6 - 522*u^5 - 10017*u^4 + 2484*u^3 - 5265*u^2 + 12150*u - 5103)1739sage: E = EllipticCurve(K,[-75295/1335852*a^5+13066735/445284*a^4+44903485/74214*a^3+17086861/24738*a^2+11373021/16492*a-1246245/2356,0])1740sage: isogenies_prime_degree_genus_plus_0_j1728(E,11)1741[Isogeny of degree 11 from Elliptic Curve defined by y^2 = x^3 + (-75295/1335852*a^5+13066735/445284*a^4+44903485/74214*a^3+17086861/24738*a^2+11373021/16492*a-1246245/2356)*x over Number Field in a with defining polynomial x^6 - 522*x^5 - 10017*x^4 + 2484*x^3 - 5265*x^2 + 12150*x - 5103 to Elliptic Curve defined by y^2 = x^3 + (9110695/1335852*a^5-1581074935/445284*a^4-5433321685/74214*a^3-3163057249/24738*a^2+1569269691/16492*a+73825125/2356)*x + (-3540460*a^3+30522492*a^2-7043652*a-5031180) over Number Field in a with defining polynomial x^6 - 522*x^5 - 10017*x^4 + 2484*x^3 - 5265*x^2 + 12150*x - 5103, Isogeny of degree 11 from Elliptic Curve defined by y^2 = x^3 + (-75295/1335852*a^5+13066735/445284*a^4+44903485/74214*a^3+17086861/24738*a^2+11373021/16492*a-1246245/2356)*x over Number Field in a with defining polynomial x^6 - 522*x^5 - 10017*x^4 + 2484*x^3 - 5265*x^2 + 12150*x - 5103 to Elliptic Curve defined by y^2 = x^3 + (9110695/1335852*a^5-1581074935/445284*a^4-5433321685/74214*a^3-3163057249/24738*a^2+1569269691/16492*a+73825125/2356)*x + (3540460*a^3-30522492*a^2+7043652*a+5031180) over Number Field in a with defining polynomial x^6 - 522*x^5 - 10017*x^4 + 2484*x^3 - 5265*x^2 + 12150*x - 5103]1742sage: i = QuadraticField(-1,'i').gen()1743sage: E = EllipticCurve([-1-2*i,0])1744sage: isogenies_prime_degree_genus_plus_0_j1728(E,17)1745[Isogeny of degree 17 from Elliptic Curve defined by y^2 = x^3 + (-2*i-1)*x over Number Field in i with defining polynomial x^2 + 1 to Elliptic Curve defined by y^2 = x^3 + (-82*i-641)*x over Number Field in i with defining polynomial x^2 + 1,1746Isogeny of degree 17 from Elliptic Curve defined by y^2 = x^3 + (-2*i-1)*x over Number Field in i with defining polynomial x^2 + 1 to Elliptic Curve defined by y^2 = x^3 + (-562*i+319)*x over Number Field in i with defining polynomial x^2 + 1]1747sage: Emin = E.global_minimal_model()1748sage: [(p,len(isogenies_prime_degree_genus_plus_0_j1728(Emin,p))) for p in [17, 29, 41]]1749[(17, 2), (29, 2), (41, 2)]1750"""1751if not l in hyperelliptic_primes:1752raise ValueError("%s must be one of %s."%(l,hyperelliptic_primes))1753F = E.base_ring()1754if E.j_invariant() != 1728:1755raise ValueError("j-invariant must be 1728.")1756if F.characteristic() in [2,3,l]:1757raise NotImplementedError("Not implemented in characteristic 2, 3 or l.")17581759Fu = PolynomialRing(F,'u')1760u = Fu.gen()1761Fuv = PolynomialRing(F,['u','v'])1762Fxuv = PolynomialRing(F,['x','u','v'])1763X = u1764data = _hyperelliptic_isogeny_data(l)1765a = Fu(data['alpha'])1766b = Fu(data['beta'])1767f = Fu(data['hyper_poly'])1768P = a1769Q = Fu((a**2 - f*b**2)/Fu(4))1770u_list = (1728**2-P*1728+Q).roots(multiplicities=False)1771c4, b2 = E.c4(), E.b2()1772kernels = []17731774if l % 4 == 1 and F(-1).is_square():1775i = F(-1).sqrt()1776endo = Fxuv(data['endo'])1777kernels += [endo(36*X+3*b2,i,-27*c4).monic(), endo(36*X+3*b2,-i,-27*c4).monic()]17781779S = []1780for u0 in u_list:1781if l % 4 == 1 and u0 == F(data['endo_u']):1782continue1783if b(u0) == 0:1784S += [[u0,v0] for v0 in (X**2-f(u0)).roots(multiplicities=False)]1785else:1786S += [[u0,(2*1728-a(u0))/b(u0)]]1787if len(S)==0 and len(kernels) == 0:1788return []1789S.sort()17901791psi = Fxuv(Psi2(l))1792for u0,v0 in S:1793A4 = Fuv(data['A4'])(u0,v0) # non-zero since j!=01794kernels += [psi((36*X+3*b2)*T,u0,v0).monic() for T in (X**2-A4/(-27*c4)).roots(multiplicities=False)]1795return [E.isogeny(ker) for ker in kernels]17961797@cached_function1798def _least_semi_primitive(p):1799"""1800Returns the smallest semi-primitive root modulo `p`, i.e., generator of the group `(\ZZ/p\ZZ)^*/\{1,-1\}`.18011802INPUT:18031804- ``p`` -- an odd prime.18051806OUTPUT:18071808the smallest semi-primitive root modulo `p`.18091810.. note::18111812This function would normally be invoked indirectly via ``E.isogenies_prime_degree_general(l)``.18131814EXAMPLES::18151816sage: from sage.schemes.elliptic_curves.isogeny_small_degree import _least_semi_primitive1817sage: _least_semi_primitive(5)181821819sage: _least_semi_primitive(13)182021821sage: _least_semi_primitive(17)182231823sage: _least_semi_primitive(73)182451825sage: _least_semi_primitive(997)182671827"""1828if not p.is_prime() or p<3:1829raise ValueError("%s is not an odd prime"%p)18301831def is_semi_primitive(a,p):1832from sage.rings.finite_rings.integer_mod_ring import Integers1833d = Integers(p)(a).multiplicative_order()1834if p%4==1:1835return (d==p-1)1836else:1837return d >= (p-1)/218381839a = 21840while not is_semi_primitive(a,p):1841a += 11842return a184318441845def isogenies_prime_degree_general(E, l):1846"""1847Returns a list of ``l`` -isogenies with domain ``E``.18481849INPUT:18501851- ``E`` -- an elliptic curve.18521853- ``l`` -- a prime.18541855OUTPUT:18561857(list) a list of all isogenies of degree l.18581859ALGORITHM:18601861This algorithm factors the ``l``-division polynomial, then1862combines its factors to otain kernels. See [KT2013]_, Chapter 3.18631864.. note::18651866This function works for any prime `l`. Normally one should use1867the function :meth:`isogenies_prime_degree` which uses special1868functions for certain small primes.18691870EXAMPLES::18711872sage: from sage.schemes.elliptic_curves.isogeny_small_degree import isogenies_prime_degree_general1873sage: E = EllipticCurve_from_j(GF(2^6,'a')(1))1874sage: isogenies_prime_degree_general(E, 7)1875[Isogeny of degree 7 from Elliptic Curve defined by y^2 + x*y = x^3 + 1 over Finite Field in a of size 2^6 to Elliptic Curve defined by y^2 + x*y = x^3 + x over Finite Field in a of size 2^6]1876sage: E = EllipticCurve_from_j(GF(3^12,'a')(2))1877sage: isogenies_prime_degree_general(E, 17)1878[Isogeny of degree 17 from Elliptic Curve defined by y^2 = x^3 + 2*x^2 + 2 over Finite Field in a of size 3^12 to Elliptic Curve defined by y^2 = x^3 + 2*x^2 + 2*x over Finite Field in a of size 3^12, Isogeny of degree 17 from Elliptic Curve defined by y^2 = x^3 + 2*x^2 + 2 over Finite Field in a of size 3^12 to Elliptic Curve defined by y^2 = x^3 + 2*x^2 + x + 2 over Finite Field in a of size 3^12]1879sage: E = EllipticCurve('50a1')1880sage: isogenies_prime_degree_general(E, 3)1881[Isogeny of degree 3 from Elliptic Curve defined by y^2 + x*y + y = x^3 - x - 2 over Rational Field to Elliptic Curve defined by y^2 + x*y + y = x^3 - 126*x - 552 over Rational Field]1882sage: isogenies_prime_degree_general(E, 5)1883[Isogeny of degree 5 from Elliptic Curve defined by y^2 + x*y + y = x^3 - x - 2 over Rational Field to Elliptic Curve defined by y^2 + x*y + y = x^3 - 76*x + 298 over Rational Field]1884sage: E = EllipticCurve([0, 0, 1, -1862, -30956])1885sage: isogenies_prime_degree_general(E, 19)1886[Isogeny of degree 19 from Elliptic Curve defined by y^2 + y = x^3 - 1862*x - 30956 over Rational Field to Elliptic Curve defined by y^2 + y = x^3 - 672182*x + 212325489 over Rational Field]1887sage: E = EllipticCurve([0, -1, 0, -6288, 211072])1888sage: isogenies_prime_degree_general(E, 37) # long time (10s)1889[Isogeny of degree 37 from Elliptic Curve defined by y^2 = x^3 - x^2 - 6288*x + 211072 over Rational Field to Elliptic Curve defined by y^2 = x^3 - x^2 - 163137088*x - 801950801728 over Rational Field]18901891sage: E = EllipticCurve([-3440, 77658])1892sage: isogenies_prime_degree_general(E, 43) # long time (16s)1893[Isogeny of degree 43 from Elliptic Curve defined by y^2 = x^3 - 3440*x + 77658 over Rational Field to Elliptic Curve defined by y^2 = x^3 - 6360560*x - 6174354606 over Rational Field]18941895Isogenies of degree equal to the characteristic are computed (but1896only the separable isogeny). In the following example we consider1897an elliptic curve which is supersingular in characteristic 21898only::18991900sage: from sage.schemes.elliptic_curves.isogeny_small_degree import isogenies_prime_degree_general1901sage: ainvs = (0,1,1,-1,-1)1902sage: for l in prime_range(50):1903....: E = EllipticCurve(GF(l),ainvs)1904....: isogenies_prime_degree_general(E,l)1905[]1906[Isogeny of degree 3 from Elliptic Curve defined by y^2 + y = x^3 + x^2 + 2*x + 2 over Finite Field of size 3 to Elliptic Curve defined by y^2 + y = x^3 + x^2 + x over Finite Field of size 3]1907[Isogeny of degree 5 from Elliptic Curve defined by y^2 + y = x^3 + x^2 + 4*x + 4 over Finite Field of size 5 to Elliptic Curve defined by y^2 + y = x^3 + x^2 + 4*x + 4 over Finite Field of size 5]1908[Isogeny of degree 7 from Elliptic Curve defined by y^2 + y = x^3 + x^2 + 6*x + 6 over Finite Field of size 7 to Elliptic Curve defined by y^2 + y = x^3 + x^2 + 4 over Finite Field of size 7]1909[Isogeny of degree 11 from Elliptic Curve defined by y^2 + y = x^3 + x^2 + 10*x + 10 over Finite Field of size 11 to Elliptic Curve defined by y^2 + y = x^3 + x^2 + x + 1 over Finite Field of size 11]1910[Isogeny of degree 13 from Elliptic Curve defined by y^2 + y = x^3 + x^2 + 12*x + 12 over Finite Field of size 13 to Elliptic Curve defined by y^2 + y = x^3 + x^2 + 12*x + 12 over Finite Field of size 13]1911[Isogeny of degree 17 from Elliptic Curve defined by y^2 + y = x^3 + x^2 + 16*x + 16 over Finite Field of size 17 to Elliptic Curve defined by y^2 + y = x^3 + x^2 + 15 over Finite Field of size 17]1912[Isogeny of degree 19 from Elliptic Curve defined by y^2 + y = x^3 + x^2 + 18*x + 18 over Finite Field of size 19 to Elliptic Curve defined by y^2 + y = x^3 + x^2 + 3*x + 12 over Finite Field of size 19]1913[Isogeny of degree 23 from Elliptic Curve defined by y^2 + y = x^3 + x^2 + 22*x + 22 over Finite Field of size 23 to Elliptic Curve defined by y^2 + y = x^3 + x^2 + 22*x + 22 over Finite Field of size 23]1914[Isogeny of degree 29 from Elliptic Curve defined by y^2 + y = x^3 + x^2 + 28*x + 28 over Finite Field of size 29 to Elliptic Curve defined by y^2 + y = x^3 + x^2 + 7*x + 27 over Finite Field of size 29]1915[Isogeny of degree 31 from Elliptic Curve defined by y^2 + y = x^3 + x^2 + 30*x + 30 over Finite Field of size 31 to Elliptic Curve defined by y^2 + y = x^3 + x^2 + 15*x + 16 over Finite Field of size 31]1916[Isogeny of degree 37 from Elliptic Curve defined by y^2 + y = x^3 + x^2 + 36*x + 36 over Finite Field of size 37 to Elliptic Curve defined by y^2 + y = x^3 + x^2 + 16*x + 17 over Finite Field of size 37]1917[Isogeny of degree 41 from Elliptic Curve defined by y^2 + y = x^3 + x^2 + 40*x + 40 over Finite Field of size 41 to Elliptic Curve defined by y^2 + y = x^3 + x^2 + 33*x + 16 over Finite Field of size 41]1918[Isogeny of degree 43 from Elliptic Curve defined by y^2 + y = x^3 + x^2 + 42*x + 42 over Finite Field of size 43 to Elliptic Curve defined by y^2 + y = x^3 + x^2 + 36 over Finite Field of size 43]1919[Isogeny of degree 47 from Elliptic Curve defined by y^2 + y = x^3 + x^2 + 46*x + 46 over Finite Field of size 47 to Elliptic Curve defined by y^2 + y = x^3 + x^2 + 42*x + 34 over Finite Field of size 47]19201921"""1922if not l.is_prime():1923raise ValueError("%s is not prime."%l)1924if l==2:1925return isogenies_2(E)1926if l==3:1927return isogenies_3(E)19281929psi_l = E.division_polynomial(l)1930factors = [h for h,e in psi_l.factor() if (l-1)/2 % h.degree() == 0]1931a = _least_semi_primitive(l)1932m = E.multiplication_by_m(a, x_only=True)1933F = psi_l.parent()1934x = F.gen()1935ker = []1936from sage.rings.arith import gcd1937from sage.misc.all import prod1938def mult(f):1939return gcd(F(f(m(x)).numerator()),psi_l).monic()1940while len(factors) > 0:1941f = factors[0]1942factors.remove(f)1943d = f.degree()1944S = [f]1945for i in range((l-1)/(2*d)-1):1946g = mult(S[i])1947S.append(g)1948if g in factors:1949factors.remove(g)1950if mult(S[-1]) == f:1951ker.append(prod(S))1952return [E.isogeny(k) for k in ker]19531954def isogenies_prime_degree(E, l):1955"""1956Returns a list of ``l`` -isogenies with domain ``E``.19571958INPUT:19591960- ``E`` -- an elliptic curve.19611962- ``l`` -- a prime.19631964OUTPUT:19651966(list) a list of all isogenies of degree `l`. If the1967characteristic is `l` then only separable isogenies are1968constructed.19691970EXAMPLES::19711972sage: from sage.schemes.elliptic_curves.isogeny_small_degree import isogenies_prime_degree1973sage: E = EllipticCurve_from_j(GF(2^6,'a')(1))1974sage: isogenies_prime_degree(E, 7)1975[Isogeny of degree 7 from Elliptic Curve defined by y^2 + x*y = x^3 + 1 over Finite Field in a of size 2^6 to Elliptic Curve defined by y^2 + x*y = x^3 + x over Finite Field in a of size 2^6]1976sage: E = EllipticCurve_from_j(GF(3^12,'a')(2))1977sage: isogenies_prime_degree(E, 17)1978[Isogeny of degree 17 from Elliptic Curve defined by y^2 = x^3 + 2*x^2 + 2 over Finite Field in a of size 3^12 to Elliptic Curve defined by y^2 = x^3 + 2*x^2 + 2*x over Finite Field in a of size 3^12, Isogeny of degree 17 from Elliptic Curve defined by y^2 = x^3 + 2*x^2 + 2 over Finite Field in a of size 3^12 to Elliptic Curve defined by y^2 = x^3 + 2*x^2 + x + 2 over Finite Field in a of size 3^12]1979sage: E = EllipticCurve('50a1')1980sage: isogenies_prime_degree(E, 3)1981[Isogeny of degree 3 from Elliptic Curve defined by y^2 + x*y + y = x^3 - x - 2 over Rational Field to Elliptic Curve defined by y^2 + x*y + y = x^3 - 126*x - 552 over Rational Field]1982sage: isogenies_prime_degree(E, 5)1983[Isogeny of degree 5 from Elliptic Curve defined by y^2 + x*y + y = x^3 - x - 2 over Rational Field to Elliptic Curve defined by y^2 + x*y + y = x^3 - 76*x + 298 over Rational Field]1984sage: E = EllipticCurve([0, 0, 1, -1862, -30956])1985sage: isogenies_prime_degree(E, 19)1986[Isogeny of degree 19 from Elliptic Curve defined by y^2 + y = x^3 - 1862*x - 30956 over Rational Field to Elliptic Curve defined by y^2 + y = x^3 - 672182*x + 212325489 over Rational Field]1987sage: E = EllipticCurve([0, -1, 0, -6288, 211072])1988sage: isogenies_prime_degree(E, 37)1989[Isogeny of degree 37 from Elliptic Curve defined by y^2 = x^3 - x^2 - 6288*x + 211072 over Rational Field to Elliptic Curve defined by y^2 = x^3 - x^2 - 163137088*x - 801950801728 over Rational Field]199019911992Isogenies of degree equal to the characteristic are computed (but1993only the separable isogeny). In the following example we consider1994an elliptic curve which is supersingular in characteristic 21995only::19961997sage: from sage.schemes.elliptic_curves.isogeny_small_degree import isogenies_prime_degree1998sage: ainvs = (0,1,1,-1,-1)1999sage: for l in prime_range(50):2000....: E = EllipticCurve(GF(l),ainvs)2001....: isogenies_prime_degree(E,l)2002[]2003[Isogeny of degree 3 from Elliptic Curve defined by y^2 + y = x^3 + x^2 + 2*x + 2 over Finite Field of size 3 to Elliptic Curve defined by y^2 + y = x^3 + x^2 + x over Finite Field of size 3]2004[Isogeny of degree 5 from Elliptic Curve defined by y^2 + y = x^3 + x^2 + 4*x + 4 over Finite Field of size 5 to Elliptic Curve defined by y^2 + y = x^3 + x^2 + 4*x + 4 over Finite Field of size 5]2005[Isogeny of degree 7 from Elliptic Curve defined by y^2 + y = x^3 + x^2 + 6*x + 6 over Finite Field of size 7 to Elliptic Curve defined by y^2 + y = x^3 + x^2 + 4 over Finite Field of size 7]2006[Isogeny of degree 11 from Elliptic Curve defined by y^2 + y = x^3 + x^2 + 10*x + 10 over Finite Field of size 11 to Elliptic Curve defined by y^2 + y = x^3 + x^2 + x + 1 over Finite Field of size 11]2007[Isogeny of degree 13 from Elliptic Curve defined by y^2 + y = x^3 + x^2 + 12*x + 12 over Finite Field of size 13 to Elliptic Curve defined by y^2 + y = x^3 + x^2 + 12*x + 12 over Finite Field of size 13]2008[Isogeny of degree 17 from Elliptic Curve defined by y^2 + y = x^3 + x^2 + 16*x + 16 over Finite Field of size 17 to Elliptic Curve defined by y^2 + y = x^3 + x^2 + 15 over Finite Field of size 17]2009[Isogeny of degree 19 from Elliptic Curve defined by y^2 + y = x^3 + x^2 + 18*x + 18 over Finite Field of size 19 to Elliptic Curve defined by y^2 + y = x^3 + x^2 + 3*x + 12 over Finite Field of size 19]2010[Isogeny of degree 23 from Elliptic Curve defined by y^2 + y = x^3 + x^2 + 22*x + 22 over Finite Field of size 23 to Elliptic Curve defined by y^2 + y = x^3 + x^2 + 22*x + 22 over Finite Field of size 23]2011[Isogeny of degree 29 from Elliptic Curve defined by y^2 + y = x^3 + x^2 + 28*x + 28 over Finite Field of size 29 to Elliptic Curve defined by y^2 + y = x^3 + x^2 + 7*x + 27 over Finite Field of size 29]2012[Isogeny of degree 31 from Elliptic Curve defined by y^2 + y = x^3 + x^2 + 30*x + 30 over Finite Field of size 31 to Elliptic Curve defined by y^2 + y = x^3 + x^2 + 15*x + 16 over Finite Field of size 31]2013[Isogeny of degree 37 from Elliptic Curve defined by y^2 + y = x^3 + x^2 + 36*x + 36 over Finite Field of size 37 to Elliptic Curve defined by y^2 + y = x^3 + x^2 + 16*x + 17 over Finite Field of size 37]2014[Isogeny of degree 41 from Elliptic Curve defined by y^2 + y = x^3 + x^2 + 40*x + 40 over Finite Field of size 41 to Elliptic Curve defined by y^2 + y = x^3 + x^2 + 33*x + 16 over Finite Field of size 41]2015[Isogeny of degree 43 from Elliptic Curve defined by y^2 + y = x^3 + x^2 + 42*x + 42 over Finite Field of size 43 to Elliptic Curve defined by y^2 + y = x^3 + x^2 + 36 over Finite Field of size 43]2016[Isogeny of degree 47 from Elliptic Curve defined by y^2 + y = x^3 + x^2 + 46*x + 46 over Finite Field of size 47 to Elliptic Curve defined by y^2 + y = x^3 + x^2 + 42*x + 34 over Finite Field of size 47]20172018Note that the computation is faster for degrees equal to one of2019the genus 0 primes (2, 3, 5, 7, 13) or one of the hyperelliptic2020primes (11, 17, 19, 23, 29, 31, 41, 47, 59, 71) than when the2021generic code must be used::20222023sage: E = EllipticCurve(GF(101), [-3440, 77658])2024sage: E.isogenies_prime_degree(71) # fast2025[]2026sage: E.isogenies_prime_degree(73) # not tested (very long time: 32s)2027[]202820292030"""2031if not l.is_prime():2032raise ValueError("%s is not prime."%l)2033if l==2:2034return isogenies_2(E)2035if l==3:2036return isogenies_3(E)20372038p = E.base_ring().characteristic()2039if l==p:2040return isogenies_prime_degree_general(E,l)20412042if l in [5,7,13] and not p in [2,3]:2043return isogenies_prime_degree_genus_0(E,l)20442045if l in hyperelliptic_primes and not p in [2,3]:2046return isogenies_prime_degree_genus_plus_0(E,l)20472048return isogenies_prime_degree_general(E,l)204920502051