An isogeny between two elliptic curves
and
is a morphism of curves that sends the origin of
to the
origin of
. Such a morphism is automatically a morphism of group
schemes and the kernel is a finite subgroup scheme of
. Such a
subscheme can either be given by a list of generators, which have to
be torsion points, or by a polynomial in the coordinate
of the
Weierstrass equation of
.
The usual way to create and work with isogenies is illustrated with the following example:
sage: k = GF(11)
sage: E = EllipticCurve(k,[1,1])
sage: Q = E(6,5)
sage: phi = E.isogeny(Q)
sage: phi
Isogeny of degree 7 from Elliptic Curve defined by y^2 = x^3 + x + 1 over Finite Field of size 11 to Elliptic Curve defined by y^2 = x^3 + 7*x + 8 over Finite Field of size 11
sage: P = E(4,5)
sage: phi(P)
(10 : 0 : 1)
sage: phi.codomain()
Elliptic Curve defined by y^2 = x^3 + 7*x + 8 over Finite Field of size 11
sage: phi.rational_maps()
((x^7 + 4*x^6 - 3*x^5 - 2*x^4 - 3*x^3 + 3*x^2 + x - 2)/(x^6 + 4*x^5 - 4*x^4 - 5*x^3 + 5*x^2), (x^9*y - 5*x^8*y - x^7*y + x^5*y - x^4*y - 5*x^3*y - 5*x^2*y - 2*x*y - 5*y)/(x^9 - 5*x^8 + 4*x^6 - 3*x^4 + 2*x^3))
The functions directly accessible from an elliptic curve E over a field are isogeny and isogeny_codomain.
The most useful functions that apply to isogenies are
Warning
Only cyclic, separable isogenies are implemented (except for [2]). Some algorithms may need the isogeny to be normalized.
AUTHORS:
Bases: sage.categories.morphism.Morphism
Class Implementing Isogenies of Elliptic Curves
This class implements cyclic, separable, normalized isogenies of elliptic curves.
Several different algorithms for computing isogenies are available. These include:
INPUT:
EXAMPLES:
A simple example of creating an isogeny of a field of small characteristic:
sage: E = EllipticCurve(GF(7), [0,0,0,1,0])
sage: phi = EllipticCurveIsogeny(E, E((0,0)) ); phi
Isogeny of degree 2 from Elliptic Curve defined by y^2 = x^3 + x over Finite Field of size 7 to Elliptic Curve defined by y^2 = x^3 + 3*x over Finite Field of size 7
sage: phi.degree() == 2
True
sage: phi.kernel_polynomial()
x
sage: phi.rational_maps()
((x^2 + 1)/x, (x^2*y - y)/x^2)
sage: phi == loads(dumps(phi)) # known bug
True
A more complicated example of a characteristic 2 field:
sage: E = EllipticCurve(GF(2^4,'alpha'), [0,0,1,0,1])
sage: P = E((1,1))
sage: phi_v = EllipticCurveIsogeny(E, P); phi_v
Isogeny of degree 3 from Elliptic Curve defined by y^2 + y = x^3 + 1 over Finite Field in alpha of size 2^4 to Elliptic Curve defined by y^2 + y = x^3 over Finite Field in alpha of size 2^4
sage: phi_ker_poly = phi_v.kernel_polynomial()
sage: phi_ker_poly
x + 1
sage: ker_poly_list = phi_ker_poly.list()
sage: phi_k = EllipticCurveIsogeny(E, ker_poly_list)
sage: phi_k == phi_v
True
sage: phi_k.rational_maps()
((x^3 + x + 1)/(x^2 + 1), (x^3*y + x^2*y + x*y + x + y)/(x^3 + x^2 + x + 1))
sage: phi_v.rational_maps()
((x^3 + x + 1)/(x^2 + 1), (x^3*y + x^2*y + x*y + x + y)/(x^3 + x^2 + x + 1))
sage: phi_k.degree() == phi_v.degree() == 3
True
sage: phi_k.is_separable()
True
sage: phi_v(E(0))
(0 : 1 : 0)
sage: alpha = E.base_field().gen()
sage: Q = E((0, alpha*(alpha + 1)))
sage: phi_v(Q)
(1 : alpha^2 + alpha : 1)
sage: phi_v(P) == phi_k(P)
True
sage: phi_k(P) == phi_v.codomain()(0)
True
We can create an isogeny that has kernel equal to the full 2 torsion:
sage: E = EllipticCurve(GF(3), [0,0,0,1,1])
sage: ker_list = E.division_polynomial(2).list()
sage: phi = EllipticCurveIsogeny(E, ker_list); phi
Isogeny of degree 4 from Elliptic Curve defined by y^2 = x^3 + x + 1 over Finite Field of size 3 to Elliptic Curve defined by y^2 = x^3 + x + 1 over Finite Field of size 3
sage: phi(E(0))
(0 : 1 : 0)
sage: phi(E((0,1)))
(1 : 0 : 1)
sage: phi(E((0,2)))
(1 : 0 : 1)
sage: phi(E((1,0)))
(0 : 1 : 0)
sage: phi.degree()
4
We can also create trivial isogenies with the trivial kernel:
sage: E = EllipticCurve(GF(17), [11, 11, 4, 12, 10])
sage: phi_v = EllipticCurveIsogeny(E, E(0))
sage: phi_v.degree()
1
sage: phi_v.rational_maps()
(x, y)
sage: E == phi_v.codomain()
True
sage: P = E.random_point()
sage: phi_v(P) == P
True
sage: E = EllipticCurve(GF(31), [23, 1, 22, 7, 18])
sage: phi_k = EllipticCurveIsogeny(E, [1]); phi_k
Isogeny of degree 1 from Elliptic Curve defined by y^2 + 23*x*y + 22*y = x^3 + x^2 + 7*x + 18 over Finite Field of size 31 to Elliptic Curve defined by y^2 + 23*x*y + 22*y = x^3 + x^2 + 7*x + 18 over Finite Field of size 31
sage: phi_k.degree()
1
sage: phi_k.rational_maps()
(x, y)
sage: phi_k.codomain() == E
True
sage: phi_k.kernel_polynomial()
1
sage: P = E.random_point(); P == phi_k(P)
True
Velu and Kohel also work in characteristic 0:
sage: E = EllipticCurve(QQ, [0,0,0,3,4])
sage: P_list = E.torsion_points()
sage: phi = EllipticCurveIsogeny(E, P_list); phi
Isogeny of degree 2 from Elliptic Curve defined by y^2 = x^3 + 3*x + 4 over Rational Field to Elliptic Curve defined by y^2 = x^3 - 27*x + 46 over Rational Field
sage: P = E((0,2))
sage: phi(P)
(6 : -10 : 1)
sage: phi_ker_poly = phi.kernel_polynomial()
sage: phi_ker_poly
x + 1
sage: ker_poly_list = phi_ker_poly.list()
sage: phi_k = EllipticCurveIsogeny(E, ker_poly_list); phi_k
Isogeny of degree 2 from Elliptic Curve defined by y^2 = x^3 + 3*x + 4 over Rational Field to Elliptic Curve defined by y^2 = x^3 - 27*x + 46 over Rational Field
sage: phi_k(P) == phi(P)
True
sage: phi_k == phi
True
sage: phi_k.degree()
2
sage: phi_k.is_separable()
True
A more complicated example over the rationals (of odd degree):
sage: E = EllipticCurve('11a1')
sage: P_list = E.torsion_points()
sage: phi_v = EllipticCurveIsogeny(E, P_list); phi_v
Isogeny of degree 5 from Elliptic Curve defined by y^2 + y = x^3 - x^2 - 10*x - 20 over Rational Field to Elliptic Curve defined by y^2 + y = x^3 - x^2 - 7820*x - 263580 over Rational Field
sage: P = E((16,-61))
sage: phi_v(P)
(0 : 1 : 0)
sage: ker_poly = phi_v.kernel_polynomial(); ker_poly
x^2 - 21*x + 80
sage: ker_poly_list = ker_poly.list()
sage: phi_k = EllipticCurveIsogeny(E, ker_poly_list); phi_k
Isogeny of degree 5 from Elliptic Curve defined by y^2 + y = x^3 - x^2 - 10*x - 20 over Rational Field to Elliptic Curve defined by y^2 + y = x^3 - x^2 - 7820*x - 263580 over Rational Field
sage: phi_k == phi_v
True
sage: phi_v(P) == phi_k(P)
True
sage: phi_k.is_separable()
True
We can also do this same example over the number field defined by
the irreducible two torsion polynomial of :
sage: E = EllipticCurve('11a1')
sage: P_list = E.torsion_points()
sage: K.<alpha> = NumberField(x^3 - 2* x^2 - 40*x - 158)
sage: EK = E.change_ring(K)
sage: P_list = [EK(P) for P in P_list]
sage: phi_v = EllipticCurveIsogeny(EK, P_list); phi_v
Isogeny of degree 5 from Elliptic Curve defined by y^2 + y = x^3 + (-1)*x^2 + (-10)*x + (-20) over Number Field in alpha with defining polynomial x^3 - 2*x^2 - 40*x - 158 to Elliptic Curve defined by y^2 + y = x^3 + (-1)*x^2 + (-7820)*x + (-263580) over Number Field in alpha with defining polynomial x^3 - 2*x^2 - 40*x - 158
sage: P = EK((alpha/2,-1/2))
sage: phi_v(P)
(122/121*alpha^2 + 1633/242*alpha - 3920/121 : -1/2 : 1)
sage: ker_poly = phi_v.kernel_polynomial()
sage: ker_poly
x^2 - 21*x + 80
sage: ker_poly_list = ker_poly.list()
sage: phi_k = EllipticCurveIsogeny(EK, ker_poly_list)
sage: phi_k
Isogeny of degree 5 from Elliptic Curve defined by y^2 + y = x^3 + (-1)*x^2 + (-10)*x + (-20) over Number Field in alpha with defining polynomial x^3 - 2*x^2 - 40*x - 158 to Elliptic Curve defined by y^2 + y = x^3 + (-1)*x^2 + (-7820)*x + (-263580) over Number Field in alpha with defining polynomial x^3 - 2*x^2 - 40*x - 158
sage: phi_v == phi_k
True
sage: phi_k(P) == phi_v(P)
True
sage: phi_k == phi_v
True
sage: phi_k.degree()
5
sage: phi_v.is_separable()
True
The following example shows how to specify an isogeny from domain and codomain:
sage: E = EllipticCurve('11a1')
sage: R.<x> = QQ[]
sage: f = x^2 - 21*x + 80
sage: phi = E.isogeny(f)
sage: E2 = phi.codomain()
sage: phi_s = EllipticCurveIsogeny(E, None, E2, 5)
sage: phi_s
Isogeny of degree 5 from Elliptic Curve defined by y^2 + y = x^3 - x^2 - 10*x - 20 over Rational Field to Elliptic Curve defined by y^2 + y = x^3 - x^2 - 7820*x - 263580 over Rational Field
sage: phi_s == phi
True
sage: phi_s.rational_maps() == phi.rational_maps()
True
However only cyclic normalized isogenies can be constructed this way. So it won’t find the isogeny [3]:
sage: E.isogeny(None, codomain=E,degree=9)
Traceback (most recent call last):
...
ValueError: The two curves are not linked by a cyclic normalized isogeny of degree 9
Also the presumed isogeny between the domain and codomain must be normalized:
sage: E2.isogeny(None,codomain=E,degree=5)
Traceback (most recent call last):
...
ValueError: The two curves are not linked by a cyclic normalized isogeny of degree 5
sage: phihat = phi.dual(); phihat
Isogeny of degree 5 from Elliptic Curve defined by y^2 + y = x^3 - x^2 - 7820*x - 263580 over Rational Field to Elliptic Curve defined by y^2 + y = x^3 - x^2 - 10*x - 20 over Rational Field
sage: phihat.is_normalized()
False
Here an example of a construction of a endomorphisms with cyclic kernel on a CM-curve:
sage: K.<i> = NumberField(x^2+1)
sage: E = EllipticCurve(K, [1,0])
sage: RK.<X> = K[]
sage: f = X^2 - 2/5*i + 1/5
sage: phi= E.isogeny(f)
sage: isom = phi.codomain().isomorphism_to(E)
sage: phi.set_post_isomorphism(isom)
sage: phi.codomain() == phi.domain()
True
sage: phi.rational_maps()
(((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)), ((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)))
Domain and codomain tests (see trac ticket #12880):
sage: E = EllipticCurve(QQ, [0,0,0,1,0])
sage: phi = EllipticCurveIsogeny(E, E(0,0))
sage: phi.domain() == E
True
sage: phi.codomain()
Elliptic Curve defined by y^2 = x^3 - 4*x over Rational Field
sage: E = EllipticCurve(GF(31), [1,0,0,1,2])
sage: phi = EllipticCurveIsogeny(E, [17, 1])
sage: phi.domain()
Elliptic Curve defined by y^2 + x*y = x^3 + x + 2 over Finite Field of size 31
sage: phi.codomain()
Elliptic Curve defined by y^2 + x*y = x^3 + 24*x + 6 over Finite Field of size 31
Composition tests (see trac ticket #16245):
sage: E = EllipticCurve(j=GF(7)(0))
sage: phi = E.isogeny([E(0), E((0,1)), E((0,-1))]); phi
Isogeny of degree 3 from Elliptic Curve defined by y^2 = x^3 + 1 over Finite Field of size 7 to Elliptic Curve defined by y^2 = x^3 + 1 over Finite Field of size 7
sage: phi2 = phi * phi; phi2
Composite map:
From: Elliptic Curve defined by y^2 = x^3 + 1 over Finite Field of size 7
To: Elliptic Curve defined by y^2 = x^3 + 1 over Finite Field of size 7
Defn: Isogeny of degree 3 from Elliptic Curve defined by y^2 = x^3 + 1 over Finite Field of size 7 to Elliptic Curve defined by y^2 = x^3 + 1 over Finite Field of size 7
then
Isogeny of degree 3 from Elliptic Curve defined by y^2 = x^3 + 1 over Finite Field of size 7 to Elliptic Curve defined by y^2 = x^3 + 1 over Finite Field of size 7
Examples over relative number fields used not to work (see trac ticket #16779):
sage: pol26 = hilbert_class_polynomial(-4*26)
sage: pol = NumberField(pol26,'a').optimized_representation()[0].polynomial()
sage: K.<a> = NumberField(pol)
sage: j = pol26.roots(K)[0][0]
sage: E = EllipticCurve(j=j)
sage: L.<b> = K.extension(x^2+26)
sage: EL = E.change_ring(L)
sage: iso2 = EL.isogenies_prime_degree(2); len(iso2)
1
sage: iso3 = EL.isogenies_prime_degree(3); len(iso3)
2
Examples over function fields used not to work (see trac ticket #11327):
sage: F.<t> = FunctionField(QQ)
sage: E = EllipticCurve([0,0,0,-t^2,0])
sage: isogs = E.isogenies_prime_degree(2)
sage: isogs[0]
Isogeny of degree 2 from Elliptic Curve defined by y^2 = x^3 + (-t^2)*x over Rational function field in t over Rational Field to Elliptic Curve defined by y^2 = x^3 + 4*t^2*x over Rational function field in t over Rational Field
sage: isogs[0].rational_maps()
((x^2 - t^2)/x, (x^3*y + t^2*x*y)/x^3)
sage: duals = [phi.dual() for phi in isogs]
sage: duals[0]
Isogeny of degree 2 from Elliptic Curve defined by y^2 = x^3 + 4*t^2*x over Rational function field in t over Rational Field to Elliptic Curve defined by y^2 = x^3 + (-t^2)*x over Rational function field in t over Rational Field
sage: duals[0].rational_maps()
((1/4*x^2 + t^2)/x, (1/8*x^3*y + (-1/2*t^2)*x*y)/x^3)
sage: duals[0]
Isogeny of degree 2 from Elliptic Curve defined by y^2 = x^3 + 4*t^2*x over Rational function field in t over Rational Field to Elliptic Curve defined by y^2 = x^3 + (-t^2)*x over Rational function field in t over Rational Field
Returns the degree of this isogeny.
EXAMPLES:
sage: E = EllipticCurve(QQ, [0,0,0,1,0])
sage: phi = EllipticCurveIsogeny(E, E((0,0)))
sage: phi.degree()
2
sage: phi = EllipticCurveIsogeny(E, [0,1,0,1])
sage: phi.degree()
4
sage: E = EllipticCurve(GF(31), [1,0,0,1,2])
sage: phi = EllipticCurveIsogeny(E, [17, 1])
sage: phi.degree()
3
Return the isogeny dual to this isogeny.
Note
If is the given isogeny, then the
dual is by definition the unique isogeny
such that the compositions
and
are
the multiplication
by the degree of
on
and
respectively.
EXAMPLES:
sage: E = EllipticCurve('11a1')
sage: R.<x> = QQ[]
sage: f = x^2 - 21*x + 80
sage: phi = EllipticCurveIsogeny(E, f)
sage: phi_hat = phi.dual()
sage: phi_hat.domain() == phi.codomain()
True
sage: phi_hat.codomain() == phi.domain()
True
sage: (X, Y) = phi.rational_maps()
sage: (Xhat, Yhat) = phi_hat.rational_maps()
sage: Xm = Xhat.subs(x=X, y=Y)
sage: Ym = Yhat.subs(x=X, y=Y)
sage: (Xm, Ym) == E.multiplication_by_m(5)
True
sage: E = EllipticCurve(GF(37), [0,0,0,1,8])
sage: R.<x> = GF(37)[]
sage: f = x^3 + x^2 + 28*x + 33
sage: phi = EllipticCurveIsogeny(E, f)
sage: phi_hat = phi.dual()
sage: phi_hat.codomain() == phi.domain()
True
sage: phi_hat.domain() == phi.codomain()
True
sage: (X, Y) = phi.rational_maps()
sage: (Xhat, Yhat) = phi_hat.rational_maps()
sage: Xm = Xhat.subs(x=X, y=Y)
sage: Ym = Yhat.subs(x=X, y=Y)
sage: (Xm, Ym) == E.multiplication_by_m(7)
True
sage: E = EllipticCurve(GF(31), [0,0,0,1,8])
sage: R.<x> = GF(31)[]
sage: f = x^2 + 17*x + 29
sage: phi = EllipticCurveIsogeny(E, f)
sage: phi_hat = phi.dual()
sage: phi_hat.codomain() == phi.domain()
True
sage: phi_hat.domain() == phi.codomain()
True
sage: (X, Y) = phi.rational_maps()
sage: (Xhat, Yhat) = phi_hat.rational_maps()
sage: Xm = Xhat.subs(x=X, y=Y)
sage: Ym = Yhat.subs(x=X, y=Y)
sage: (Xm, Ym) == E.multiplication_by_m(5)
True
Test (for trac ticket 7096):
sage: E = EllipticCurve('11a1')
sage: phi = E.isogeny(E(5,5))
sage: phi.dual().dual() == phi
True
sage: k = GF(103)
sage: E = EllipticCurve(k,[11,11])
sage: phi = E.isogeny(E(4,4))
sage: phi
Isogeny of degree 5 from Elliptic Curve defined by y^2 = x^3 + 11*x + 11 over Finite Field of size 103 to Elliptic Curve defined by y^2 = x^3 + 25*x + 80 over Finite Field of size 103
sage: from sage.schemes.elliptic_curves.weierstrass_morphism import WeierstrassIsomorphism
sage: phi.set_post_isomorphism(WeierstrassIsomorphism(phi.codomain(),(5,0,1,2)))
sage: phi.dual().dual() == phi
True
sage: E = EllipticCurve(GF(103),[1,0,0,1,-1])
sage: phi = E.isogeny(E(60,85))
sage: phi.dual()
Isogeny of degree 7 from Elliptic Curve defined by y^2 + x*y = x^3 + 84*x + 34 over Finite Field of size 103 to Elliptic Curve defined by y^2 + x*y = x^3 + x + 102 over Finite Field of size 103
Return the formal isogeny as a power series in the variable
on the domain curve.
INPUT:
EXAMPLES:
sage: E = EllipticCurve(GF(13),[1,7])
sage: phi = E.isogeny(E(10,4))
sage: phi.formal()
t + 12*t^13 + 2*t^17 + 8*t^19 + 2*t^21 + O(t^23)
sage: E = EllipticCurve([0,1])
sage: phi = E.isogeny(E(2,3))
sage: phi.formal(prec=10)
t + 54*t^5 + 255*t^7 + 2430*t^9 + 19278*t^11 + O(t^13)
sage: E = EllipticCurve('11a2')
sage: R.<x> = QQ[]
sage: phi = E.isogeny(x^2 + 101*x + 12751/5)
sage: phi.formal(prec=7)
t - 2724/5*t^5 + 209046/5*t^7 - 4767/5*t^8 + 29200946/5*t^9 + O(t^10)
Return the post-isomorphism of this isogeny, or None.
EXAMPLES:
sage: E = EllipticCurve(j=GF(31)(0))
sage: R.<x> = GF(31)[]
sage: phi = EllipticCurveIsogeny(E, x+18)
sage: phi.get_post_isomorphism()
sage: from sage.schemes.elliptic_curves.weierstrass_morphism import WeierstrassIsomorphism
sage: isom = WeierstrassIsomorphism(phi.codomain(), (6,8,10,12))
sage: phi.set_post_isomorphism(isom)
sage: isom == phi.get_post_isomorphism()
True
sage: E = EllipticCurve(GF(83), [1,0,1,1,0])
sage: R.<x> = GF(83)[]; f = x+24
sage: phi = EllipticCurveIsogeny(E, f)
sage: E2 = phi.codomain()
sage: phi2 = EllipticCurveIsogeny(E, None, E2, 2)
sage: phi2.get_post_isomorphism()
Generic morphism:
From: Abelian group of points on Elliptic Curve defined by y^2 = x^3 + 65*x + 69 over Finite Field of size 83
To: Abelian group of points on Elliptic Curve defined by y^2 + x*y + y = x^3 + 4*x + 16 over Finite Field of size 83
Via: (u,r,s,t) = (1, 7, 42, 42)
Return the pre-isomorphism of this isogeny, or None.
EXAMPLES:
sage: E = EllipticCurve(GF(31), [1,1,0,1,-1])
sage: R.<x> = GF(31)[]
sage: f = x^3 + 9*x^2 + x + 30
sage: phi = EllipticCurveIsogeny(E, f)
sage: phi.get_post_isomorphism()
sage: Epr = E.short_weierstrass_model()
sage: isom = Epr.isomorphism_to(E)
sage: phi.set_pre_isomorphism(isom)
sage: isom == phi.get_pre_isomorphism()
True
sage: E = EllipticCurve(GF(83), [1,0,1,1,0])
sage: R.<x> = GF(83)[]; f = x+24
sage: phi = EllipticCurveIsogeny(E, f)
sage: E2 = phi.codomain()
sage: phi2 = EllipticCurveIsogeny(E, None, E2, 2)
sage: phi2.get_pre_isomorphism()
Generic morphism:
From: Abelian group of points on Elliptic Curve defined by y^2 + x*y + y = x^3 + x over Finite Field of size 83
To: Abelian group of points on Elliptic Curve defined by y^2 = x^3 + 62*x + 74 over Finite Field of size 83
Via: (u,r,s,t) = (1, 76, 41, 3)
Return True if and only if this isogeny has trivial kernel.
EXAMPLES:
sage: E = EllipticCurve('11a1')
sage: R.<x> = QQ[]
sage: f = x^2 + x - 29/5
sage: phi = EllipticCurveIsogeny(E, f)
sage: phi.is_injective()
False
sage: phi = EllipticCurveIsogeny(E, R(1))
sage: phi.is_injective()
True
sage: F = GF(7)
sage: E = EllipticCurve(j=F(0))
sage: phi = EllipticCurveIsogeny(E, [ E((0,-1)), E((0,1))])
sage: phi.is_injective()
False
sage: phi = EllipticCurveIsogeny(E, E(0))
sage: phi.is_injective()
True
Return whether this isogeny is normalized.
Note
An isogeny between two given
Weierstrass equations is said to be normalized if the
constant
is
in
,
where
and
are the invariant
differentials on
and
corresponding to the given
equation.
INPUT:
EXAMPLES:
sage: from sage.schemes.elliptic_curves.weierstrass_morphism import WeierstrassIsomorphism
sage: E = EllipticCurve(GF(7), [0,0,0,1,0])
sage: R.<x> = GF(7)[]
sage: phi = EllipticCurveIsogeny(E, x)
sage: phi.is_normalized()
True
sage: isom = WeierstrassIsomorphism(phi.codomain(), (3, 0, 0, 0))
sage: phi.set_post_isomorphism(isom)
sage: phi.is_normalized()
False
sage: isom = WeierstrassIsomorphism(phi.codomain(), (5, 0, 0, 0))
sage: phi.set_post_isomorphism(isom)
sage: phi.is_normalized()
True
sage: isom = WeierstrassIsomorphism(phi.codomain(), (1, 1, 1, 1))
sage: phi.set_post_isomorphism(isom)
sage: phi.is_normalized()
True
sage: F = GF(2^5, 'alpha'); alpha = F.gen()
sage: E = EllipticCurve(F, [1,0,1,1,1])
sage: R.<x> = F[]
sage: phi = EllipticCurveIsogeny(E, x+1)
sage: isom = WeierstrassIsomorphism(phi.codomain(), (alpha, 0, 0, 0))
sage: phi.is_normalized()
True
sage: phi.set_post_isomorphism(isom)
sage: phi.is_normalized()
False
sage: isom = WeierstrassIsomorphism(phi.codomain(), (1/alpha, 0, 0, 0))
sage: phi.set_post_isomorphism(isom)
sage: phi.is_normalized()
True
sage: isom = WeierstrassIsomorphism(phi.codomain(), (1, 1, 1, 1))
sage: phi.set_post_isomorphism(isom)
sage: phi.is_normalized()
True
sage: E = EllipticCurve('11a1')
sage: R.<x> = QQ[]
sage: f = x^3 - x^2 - 10*x - 79/4
sage: phi = EllipticCurveIsogeny(E, f)
sage: isom = WeierstrassIsomorphism(phi.codomain(), (2, 0, 0, 0))
sage: phi.is_normalized()
True
sage: phi.set_post_isomorphism(isom)
sage: phi.is_normalized()
False
sage: isom = WeierstrassIsomorphism(phi.codomain(), (1/2, 0, 0, 0))
sage: phi.set_post_isomorphism(isom)
sage: phi.is_normalized()
True
sage: isom = WeierstrassIsomorphism(phi.codomain(), (1, 1, 1, 1))
sage: phi.set_post_isomorphism(isom)
sage: phi.is_normalized()
True
Return whether or not this isogeny is separable.
Note
This function always returns True as currently this class only implements separable isogenies.
EXAMPLES:
sage: E = EllipticCurve(GF(17), [0,0,0,3,0])
sage: phi = EllipticCurveIsogeny(E, E((0,0)))
sage: phi.is_separable()
True
sage: E = EllipticCurve('11a1')
sage: phi = EllipticCurveIsogeny(E, E.torsion_points())
sage: phi.is_separable()
True
Return True if and only if this isogeny is surjective.
Note
This function always returns True, as a non-constant
map of algebraic curves must be surjective, and this class
does not model the constant isogeny.
EXAMPLES:
sage: E = EllipticCurve('11a1')
sage: R.<x> = QQ[]
sage: f = x^2 + x - 29/5
sage: phi = EllipticCurveIsogeny(E, f)
sage: phi.is_surjective()
True
sage: E = EllipticCurve(GF(7), [0,0,0,1,0])
sage: phi = EllipticCurveIsogeny(E, E((0,0)))
sage: phi.is_surjective()
True
sage: F = GF(2^5, 'omega')
sage: E = EllipticCurve(j=F(0))
sage: R.<x> = F[]
sage: phi = EllipticCurveIsogeny(E, x)
sage: phi.is_surjective()
True
Return whether this isogeny is zero.
Note
Currently this class does not allow zero isogenies, so this function will always return True.
EXAMPLES:
sage: E = EllipticCurve(j=GF(7)(0))
sage: phi = EllipticCurveIsogeny(E, [ E((0,1)), E((0,-1))])
sage: phi.is_zero()
False
Return the kernel polynomial of this isogeny.
EXAMPLES:
sage: E = EllipticCurve(QQ, [0,0,0,2,0])
sage: phi = EllipticCurveIsogeny(E, E((0,0)))
sage: phi.kernel_polynomial()
x
sage: E = EllipticCurve('11a1')
sage: phi = EllipticCurveIsogeny(E, E.torsion_points())
sage: phi.kernel_polynomial()
x^2 - 21*x + 80
sage: E = EllipticCurve(GF(17), [1,-1,1,-1,1])
sage: phi = EllipticCurveIsogeny(E, [1])
sage: phi.kernel_polynomial()
1
sage: E = EllipticCurve(GF(31), [0,0,0,3,0])
sage: phi = EllipticCurveIsogeny(E, [0,3,0,1])
sage: phi.kernel_polynomial()
x^3 + 3*x
Numerical Approximation inherited from Map (through morphism), nonsensical for isogenies.
EXAMPLES:
sage: E = EllipticCurve(j=GF(7)(0))
sage: phi = EllipticCurveIsogeny(E, [ E((0,1)), E((0,-1))])
sage: phi.n()
Traceback (most recent call last):
...
NotImplementedError: Numerical approximations do not make sense for Elliptic Curve Isogenies
Return the post-composition of this isogeny with left.
EXAMPLES:
sage: E = EllipticCurve(j=GF(7)(0))
sage: phi = EllipticCurveIsogeny(E, [ E((0,1)), E((0,-1))])
sage: phi.post_compose(phi)
Traceback (most recent call last):
...
NotImplementedError: post-composition of isogenies not yet implemented
Return the pre-composition of this isogeny with right.
EXAMPLES:
sage: E = EllipticCurve(j=GF(7)(0))
sage: phi = EllipticCurveIsogeny(E, [ E((0,1)), E((0,-1))])
sage: phi.pre_compose(phi)
Traceback (most recent call last):
...
NotImplementedError: pre-composition of isogenies not yet implemented
Return the pair of rational maps defining this isogeny.
Note
Both components are returned as elements of the function
field in two variables over the base field
,
though the first only involves
. To obtain the
-coordinate function as a rational function in
,
use x_rational_map().
EXAMPLES:
sage: E = EllipticCurve(QQ, [0,2,0,1,-1])
sage: phi = EllipticCurveIsogeny(E, [1])
sage: phi.rational_maps()
(x, y)
sage: E = EllipticCurve(GF(17), [0,0,0,3,0])
sage: phi = EllipticCurveIsogeny(E, E((0,0)))
sage: phi.rational_maps()
((x^2 + 3)/x, (x^2*y - 3*y)/x^2)
Modify this isogeny by postcomposing with a Weierstrass isomorphism.
EXAMPLES:
sage: E = EllipticCurve(j=GF(31)(0))
sage: R.<x> = GF(31)[]
sage: phi = EllipticCurveIsogeny(E, x+18)
sage: from sage.schemes.elliptic_curves.weierstrass_morphism import WeierstrassIsomorphism
sage: phi.set_post_isomorphism(WeierstrassIsomorphism(phi.codomain(), (6,8,10,12)))
sage: phi
Isogeny of degree 3 from Elliptic Curve defined by y^2 = x^3 + 1 over Finite Field of size 31 to Elliptic Curve defined by y^2 + 24*x*y + 7*y = x^3 + 22*x^2 + 16*x + 20 over Finite Field of size 31
sage: E = EllipticCurve(j=GF(47)(0))
sage: f = E.torsion_polynomial(3)/3
sage: phi = EllipticCurveIsogeny(E, f)
sage: E2 = phi.codomain()
sage: post_isom = E2.isomorphism_to(E)
sage: phi.set_post_isomorphism(post_isom)
sage: phi.rational_maps() == E.multiplication_by_m(3)
False
sage: phi.switch_sign()
sage: phi.rational_maps() == E.multiplication_by_m(3)
True
Example over a number field:
sage: R.<x> = QQ[]
sage: K.<a> = NumberField(x^2 + 2)
sage: E = EllipticCurve(j=K(1728))
sage: ker_list = E.torsion_points()
sage: phi = EllipticCurveIsogeny(E, ker_list)
sage: from sage.schemes.elliptic_curves.weierstrass_morphism import WeierstrassIsomorphism
sage: post_isom = WeierstrassIsomorphism(phi.codomain(), (a,2,3,5))
sage: phi
Isogeny of degree 4 from Elliptic Curve defined by y^2 = x^3 + x over Number Field in a with defining polynomial x^2 + 2 to Elliptic Curve defined by y^2 = x^3 + (-44)*x + 112 over Number Field in a with defining polynomial x^2 + 2
Modify this isogeny by precomposing with a Weierstrass isomorphism.
EXAMPLES:
sage: E = EllipticCurve(GF(31), [1,1,0,1,-1])
sage: R.<x> = GF(31)[]
sage: f = x^3 + 9*x^2 + x + 30
sage: phi = EllipticCurveIsogeny(E, f)
sage: Epr = E.short_weierstrass_model()
sage: isom = Epr.isomorphism_to(E)
sage: phi.set_pre_isomorphism(isom)
sage: phi.rational_maps()
((-6*x^4 - 3*x^3 + 12*x^2 + 10*x - 1)/(x^3 + x - 12), (3*x^7 + x^6*y - 14*x^6 - 3*x^5 + 5*x^4*y + 7*x^4 + 8*x^3*y - 8*x^3 - 5*x^2*y + 5*x^2 - 14*x*y + 14*x - 6*y - 6)/(x^6 + 2*x^4 + 7*x^3 + x^2 + 7*x - 11))
sage: phi(Epr((0,22)))
(13 : 21 : 1)
sage: phi(Epr((3,7)))
(14 : 17 : 1)
sage: E = EllipticCurve(GF(29), [0,0,0,1,0])
sage: R.<x> = GF(29)[]
sage: f = x^2 + 5
sage: phi = EllipticCurveIsogeny(E, f)
sage: phi
Isogeny of degree 5 from Elliptic Curve defined by y^2 = x^3 + x over Finite Field of size 29 to Elliptic Curve defined by y^2 = x^3 + 20*x over Finite Field of size 29
sage: from sage.schemes.elliptic_curves.weierstrass_morphism import WeierstrassIsomorphism
sage: inv_isom = WeierstrassIsomorphism(E, (1,-2,5,10))
sage: Epr = inv_isom.codomain().codomain()
sage: isom = Epr.isomorphism_to(E)
sage: phi.set_pre_isomorphism(isom); phi
Isogeny of degree 5 from Elliptic Curve defined by y^2 + 10*x*y + 20*y = x^3 + 27*x^2 + 6 over Finite Field of size 29 to Elliptic Curve defined by y^2 = x^3 + 20*x over Finite Field of size 29
sage: phi(Epr((12,1)))
(26 : 0 : 1)
sage: phi(Epr((2,9)))
(0 : 0 : 1)
sage: phi(Epr((21,12)))
(3 : 0 : 1)
sage: phi.rational_maps()[0]
(x^5 - 10*x^4 - 6*x^3 - 7*x^2 - x + 3)/(x^4 - 8*x^3 + 5*x^2 - 14*x - 6)
sage: E = EllipticCurve('11a1')
sage: R.<x> = QQ[]
sage: f = x^2 - 21*x + 80
sage: phi = EllipticCurveIsogeny(E, f); phi
Isogeny of degree 5 from Elliptic Curve defined by y^2 + y = x^3 - x^2 - 10*x - 20 over Rational Field to Elliptic Curve defined by y^2 + y = x^3 - x^2 - 7820*x - 263580 over Rational Field
sage: from sage.schemes.elliptic_curves.weierstrass_morphism import WeierstrassIsomorphism
sage: Epr = E.short_weierstrass_model()
sage: isom = Epr.isomorphism_to(E)
sage: phi.set_pre_isomorphism(isom)
sage: phi
Isogeny of degree 5 from Elliptic Curve defined by y^2 = x^3 - 13392*x - 1080432 over Rational Field to Elliptic Curve defined by y^2 + y = x^3 - x^2 - 7820*x - 263580 over Rational Field
sage: phi(Epr((168,1188)))
(0 : 1 : 0)
Compose this isogeny with (negation).
EXAMPLES:
sage: E = EllipticCurve(GF(23), [0,0,0,1,0])
sage: f = E.torsion_polynomial(3)/3
sage: phi = EllipticCurveIsogeny(E, f, E)
sage: phi.rational_maps() == E.multiplication_by_m(3)
False
sage: phi.switch_sign()
sage: phi.rational_maps() == E.multiplication_by_m(3)
True
sage: E = EllipticCurve(GF(17), [-2, 3, -5, 7, -11])
sage: R.<x> = GF(17)[]
sage: f = x+6
sage: phi = EllipticCurveIsogeny(E, f)
sage: phi
Isogeny of degree 2 from Elliptic Curve defined by y^2 + 15*x*y + 12*y = x^3 + 3*x^2 + 7*x + 6 over Finite Field of size 17 to Elliptic Curve defined by y^2 + 15*x*y + 12*y = x^3 + 3*x^2 + 4*x + 8 over Finite Field of size 17
sage: phi.rational_maps()
((x^2 + 6*x + 4)/(x + 6), (x^2*y - 5*x*y + 8*x - 2*y)/(x^2 - 5*x + 2))
sage: phi.switch_sign()
sage: phi
Isogeny of degree 2 from Elliptic Curve defined by y^2 + 15*x*y + 12*y = x^3 + 3*x^2 + 7*x + 6 over Finite Field of size 17 to Elliptic Curve defined by y^2 + 15*x*y + 12*y = x^3 + 3*x^2 + 4*x + 8 over Finite Field of size 17
sage: phi.rational_maps()
((x^2 + 6*x + 4)/(x + 6),
(2*x^3 - x^2*y - 5*x^2 + 5*x*y - 4*x + 2*y + 7)/(x^2 - 5*x + 2))
sage: E = EllipticCurve('11a1')
sage: R.<x> = QQ[]
sage: f = x^2 - 21*x + 80
sage: phi = EllipticCurveIsogeny(E, f)
sage: (xmap1, ymap1) = phi.rational_maps()
sage: phi.switch_sign()
sage: (xmap2, ymap2) = phi.rational_maps()
sage: xmap1 == xmap2
True
sage: ymap1 == -ymap2 - E.a1()*xmap2 - E.a3()
True
sage: K.<a> = NumberField(x^2 + 1)
sage: E = EllipticCurve(K, [0,0,0,1,0])
sage: R.<x> = K[]
sage: phi = EllipticCurveIsogeny(E, x-a)
sage: phi.rational_maps()
((x^2 + (-a)*x - 2)/(x + (-a)), (x^2*y + (-2*a)*x*y + y)/(x^2 + (-2*a)*x - 1))
sage: phi.switch_sign()
sage: phi.rational_maps()
((x^2 + (-a)*x - 2)/(x + (-a)), (-x^2*y + (2*a)*x*y - y)/(x^2 + (-2*a)*x - 1))
Return the rational map giving the -coordinate of this isogeny.
Note
This function returns the -coordinate component of the
isogeny as a rational function in
, where
is the
base field. To obtain both coordiunate functions as
elements of the function field
in two variables,
use rational_maps().
EXAMPLES:
sage: E = EllipticCurve(QQ, [0,2,0,1,-1])
sage: phi = EllipticCurveIsogeny(E, [1])
sage: phi.x_rational_map()
x
sage: E = EllipticCurve(GF(17), [0,0,0,3,0])
sage: phi = EllipticCurveIsogeny(E, E((0,0)))
sage: phi.x_rational_map()
(x^2 + 3)/x
Compute the codomain curve given parameters and
(as in
Velu / Kohel / etc formulas).
INPUT:
OUTPUT:
The elliptic curve with invariants
where
.
EXAMPLES:
This formula is used by every Isogeny instantiation:
sage: E = EllipticCurve(GF(19), [1,2,3,4,5])
sage: phi = EllipticCurveIsogeny(E, E((1,2)) )
sage: phi.codomain()
Elliptic Curve defined by y^2 + x*y + 3*y = x^3 + 2*x^2 + 9*x + 13 over Finite Field of size 19
sage: from sage.schemes.elliptic_curves.ell_curve_isogeny import compute_codomain_formula
sage: v = phi._EllipticCurveIsogeny__v
sage: w = phi._EllipticCurveIsogeny__w
sage: compute_codomain_formula(E, v, w) == phi.codomain()
True
Compute the codomain from the kernel polynomial using Kohel’s formulas.
INPUT:
OUTPUT:
(elliptic curve) – the codomain elliptic curve E/kernel
EXAMPLES:
sage: from sage.schemes.elliptic_curves.ell_curve_isogeny import compute_codomain_kohel
sage: E = EllipticCurve(GF(19), [1,2,3,4,5])
sage: phi = EllipticCurveIsogeny(E, [9,1])
sage: phi.codomain() == isogeny_codomain_from_kernel(E, [9,1])
True
sage: compute_codomain_kohel(E, [9,1], 2)
Elliptic Curve defined by y^2 + x*y + 3*y = x^3 + 2*x^2 + 9*x + 8 over Finite Field of size 19
sage: R.<x> = GF(19)[]
sage: E = EllipticCurve(GF(19), [18,17,16,15,14])
sage: phi = EllipticCurveIsogeny(E, x^3 + 14*x^2 + 3*x + 11)
sage: phi.codomain() == isogeny_codomain_from_kernel(E, x^3 + 14*x^2 + 3*x + 11)
True
sage: compute_codomain_kohel(E, x^3 + 14*x^2 + 3*x + 11, 7)
Elliptic Curve defined by y^2 + 18*x*y + 16*y = x^3 + 17*x^2 + 18*x + 18 over Finite Field of size 19
sage: E = EllipticCurve(GF(19), [1,2,3,4,5])
sage: phi = EllipticCurveIsogeny(E, x^3 + 7*x^2 + 15*x + 12)
sage: isogeny_codomain_from_kernel(E, x^3 + 7*x^2 + 15*x + 12) == phi.codomain()
True
sage: compute_codomain_kohel(E, x^3 + 7*x^2 + 15*x + 12,4)
Elliptic Curve defined by y^2 + x*y + 3*y = x^3 + 2*x^2 + 3*x + 15 over Finite Field of size 19
Note
This function uses the formulas of Section 2.4 of [K96].
REFERENCES:
[K96] | Kohel, “Endomorphism Rings of Elliptic Curves over Finite Fields”, UC Berkeley PhD thesis 1996. |
Return intermediate curves and isomorphisms.
Note
This is used so we can compute functions from the short
Weierstrass model more easily.
Warning
The base field must be of characteristic not equal to 2,3.
INPUT:
OUTPUT:
tuple (pre_isomorphism, post_isomorphism, intermediate_domain, intermediate_codomain):
EXAMPLES:
sage: from sage.schemes.elliptic_curves.ell_curve_isogeny import compute_intermediate_curves
sage: E = EllipticCurve(GF(83), [1,0,1,1,0])
sage: R.<x> = GF(83)[]; f = x+24
sage: phi = EllipticCurveIsogeny(E, f)
sage: E2 = phi.codomain()
sage: compute_intermediate_curves(E, E2)
(Elliptic Curve defined by y^2 = x^3 + 62*x + 74 over Finite Field of size 83,
Elliptic Curve defined by y^2 = x^3 + 65*x + 69 over Finite Field of size 83,
Generic morphism:
From: Abelian group of points on Elliptic Curve defined by y^2 + x*y + y = x^3 + x over Finite Field of size 83
To: Abelian group of points on Elliptic Curve defined by y^2 = x^3 + 62*x + 74 over Finite Field of size 83
Via: (u,r,s,t) = (1, 76, 41, 3),
Generic morphism:
From: Abelian group of points on Elliptic Curve defined by y^2 = x^3 + 65*x + 69 over Finite Field of size 83
To: Abelian group of points on Elliptic Curve defined by y^2 + x*y + y = x^3 + 4*x + 16 over Finite Field of size 83
Via: (u,r,s,t) = (1, 7, 42, 42))
sage: R.<x> = QQ[]
sage: K.<i> = NumberField(x^2 + 1)
sage: E = EllipticCurve(K, [0,0,0,1,0])
sage: E2 = EllipticCurve(K, [0,0,0,16,0])
sage: compute_intermediate_curves(E, E2)
(Elliptic Curve defined by y^2 = x^3 + x over Number Field in i with defining polynomial x^2 + 1,
Elliptic Curve defined by y^2 = x^3 + 16*x over Number Field in i with defining polynomial x^2 + 1,
Generic endomorphism of Abelian group of points on Elliptic Curve defined by y^2 = x^3 + x over Number Field in i with defining polynomial x^2 + 1
Via: (u,r,s,t) = (1, 0, 0, 0),
Generic endomorphism of Abelian group of points on Elliptic Curve defined by y^2 = x^3 + 16*x over Number Field in i with defining polynomial x^2 + 1
Via: (u,r,s,t) = (1, 0, 0, 0))
Return the kernel polynomial of an isogeny of degree ell between E1 and E2.
INPUT:
OUTPUT:
polynomial over the field of definition of E1, E2, that is the kernel polynomial of the isogeny from E1 to E2.
Note
If there is no degree ell, cyclic, separable, normalized isogeny from E1 to E2 then an error will be raised.
EXAMPLES:
sage: from sage.schemes.elliptic_curves.ell_curve_isogeny import compute_isogeny_kernel_polynomial
sage: E = EllipticCurve(GF(37), [0,0,0,1,8])
sage: R.<x> = GF(37)[]
sage: f = (x + 14) * (x + 30)
sage: phi = EllipticCurveIsogeny(E, f)
sage: E2 = phi.codomain()
sage: compute_isogeny_kernel_polynomial(E, E2, 5)
x^2 + 7*x + 13
sage: f
x^2 + 7*x + 13
sage: R.<x> = QQ[]
sage: K.<i> = NumberField(x^2 + 1)
sage: E = EllipticCurve(K, [0,0,0,1,0])
sage: E2 = EllipticCurve(K, [0,0,0,16,0])
sage: compute_isogeny_kernel_polynomial(E, E2, 4)
x^3 + x
Return the kernel polynomials of an isogeny of degree ell between E1 and E2.
INPUT:
OUTPUT:
polynomial over the field of definition of E1, E2, that is the kernel polynomial of the isogeny from E1 to E2.
Note
There must be a degree ell, separable, normalized cyclic isogeny from E1 to E2, or an error will be raised.
ALGORITHM:
This function uses Starks Algorithm as presented in section 6.2 of [BMSS].
Note
As published in [BMSS], the algorithm is incorrect, and a correct version (with slightly different notation) can be found in [M09]. The algorithm originates in [S72].
REFERENCES:
[BMSS] | (1, 2) Boston, Morain, Salvy, Schost, “Fast Algorithms for Isogenies.” |
[M09] | Moody, “The Diffie-Hellman Problem and Generalization of Verheul’s Theorem” |
[S72] | Stark, “Class-numbers of complex quadratic fields.” |
EXAMPLES:
sage: from sage.schemes.elliptic_curves.ell_curve_isogeny import compute_isogeny_starks, compute_sequence_of_maps
sage: E = EllipticCurve(GF(97), [1,0,1,1,0])
sage: R.<x> = GF(97)[]; f = x^5 + 27*x^4 + 61*x^3 + 58*x^2 + 28*x + 21
sage: phi = EllipticCurveIsogeny(E, f)
sage: E2 = phi.codomain()
sage: (isom1, isom2, E1pr, E2pr, ker_poly) = compute_sequence_of_maps(E, E2, 11)
sage: compute_isogeny_starks(E1pr, E2pr, 11)
x^10 + 37*x^9 + 53*x^8 + 66*x^7 + 66*x^6 + 17*x^5 + 57*x^4 + 6*x^3 + 89*x^2 + 53*x + 8
sage: E = EllipticCurve(GF(37), [0,0,0,1,8])
sage: R.<x> = GF(37)[]
sage: f = (x + 14) * (x + 30)
sage: phi = EllipticCurveIsogeny(E, f)
sage: E2 = phi.codomain()
sage: compute_isogeny_starks(E, E2, 5)
x^4 + 14*x^3 + x^2 + 34*x + 21
sage: f**2
x^4 + 14*x^3 + x^2 + 34*x + 21
sage: E = EllipticCurve(QQ, [0,0,0,1,0])
sage: R.<x> = QQ[]
sage: f = x
sage: phi = EllipticCurveIsogeny(E, f)
sage: E2 = phi.codomain()
sage: compute_isogeny_starks(E, E2, 2)
x
Return intermediate curves, isomorphisms and kernel polynomial.
INPUT:
OUTPUT:
(pre_isom, post_isom, E1pr, E2pr, ker_poly) where:
EXAMPLES:
sage: from sage.schemes.elliptic_curves.ell_curve_isogeny import compute_sequence_of_maps
sage: E = EllipticCurve('11a1')
sage: R.<x> = QQ[]; f = x^2 - 21*x + 80
sage: phi = EllipticCurveIsogeny(E, f)
sage: E2 = phi.codomain()
sage: compute_sequence_of_maps(E, E2, 5)
(Generic morphism:
From: Abelian group of points on Elliptic Curve defined by y^2 + y = x^3 - x^2 - 10*x - 20 over Rational Field
To: Abelian group of points on Elliptic Curve defined by y^2 = x^3 - 31/3*x - 2501/108 over Rational Field
Via: (u,r,s,t) = (1, 1/3, 0, -1/2),
Generic morphism:
From: Abelian group of points on Elliptic Curve defined by y^2 = x^3 - 23461/3*x - 28748141/108 over Rational Field
To: Abelian group of points on Elliptic Curve defined by y^2 + y = x^3 - x^2 - 7820*x - 263580 over Rational Field
Via: (u,r,s,t) = (1, -1/3, 0, 1/2),
Elliptic Curve defined by y^2 = x^3 - 31/3*x - 2501/108 over Rational Field,
Elliptic Curve defined by y^2 = x^3 - 23461/3*x - 28748141/108 over Rational Field,
x^2 - 61/3*x + 658/9)
sage: K.<i> = NumberField(x^2 + 1)
sage: E = EllipticCurve(K, [0,0,0,1,0])
sage: E2 = EllipticCurve(K, [0,0,0,16,0])
sage: compute_sequence_of_maps(E, E2, 4)
(Generic endomorphism of Abelian group of points on Elliptic Curve defined by y^2 = x^3 + x over Number Field in i with defining polynomial x^2 + 1
Via: (u,r,s,t) = (1, 0, 0, 0),
Generic endomorphism of Abelian group of points on Elliptic Curve defined by y^2 = x^3 + 16*x over Number Field in i with defining polynomial x^2 + 1
Via: (u,r,s,t) = (1, 0, 0, 0),
Elliptic Curve defined by y^2 = x^3 + x over Number Field in i with defining polynomial x^2 + 1,
Elliptic Curve defined by y^2 = x^3 + 16*x over Number Field in i with defining polynomial x^2 + 1,
x^3 + x)
sage: E = EllipticCurve(GF(97), [1,0,1,1,0])
sage: R.<x> = GF(97)[]; f = x^5 + 27*x^4 + 61*x^3 + 58*x^2 + 28*x + 21
sage: phi = EllipticCurveIsogeny(E, f)
sage: E2 = phi.codomain()
sage: compute_sequence_of_maps(E, E2, 11)
(Generic morphism:
From: Abelian group of points on Elliptic Curve defined by y^2 + x*y + y = x^3 + x over Finite Field of size 97
To: Abelian group of points on Elliptic Curve defined by y^2 = x^3 + 52*x + 31 over Finite Field of size 97
Via: (u,r,s,t) = (1, 8, 48, 44),
Generic morphism:
From: Abelian group of points on Elliptic Curve defined by y^2 = x^3 + 41*x + 66 over Finite Field of size 97
To: Abelian group of points on Elliptic Curve defined by y^2 + x*y + y = x^3 + 87*x + 26 over Finite Field of size 97
Via: (u,r,s,t) = (1, 89, 49, 49),
Elliptic Curve defined by y^2 = x^3 + 52*x + 31 over Finite Field of size 97,
Elliptic Curve defined by y^2 = x^3 + 41*x + 66 over Finite Field of size 97,
x^5 + 67*x^4 + 13*x^3 + 35*x^2 + 77*x + 69)
Compute Velu’s (v,w) using Kohel’s formulas for isogenies of degree exactly divisible by 2.
INPUT:
OUTPUT:
(tuple) Velu’s isogeny parameters (v,w).
EXAMPLES:
This function will be implicitly called by the following example:
sage: E = EllipticCurve(GF(19), [1,2,3,4,5])
sage: phi = EllipticCurveIsogeny(E, [9,1]); phi
Isogeny of degree 2 from Elliptic Curve defined by y^2 + x*y + 3*y = x^3 + 2*x^2 + 4*x + 5 over Finite Field of size 19 to Elliptic Curve defined by y^2 + x*y + 3*y = x^3 + 2*x^2 + 9*x + 8 over Finite Field of size 19
sage: from sage.schemes.elliptic_curves.ell_curve_isogeny import compute_vw_kohel_even_deg1
sage: a1,a2,a3,a4,a6 = E.ainvs()
sage: x0 = -9
sage: y0 = -(a1*x0 + a3)/2
sage: compute_vw_kohel_even_deg1(x0, y0, a1, a2, a4)
(18, 9)
Compute Velu’s (v,w) using Kohel’s formulas for isogenies of degree divisible by 4.
INPUT:
OUTPUT:
(tuple) Velu’s isogeny parameters (v,w).
EXAMPLES:
This function will be implicitly called by the following example:
sage: E = EllipticCurve(GF(19), [1,2,3,4,5])
sage: R.<x> = GF(19)[]
sage: phi = EllipticCurveIsogeny(E, x^3 + 7*x^2 + 15*x + 12); phi
Isogeny of degree 4 from Elliptic Curve defined by y^2 + x*y + 3*y = x^3 + 2*x^2 + 4*x + 5 over Finite Field of size 19 to Elliptic Curve defined by y^2 + x*y + 3*y = x^3 + 2*x^2 + 3*x + 15 over Finite Field of size 19
sage: from sage.schemes.elliptic_curves.ell_curve_isogeny import compute_vw_kohel_even_deg3
sage: (b2,b4) = (E.b2(), E.b4())
sage: (s1, s2, s3) = (-7, 15, -12)
sage: compute_vw_kohel_even_deg3(b2, b4, s1, s2, s3)
(4, 7)
Compute Velu’s (v,w) using Kohel’s formulas for isogenies of odd degree.
INPUT:
OUTPUT:
(tuple) Velu’s isogeny parameters (v,w).
EXAMPLES:
This function will be implicitly called by the following example:
sage: E = EllipticCurve(GF(19), [18,17,16,15,14])
sage: R.<x> = GF(19)[]
sage: phi = EllipticCurveIsogeny(E, x^3 + 14*x^2 + 3*x + 11); phi
Isogeny of degree 7 from Elliptic Curve defined by y^2 + 18*x*y + 16*y = x^3 + 17*x^2 + 15*x + 14 over Finite Field of size 19 to Elliptic Curve defined by y^2 + 18*x*y + 16*y = x^3 + 17*x^2 + 18*x + 18 over Finite Field of size 19
sage: from sage.schemes.elliptic_curves.ell_curve_isogeny import compute_vw_kohel_odd
sage: (b2,b4,b6) = (E.b2(), E.b4(), E.b6())
sage: (s1,s2,s3) = (-14,3,-11)
sage: compute_vw_kohel_odd(b2,b4,b6,s1,s2,s3,3)
(7, 1)
Returns a filled isogeny matrix giving all degrees from one giving only prime degrees.
INPUT:
OUTPUT:
(matrix) a square matrix with entries on the diagonal, and in
general the
,
entry is
if
is the minimal degree
of an isogeny from the
‘th to the
‘th curve,
EXAMPLES:
sage: M = Matrix([[0, 2, 3, 3, 0, 0], [2, 0, 0, 0, 3, 3], [3, 0, 0, 0, 2, 0], [3, 0, 0, 0, 0, 2], [0, 3, 2, 0, 0, 0], [0, 3, 0, 2, 0, 0]]); M
[0 2 3 3 0 0]
[2 0 0 0 3 3]
[3 0 0 0 2 0]
[3 0 0 0 0 2]
[0 3 2 0 0 0]
[0 3 0 2 0 0]
sage: from sage.schemes.elliptic_curves.ell_curve_isogeny import fill_isogeny_matrix
sage: fill_isogeny_matrix(M)
[ 1 2 3 3 6 6]
[ 2 1 6 6 3 3]
[ 3 6 1 9 2 18]
[ 3 6 9 1 18 2]
[ 6 3 2 18 1 9]
[ 6 3 18 2 9 1]
Compute the isogeny codomain given a kernel.
INPUT:
E - The domain elliptic curve.
kernel polynomial (specified as a either a univariate polynomial or a coefficient list.)
of the kernel.
OUTPUT:
(elliptic curve) the codomain of the separable normalized isogeny from this kernel
EXAMPLES:
sage: from sage.schemes.elliptic_curves.ell_curve_isogeny import isogeny_codomain_from_kernel
sage: E = EllipticCurve(GF(7), [1,0,1,0,1])
sage: R.<x> = GF(7)[]
sage: isogeny_codomain_from_kernel(E, [4,1], degree=3)
Elliptic Curve defined by y^2 + x*y + y = x^3 + 4*x + 6 over Finite Field of size 7
sage: EllipticCurveIsogeny(E, [4,1]).codomain() == isogeny_codomain_from_kernel(E, [4,1], degree=3)
True
sage: isogeny_codomain_from_kernel(E, x^3 + x^2 + 4*x + 3)
Elliptic Curve defined by y^2 + x*y + y = x^3 + 4*x + 6 over Finite Field of size 7
sage: isogeny_codomain_from_kernel(E, x^3 + 2*x^2 + 4*x + 3)
Elliptic Curve defined by y^2 + x*y + y = x^3 + 5*x + 2 over Finite Field of size 7
sage: E = EllipticCurve(GF(19), [1,2,3,4,5])
sage: kernel_list = [E((15,10)), E((10,3)),E((6,5))]
sage: isogeny_codomain_from_kernel(E, kernel_list)
Elliptic Curve defined by y^2 + x*y + 3*y = x^3 + 2*x^2 + 3*x + 15 over Finite Field of size 19
Helper function that allows the various isogeny functions to infer the algorithm type from the parameters passed in.
INPUT:
OUTPUT:
(string) either ‘velu’ or ‘kohel’
If kernel is a list of points on the EllipticCurve , then
we will try to use Velu’s algorithm.
If kernel is a list of coefficients or a univariate polynomial, we will try to use the Kohel’s algorithms.
EXAMPLES:
This helper function will be implicitly called by the following examples:
sage: R.<x> = GF(5)[]
sage: E = EllipticCurve(GF(5), [0,0,0,1,0])
We can construct the same isogeny from a kernel polynomial:
sage: phi = EllipticCurveIsogeny(E, x+3)
or from a list of coefficients of a kernel polynomial:
sage: phi == EllipticCurveIsogeny(E, [3,1])
True
or from a rational point which generates the kernel:
sage: phi == EllipticCurveIsogeny(E, E((2,0)) )
True
In the first two cases, Kohel’s algorithm will be used, while in the third case it is Velu:
sage: from sage.schemes.elliptic_curves.ell_curve_isogeny import isogeny_determine_algorithm
sage: isogeny_determine_algorithm(E, x+3)
'kohel'
sage: isogeny_determine_algorithm(E, [3, 1])
'kohel'
sage: isogeny_determine_algorithm(E, E((2,0)))
'velu'
Internal helper function for compute_isogeny_kernel_polynomial.
INPUT:
OUTPUT:
The maximum separable divisor of poly. If the input is a full
kernel polynomial where the roots which are -coordinates of
points of order greater than 2 have multiplicity 1, the output
will be a polynomial with the same roots, all of multiplicity 1.
EXAMPLES:
The following example implicitly exercises this function:
sage: E = EllipticCurve(GF(37), [0,0,0,1,8])
sage: R.<x> = GF(37)[]
sage: f = (x + 10) * (x + 12) * (x + 16)
sage: phi = EllipticCurveIsogeny(E, f)
sage: E2 = phi.codomain()
sage: from sage.schemes.elliptic_curves.ell_curve_isogeny import compute_isogeny_starks
sage: from sage.schemes.elliptic_curves.ell_curve_isogeny import split_kernel_polynomial
sage: ker_poly = compute_isogeny_starks(E, E2, 7); ker_poly
x^6 + 2*x^5 + 20*x^4 + 11*x^3 + 36*x^2 + 35*x + 16
sage: ker_poly.factor()
(x + 10)^2 * (x + 12)^2 * (x + 16)^2
sage: poly = split_kernel_polynomial(ker_poly); poly
x^3 + x^2 + 28*x + 33
sage: poly.factor()
(x + 10) * (x + 12) * (x + 16)
Returns the greatest common divisor of psi and the 2 torsion
polynomial of .
INPUT:
OUTPUT:
(polynomial) the gcd of psi and the 2-torsion polynomial of E.
EXAMPLES:
Every function that computes the kernel polynomial via Kohel’s formulas will call this function:
sage: E = EllipticCurve(GF(19), [1,2,3,4,5])
sage: R.<x> = GF(19)[]
sage: phi = EllipticCurveIsogeny(E, x + 13)
sage: isogeny_codomain_from_kernel(E, x + 13) == phi.codomain()
True
sage: from sage.schemes.elliptic_curves.ell_curve_isogeny import two_torsion_part
sage: two_torsion_part(E, x+13)
x + 13
Reverses the action of fill_isogeny_matrix.
INPUT:
OUTPUT:
(matrix) a square symmetric matrix obtained from M by
replacing non-prime entries with .
EXAMPLES:
sage: M = Matrix([[0, 2, 3, 3, 0, 0], [2, 0, 0, 0, 3, 3], [3, 0, 0, 0, 2, 0], [3, 0, 0, 0, 0, 2], [0, 3, 2, 0, 0, 0], [0, 3, 0, 2, 0, 0]]); M
[0 2 3 3 0 0]
[2 0 0 0 3 3]
[3 0 0 0 2 0]
[3 0 0 0 0 2]
[0 3 2 0 0 0]
[0 3 0 2 0 0]
sage: from sage.schemes.elliptic_curves.ell_curve_isogeny import fill_isogeny_matrix, unfill_isogeny_matrix
sage: M1 = fill_isogeny_matrix(M); M1
[ 1 2 3 3 6 6]
[ 2 1 6 6 3 3]
[ 3 6 1 9 2 18]
[ 3 6 9 1 18 2]
[ 6 3 2 18 1 9]
[ 6 3 18 2 9 1]
sage: unfill_isogeny_matrix(M1)
[0 2 3 3 0 0]
[2 0 0 0 3 3]
[3 0 0 0 2 0]
[3 0 0 0 0 2]
[0 3 2 0 0 0]
[0 3 0 2 0 0]
sage: unfill_isogeny_matrix(M1) == M
True