AUTHOR:
-Daniel Shumow <shumow@gmail.com>: 2009-04-19: initial version
Class Implementing Isogenies of Elliptic Curves
This class implements separable, normalized isogenies of elliptic curves.
Several different algorithms for computing isogenies are available. These include:
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)), E(0)], algorithm="velu"); 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))
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, 2*P, 3*P], algorithm="velu"); 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.univariate_polynomial().list()
sage: phi_k = EllipticCurveIsogeny(E, ker_poly_list, algorithm="kohel")
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()
True
sage: phi_k.degree()
3
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, algorithm="kohel"); 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)], algorithm="velu")
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], algorithm="kohel")
sage: 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, algorithm="velu")
sage: 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.univariate_polynomial().list()
sage: phi_k = EllipticCurveIsogeny(E, ker_poly_list, algorithm="kohel"); 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 (with odd degree isogeny):
sage: E = EllipticCurve('11a1')
sage: P_list = E.torsion_points()
sage: phi_v = EllipticCurveIsogeny(E, P_list, algorithm="velu"); 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.univariate_polynomial().list()
sage: phi_k = EllipticCurveIsogeny(E, ker_poly_list, algorithm="kohel"); 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 - x^2 - 10*x - 79/4)
sage: EK = E.change_ring(K)
sage: phi_v = EllipticCurveIsogeny(EK, P_list, algorithm="velu"); 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 - x^2 - 10*x - 79/4 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 - x^2 - 10*x - 79/4
sage: P = EK((alpha,-1/2))
sage: phi_v(P)
(488/121*alpha^2 + 1633/121*alpha - 3920/121 : -1/2 : 1)
sage: ker_poly = phi_v.kernel_polynomial().univariate_polynomial()
sage: ker_poly
x^2 - 21*x + 80
sage: ker_poly_list = ker_poly.list()
sage: phi_k = EllipticCurveIsogeny(EK, ker_poly_list, algorithm="kohel")
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 - x^2 - 10*x - 79/4 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 - x^2 - 10*x - 79/4
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
A private function to clear the hash if the codomain has been modified by a pre or post isomorphism.
EXAMPLES:
sage: F = GF(7);
sage: E = EllipticCurve(j=F(0))
sage: phi = EllipticCurveIsogeny(E, [E(0), E((0,-1)), E((0,1))])
sage: old_hash = hash(phi)
sage: from sage.schemes.elliptic_curves.weierstrass_morphism import WeierstrassIsomorphism
sage: phi.set_post_isomorphism(WeierstrassIsomorphism(phi.codomain(), (-1,2,-3,4)))
sage: hash(phi) == old_hash
False
sage: R.<x> = QQ[]
sage: E = EllipticCurve(QQ, [0,0,0,1,0])
sage: phi = EllipticCurveIsogeny(E, x)
sage: old_ratl_maps = phi.rational_maps()
sage: from sage.schemes.elliptic_curves.weierstrass_morphism import WeierstrassIsomorphism
sage: phi.set_post_isomorphism(WeierstrassIsomorphism(phi.codomain(), (-1,0,0,0)))
sage: old_ratl_maps == phi.rational_maps()
False
sage: old_ratl_maps[1] == -phi.rational_maps()[1]
True
sage: F = GF(127); R.<x> = F[]
sage: E = EllipticCurve(j=F(1728))
sage: f = x^5 + 43*x^4 + 97*x^3 + 81*x^2 + 42*x + 82
sage: phi = EllipticCurveIsogeny(E, f)
sage: old_hash = hash(phi)
sage: old_ratl_maps = phi.rational_maps()
sage: from sage.schemes.elliptic_curves.weierstrass_morphism import WeierstrassIsomorphism
sage: phi.set_post_isomorphism(WeierstrassIsomorphism(phi.codomain(), (-13,13,-13,13)))
sage: old_hash == hash(phi)
False
sage: old_ratl_maps == phi.rational_maps()
False
sage: phi._EllipticCurveIsogeny__clear_cached_values()
Private function that computes and sets the isogeny codomain.
EXAMPLES:
These examples inherently exercise this function:
sage: E = EllipticCurve(j=GF(7)(1728))
sage: phi = EllipticCurveIsogeny(E, [E(0), E((0,0))])
sage: phi.codomain()
Elliptic Curve defined by y^2 = x^3 + 3*x over Finite Field of size 7
sage: phi._EllipticCurveIsogeny__compute_E2()
sage: R.<x> = GF(7)[]
sage: phi = EllipticCurveIsogeny(E, x)
sage: phi.codomain()
Elliptic Curve defined by y^2 = x^3 + 3*x over Finite Field of size 7
sage: phi._EllipticCurveIsogeny__compute_E2()
Private function that computes and initializes the codomain of the isogeny (via Kohel’s.)
EXAMPLES:
These examples inherently exercise this private function:
sage: R.<x> = GF(7)[]
sage: E = EllipticCurve(GF(7), [0,-1,0,0,1])
sage: phi = EllipticCurveIsogeny(E, x+6, degree=3)
sage: phi.codomain()
Elliptic Curve defined by y^2 = x^3 + 6*x^2 + 4*x + 2 over Finite Field of size 7
sage: phi._EllipticCurveIsogeny__compute_E2_via_kohel()
Elliptic Curve defined by y^2 = x^3 + 6*x^2 + 4*x + 2 over Finite Field of size 7
Private function that computes the codomain via Velu’s formulas.
EXAMPLES:
The following example inherently exercises this function:
sage: E = EllipticCurve(GF(7), [0,0,0,-1,0])
sage: P = E((4,2))
sage: phi = EllipticCurveIsogeny(E, [P, 2*P, 3*P, 4*P]);
sage: phi.codomain()
Elliptic Curve defined by y^2 = x^3 + 2*x over Finite Field of size 7
sage: phi._EllipticCurveIsogeny__compute_E2_via_velu()
Elliptic Curve defined by y^2 = x^3 + 2*x over Finite Field of size 7
Private function that initializes the omega polynomial (from Kohel’s formulas) in the case that the characteristic of the underlying field is not 2.
EXAMPLES:
These examples inherently exercise this private function:
sage: R.<x> = GF(7)[]
sage: E = EllipticCurve(GF(7), [0,-1,0,0,1])
sage: phi = EllipticCurveIsogeny(E, x+6, degree=3); phi
Isogeny of degree 3 from Elliptic Curve defined by y^2 = x^3 + 6*x^2 + 1 over Finite Field of size 7 to Elliptic Curve defined by y^2 = x^3 + 6*x^2 + 4*x + 2 over Finite Field of size 7
sage: R.<x,y> = GF(7)[]
sage: psi = phi._EllipticCurveIsogeny__psi
sage: psi_pr = psi.derivative(x)
sage: fi = phi._EllipticCurveIsogeny__phi
sage: fi_pr = fi.derivative(x)
sage: phi._EllipticCurveIsogeny__compute_omega_fast(E, psi, psi_pr, fi, fi_pr)
x^3*y - 3*x^2*y + x*y
Private function that initializes the omega polynomial (from Kohel’s formulas) in the case of general characteristic of the underlying field.
EXAMPLES:
These examples inherently exercise this private function:
sage: F = GF(2^4, 'alpha'); R.<x> = F[]
sage: alpha = F.gen()
sage: E = EllipticCurve(F, [1,1,F.gen(),F.gen()^2+1,1])
sage: f = x + alpha^2 + 1
sage: phi = EllipticCurveIsogeny(E, f); phi
Isogeny of degree 3 from Elliptic Curve defined by y^2 + x*y + alpha*y = x^3 + x^2 + (alpha^2+1)*x + 1 over Finite Field in alpha of size 2^4 to Elliptic Curve defined by y^2 + x*y + alpha*y = x^3 + x^2 + alpha*x + alpha^3 over Finite Field in alpha of size 2^4
sage: R.<x,y> = F[]
sage: psi = phi._EllipticCurveIsogeny__psi
sage: psi_pr = psi.derivative(x)
sage: fi = phi._EllipticCurveIsogeny__phi
sage: fi_pr = fi.derivative(x)
sage: phi._EllipticCurveIsogeny__compute_omega_general(E, psi, psi_pr, fi, fi_pr)
x^3*y + (alpha^2 + 1)*x^2*y + (alpha^2 + alpha + 1)*x^2 + (alpha^2 + 1)*x*y + (alpha^2 + alpha)*x + (alpha)*y + (alpha)
Private function that applies Kohel’s formulas.
EXAMPLES:
These examples inherently exercise this private function:
sage: R.<x> = GF(7)[]
sage: E = EllipticCurve(GF(7), [0,-1,0,0,1])
sage: phi = EllipticCurveIsogeny(E, x+6, degree=3)
sage: P = E((0,1)); phi(P)
(2 : 0 : 1)
sage: phi.rational_maps()
((x^3 - 2*x^2 + 3*x + 2)/(x^2 - 2*x + 1),
(x^3*y - 3*x^2*y + x*y)/(x^3 - 3*x^2 + 3*x - 1))
sage: phi._EllipticCurveIsogeny__compute_via_kohel(0,1)
(2, 0)
sage: R.<x,y> = GF(7)[]
sage: phi._EllipticCurveIsogeny__compute_via_kohel(x,y)
((x^3 - 2*x^2 + 3*x + 2)/(x^2 - 2*x + 1),
(x^3*y - 3*x^2*y + x*y)/(x^3 - 3*x^2 + 3*x - 1))
Private function that computes a numeric result of this isogeny (via Kohel’s formulas.)
EXAMPLES:
These examples inherently exercise this private function:
sage: R.<x> = GF(7)[]
sage: E = EllipticCurve(GF(7), [0,-1,0,0,1])
sage: phi = EllipticCurveIsogeny(E, x+6, degree=3)
sage: P = E((0,1)); phi(P)
(2 : 0 : 1)
sage: P = E((1,1)); phi(P)
(0 : 1 : 0)
sage: phi._EllipticCurveIsogeny__compute_via_kohel_numeric(0, 1)
(2, 0)
sage: phi._EllipticCurveIsogeny__compute_via_kohel_numeric(1, 1)
(0 : 1 : 0)
Private function for Velu’s formulas, to perform the summation.
EXAMPLES:
The following example inherently exercises this function:
sage: E = EllipticCurve(GF(7), [0,0,0,-1,0])
sage: P = E((4,2))
sage: phi = EllipticCurveIsogeny(E, [P, 2*P, 3*P, 4*P]);
sage: Q = E((0,0)); phi(Q)
(0 : 0 : 1)
sage: phi.rational_maps()
((x^4 - 2*x^3 + x^2 - 3*x)/(x^3 - 2*x^2 + 3*x - 2),
(x^5*y - 2*x^3*y - x^2*y - 2*x*y + 2*y)/(x^5 + 3*x^3 + 3*x^2 + x - 1))
sage: phi._EllipticCurveIsogeny__compute_via_velu(0, 0)
(0, 0)
sage: R.<x,y> = GF(7)[]
sage: phi._EllipticCurveIsogeny__compute_via_velu(x, y)
((x^4 - 2*x^3 + x^2 - 3*x)/(x^3 - 2*x^2 + 3*x - 2),
(x^5*y - 2*x^3*y - x^2*y - 2*x*y + 2*y)/(x^5 + 3*x^3 + 3*x^2 + x - 1))
Private function that sorts the list of points in the kernel (For Velu’s formulas.) sorts out the 2 torsion points, and puts them in a diction
EXAMPLES:
The following example inherently exercises this function:
sage: E = EllipticCurve(GF(7), [0,0,0,-1,0])
sage: P = E((4,2))
sage: phi = EllipticCurveIsogeny(E, [P, 2*P, 3*P, 4*P]);
sage: Q = E((0,0)); phi(Q)
(0 : 0 : 1)
sage: Q = E((-1,0)); phi(Q)
(0 : 0 : 1)
sage: phi._EllipticCurveIsogeny__compute_via_velu_numeric(0, 0)
(0, 0)
sage: phi._EllipticCurveIsogeny__compute_via_velu_numeric(-1, 0)
(0, 0)
An internal function for EllipticCurveIsogeny objects that sets up the member variables necessary for algebra.
EXAMPLES:
The following tests inherently exercise this function:
sage: E = EllipticCurve(j=GF(17)(0))
sage: phi = EllipticCurveIsogeny(E, [E(0), E((-1,0))])
sage: phi._EllipticCurveIsogeny__init_algebraic_structs(E)
sage: E = EllipticCurve(QQ, [0,0,0,1,0])
sage: phi = EllipticCurveIsogeny(E, [E(0), E((0,0))])
sage: phi._EllipticCurveIsogeny__init_algebraic_structs(E)
sage: F = GF(19); R.<x> = F[]
sage: E = EllipticCurve(j=GF(19)(0))
sage: phi = EllipticCurveIsogeny(E, x)
sage: phi._EllipticCurveIsogeny__init_algebraic_structs(E)
Private function that initializes the isogeny from a kernel polynomial, for Kohel’s algorithm in the even degree case.
EXAMPLES:
These examples inherently exercise this private function:
sage: R.<x> = GF(7)[]
sage: E = EllipticCurve(GF(7), [0,0,0,-1,0])
sage: phi = EllipticCurveIsogeny(E, x);phi
Isogeny of degree 2 from Elliptic Curve defined by y^2 = x^3 + 6*x over Finite Field of size 7 to Elliptic Curve defined by y^2 = x^3 + 4*x over Finite Field of size 7
sage: phi._EllipticCurveIsogeny__init_even_kernel_polynomial(E, x)
(x^3 - x, x^3*y + x*y, 6, 0, 1, 2)
sage: F = GF(2^4, 'alpha'); R.<x> = F[]
sage: E = EllipticCurve(F, [1,1,0,1,0])
sage: phi = EllipticCurveIsogeny(E, x); phi
Isogeny of degree 2 from Elliptic Curve defined by y^2 + x*y = x^3 + x^2 + x over Finite Field in alpha of size 2^4 to Elliptic Curve defined by y^2 + x*y = x^3 + x^2 + 1 over Finite Field in alpha of size 2^4
sage: phi._EllipticCurveIsogeny__init_even_kernel_polynomial(E, x)
(x^3 + x, x^3*y + x^2 + x*y, 1, 0, 1, 2)
sage: E = EllipticCurve(GF(7), [0,-1,0,0,1])
sage: R.<x> = GF(7)[]
sage: f = x^3 + 6*x^2 + 1
sage: phi = EllipticCurveIsogeny(E, f); phi
Isogeny of degree 4 from Elliptic Curve defined by y^2 = x^3 + 6*x^2 + 1 over Finite Field of size 7 to Elliptic Curve defined by y^2 = x^3 + 6*x^2 + 2*x + 5 over Finite Field of size 7
sage: P_list = E.points()
sage: Q = E(0)
sage: Q = phi.codomain()(0)
sage: for P in P_list: Q += phi(P)
sage: Q
(0 : 1 : 0)
sage: R.<x,y> = GF(7)[]
sage: f = 1 + x + 2*x^2 + 3*x^3
sage: phi._EllipticCurveIsogeny__init_even_kernel_polynomial(E, f)
(-3*x^7 + 3*x^6 - 2*x^5 - 2*x^4 + 3*x^3 - x^2 + 2,
-2*x^9*y + 3*x^8*y + 3*x^7*y - x^6*y - x^5*y - x^4*y - x^3*y + x^2*y + 3*x*y + 3*y,
3,
2,
3,
4)
Private function that initializes the isogeny from a list of points in the kernel (For Velu’s formulas.)
EXAMPLES:
The following example inherently exercises this function:
sage: E = EllipticCurve(GF(7), [0,0,0,-1,0])
sage: phi = EllipticCurveIsogeny(E, [E(0), E((0,0))]); phi
Isogeny of degree 2 from Elliptic Curve defined by y^2 = x^3 + 6*x over Finite Field of size 7 to Elliptic Curve defined by y^2 = x^3 + 4*x over Finite Field of size 7
sage: phi._EllipticCurveIsogeny__init_from_kernel_list([E(0), E((0,0))])
Private function that initializes the isogeny from a kernel polynomial.
EXAMPLES:
These examples inherently exercise this private function:
sage: R.<x> = GF(7)[]
sage: E = EllipticCurve(GF(7), [0,0,0,-1,0])
sage: phi = EllipticCurveIsogeny(E, x);phi
Isogeny of degree 2 from Elliptic Curve defined by y^2 = x^3 + 6*x over Finite Field of size 7 to Elliptic Curve defined by y^2 = x^3 + 4*x over Finite Field of size 7
sage: phi._EllipticCurveIsogeny__init_from_kernel_polynomial(x)
sage: E = EllipticCurve(GF(7), [0,-1,0,0,1])
sage: phi = EllipticCurveIsogeny(E, x+6, degree=3); phi
Isogeny of degree 3 from Elliptic Curve defined by y^2 = x^3 + 6*x^2 + 1 over Finite Field of size 7 to Elliptic Curve defined by y^2 = x^3 + 6*x^2 + 4*x + 2 over Finite Field of size 7
sage: phi._EllipticCurveIsogeny__init_from_kernel_polynomial(x+6, degree=3)
Private function that initializes the kernel polynomial (if the algorithm does not take it as a parameter.)
EXAMPLES:
The following examples inherently exercise this function:
sage: E = EllipticCurve(j=GF(7)(1728))
sage: phi = EllipticCurveIsogeny(E, [E(0), E((0,0))])
sage: phi.kernel_polynomial()
x
sage: phi._EllipticCurveIsogeny__init_kernel_polynomial()
[0, 1]
Private function for Velu’s formulas, helper function to initialize the rational maps.
EXAMPLES:
The following example inherently exercises this function:
sage: E = EllipticCurve(GF(7), [0,0,0,-1,0])
sage: P = E((4,2))
sage: phi = EllipticCurveIsogeny(E, [P, 2*P, 3*P, 4*P]);
sage: phi.kernel_polynomial()
x^2 + 2*x - 3
sage: phi._EllipticCurveIsogeny__init_kernel_polynomial_velu()
[4, 2, 1]
Private function that initializes the isogeny from a kernel polynomial.
EXAMPLES:
These examples inherently exercise this private function:
sage: R.<x> = GF(7)[]
sage: E = EllipticCurve(GF(7), [0,-1,0,0,1])
sage: phi = EllipticCurveIsogeny(E, x+6, degree=3); phi
Isogeny of degree 3 from Elliptic Curve defined by y^2 = x^3 + 6*x^2 + 1 over Finite Field of size 7 to Elliptic Curve defined by y^2 = x^3 + 6*x^2 + 4*x + 2 over Finite Field of size 7
sage: R.<x,y> = GF(7)[]
sage: phi._EllipticCurveIsogeny__init_odd_kernel_polynomial(E, x+6)
(x^3 - 2*x^2 + 3*x + 2, x^3*y - 3*x^2*y + x*y, 2, 6, 1, 3)
sage: F = GF(2^4, 'alpha'); R.<x> = F[]
sage: alpha = F.gen()
sage: E = EllipticCurve(F, [1,1,F.gen(),F.gen()^2+1,1])
sage: f = x + alpha^2 + 1
sage: phi = EllipticCurveIsogeny(E, f); phi
Isogeny of degree 3 from Elliptic Curve defined by y^2 + x*y + alpha*y = x^3 + x^2 + (alpha^2+1)*x + 1 over Finite Field in alpha of size 2^4 to Elliptic Curve defined by y^2 + x*y + alpha*y = x^3 + x^2 + alpha*x + alpha^3 over Finite Field in alpha of size 2^4
sage: R.<x,y> = F[]
sage: f = x + alpha^2 + 1
sage: phi._EllipticCurveIsogeny__init_odd_kernel_polynomial(E, f)
(x^3 + (alpha^2 + 1)*x + (alpha^3 + alpha^2 + alpha),
x^3*y + (alpha^2 + 1)*x^2*y + (alpha^2 + alpha + 1)*x^2 + (alpha^2 + 1)*x*y + (alpha^2 + alpha)*x + (alpha)*y + (alpha),
alpha^2 + alpha + 1,
alpha^3 + alpha^2 + alpha,
1,
3)
Private function that computes and initializes the rational maps.
EXAMPLES:
The following examples inherently exercise this function:
sage: E = EllipticCurve(j=GF(7)(1728))
sage: phi = EllipticCurveIsogeny(E, [E(0), E((0,0))])
sage: phi._EllipticCurveIsogeny__initialize_rational_maps()
sage: phi.rational_maps()
((x^2 + 1)/x, (x^2*y - y)/x^2)
sage: R.<x> = GF(7)[]
sage: phi = EllipticCurveIsogeny(E, x)
sage: phi = EllipticCurveIsogeny(E, x)
sage: phi.rational_maps()
((x^2 + 1)/x, (x^2*y - y)/x^2)
sage: phi._EllipticCurveIsogeny__initialize_rational_maps()
Private function that computes and initializes the rational maps of this isogeny.
EXAMPLES:
These examples inherently exercise this private function:
sage: R.<x> = GF(7)[]
sage: E = EllipticCurve(GF(7), [0,-1,0,0,1])
sage: phi = EllipticCurveIsogeny(E, x+6, degree=3)
sage: phi.rational_maps()
((x^3 - 2*x^2 + 3*x + 2)/(x^2 - 2*x + 1),
(x^3*y - 3*x^2*y + x*y)/(x^3 - 3*x^2 + 3*x - 1))
sage: phi._EllipticCurveIsogeny__initialize_rational_maps_via_kohel()
((x^3 - 2*x^2 + 3*x + 2)/(x^2 - 2*x + 1),
(x^3*y - 3*x^2*y + x*y)/(x^3 - 3*x^2 + 3*x - 1))
Private function for Velu’s formulas, helper function to initialize the rational maps.
EXAMPLES:
The following example inherently exercises this function:
sage: E = EllipticCurve(GF(7), [0,0,0,-1,0])
sage: P = E((4,2))
sage: phi = EllipticCurveIsogeny(E, [P, 2*P, 3*P, 4*P]);
sage: phi.rational_maps()
((x^4 - 2*x^3 + x^2 - 3*x)/(x^3 - 2*x^2 + 3*x - 2),
(x^5*y - 2*x^3*y - x^2*y - 2*x*y + 2*y)/(x^5 + 3*x^3 + 3*x^2 + x - 1))
sage: phi._EllipticCurveIsogeny__initialize_rational_maps_via_velu()
((x^4 - 2*x^3 + x^2 - 3*x)/(x^3 - 2*x^2 + 3*x - 2),
(x^5*y - 2*x^3*y - x^2*y - 2*x*y + 2*y)/(x^5 + 3*x^3 + 3*x^2 + x - 1))
Internal helper function, sets values on the super classes of this class.
EXAMPLES:
The following examples will implicitly exercise this function:
sage: E = EllipticCurve(GF(43), [2,3,5,7,11])
sage: R.<x> = GF(43)[]; f = x + 42
sage: phi = EllipticCurveIsogeny(E, f)
sage: phi._EllipticCurveIsogeny__perform_inheritance_housekeeping()
sage: from sage.schemes.elliptic_curves.weierstrass_morphism import WeierstrassIsomorphism
sage: E2 = phi.codomain()
sage: post_isom = WeierstrassIsomorphism(E2, (41, 37, 31, 29))
sage: phi.set_post_isomorphism(post_isom)
sage: E1pr = WeierstrassIsomorphism(E, (-1, 2, -3, 4)).codomain().codomain()
sage: pre_isom = E1pr.isomorphism_to(E)
sage: phi.set_pre_isomorphism(pre_isom)
Private function to set the post isomorphism and codomain (and keep track of the codomain of the isogeny.)
EXAMPLES:
The following examples inherently exercise this function:
sage: E = EllipticCurve(j=GF(7)(1728))
sage: phi = EllipticCurveIsogeny(E, [E(0), E((0,0))])
sage: from sage.schemes.elliptic_curves.weierstrass_morphism import WeierstrassIsomorphism
sage: E2 = phi.codomain()
sage: isom = WeierstrassIsomorphism(E2, (-1,2,-3,4))
sage: phi.set_post_isomorphism(isom)
sage: phi._EllipticCurveIsogeny__set_post_isomorphism(E2, WeierstrassIsomorphism(phi.codomain(), (1,-2,3,-4)))
sage: E2 == phi.codomain()
True
Private function to set the pre isomorphism and domain (and keep track of the domain of the isogeny.)
EXAMPLES:
sage: E = EllipticCurve(GF(43), [2,3,5,7,11])
sage: R.<x> = GF(43)[]; f = x + 42
sage: phi = EllipticCurveIsogeny(E, f)
sage: phi._EllipticCurveIsogeny__perform_inheritance_housekeeping()
sage: from sage.schemes.elliptic_curves.weierstrass_morphism import WeierstrassIsomorphism
sage: E1pr = WeierstrassIsomorphism(E, (-1, 2, -3, 4)).codomain().codomain()
sage: pre_isom = E1pr.isomorphism_to(E)
sage: phi.set_pre_isomorphism(pre_isom)
sage: phi._EllipticCurveIsogeny__set_pre_isomorphism(E, WeierstrassIsomorphism(E, (-1, 3, -3, 4)))
sage: E == phi.domain()
True
Private function to set up the post isomorphism given the codomain.
EXAMPLES:
The following examples inherently exercise this function:
sage: E = EllipticCurve(j=GF(7)(1728))
sage: E2 = EllipticCurve(GF(7), [0,0,0,5,0])
sage: phi = EllipticCurveIsogeny(E, [E(0), E((0,0))], E2); 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 + 5*x over Finite Field of size 7
sage: E3 = EllipticCurve(GF(7), [0,0,0,6,0])
sage: phi._EllipticCurveIsogeny__setup_post_isomorphism(E3, None)
sage: 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 + 6*x over Finite Field of size 7
sage: R.<x> = QQ[]
sage: E = EllipticCurve(j=1728)
sage: f = x^3 - x
sage: phi = EllipticCurveIsogeny(E, f, model='minimal'); phi
Isogeny of degree 4 from Elliptic Curve defined by y^2 = x^3 - x over Rational Field to Elliptic Curve defined by y^2 = x^3 - x over Rational Field
sage: phi = EllipticCurveIsogeny(E, f, model=None)
sage: phi._EllipticCurveIsogeny__setup_post_isomorphism(None, 'minimal')
sage: phi
Isogeny of degree 4 from Elliptic Curve defined by y^2 = x^3 - x over Rational Field to Elliptic Curve defined by y^2 = x^3 - x over Rational Field
Private function that sorts the list of points in the kernel (For Velu’s formulas.) sorts out the 2 torsion points, and puts them in a dictionary.
EXAMPLES:
The following example inherently exercises this function:
sage: E = EllipticCurve(GF(7), [0,0,0,-1,0])
sage: P = E((4,2))
sage: phi = EllipticCurveIsogeny(E, [P, 2*P, 3*P, 4*P]); phi
Isogeny of degree 4 from Elliptic Curve defined by y^2 = x^3 + 6*x over Finite Field of size 7 to Elliptic Curve defined by y^2 = x^3 + 2*x over Finite Field of size 7
sage: phi._EllipticCurveIsogeny__kernel_2tor = {}
sage: phi._EllipticCurveIsogeny__kernel_non2tor = {}
sage: phi._EllipticCurveIsogeny__sort_kernel_list()
Private function for Velu’s formulas, helper function to help perform the summation.
EXAMPLES:
The following example inherently exercises this function:
sage: E = EllipticCurve(GF(7), [0,0,0,-1,0])
sage: P = E((4,2))
sage: phi = EllipticCurveIsogeny(E, [P, 2*P, 3*P, 4*P]);
sage: Q = E((0,0)); phi(Q)
(0 : 0 : 1)
sage: phi.rational_maps()
((x^4 - 2*x^3 + x^2 - 3*x)/(x^3 - 2*x^2 + 3*x - 2),
(x^5*y - 2*x^3*y - x^2*y - 2*x*y + 2*y)/(x^5 + 3*x^3 + 3*x^2 + x - 1))
sage: E = EllipticCurve(GF(7), [0,0,0,1,0])
sage: phi = EllipticCurveIsogeny(E, [E(0), E((0,0))])
sage: Qvals = phi._EllipticCurveIsogeny__kernel_2tor[0]
sage: phi._EllipticCurveIsogeny__velu_sum_helper(Qvals, 0, 0, 5, 5)
(3, 3)
sage: R.<x,y> = GF(7)[]
sage: phi._EllipticCurveIsogeny__velu_sum_helper(Qvals, 0, 0, x, y)
(1/x, y/x^2)
Function that implements the call-ability of elliptic curve isogenies.
EXAMPLES:
sage: E = EllipticCurve(GF(17), [1, 9, 5, 4, 3])
sage: phi = EllipticCurveIsogeny(E, [6,13,1], algorithm="kohel")
sage: phi(E((1,0)))
(15 : 13 : 1)
sage: E = EllipticCurve(GF(23), [0,0,0,1,0])
sage: phi = EllipticCurveIsogeny(E, [E((0,0)), E(0)], algorithm="velu")
sage: phi(E((1,5)))
(2 : 0 : 1)
sage: phi(E(15,20), output_base_ring=GF(23^2,'alpha'))
(12 : 1 : 1)
sage: E = EllipticCurve(QQ, [0,0,0,3,0])
sage: P = E((1,2))
sage: phi = EllipticCurveIsogeny(E, [0,1], algorithm="kohel")
sage: phi(P)
(4 : -4 : 1)
sage: phi(-P)
(4 : 4 : 1)
Function that implements comparisons between isogeny objects.
This function works by comparing the underlying kernel objects.
EXAMPLES:
sage: E = EllipticCurve(QQ, [0,0,0,1,0])
sage: phi_v = EllipticCurveIsogeny(E, [E(0), E((0,0))], algorithm="velu")
sage: phi_k = EllipticCurveIsogeny(E, [0,1], algorithm="kohel")
sage: phi_k == phi_v
True
sage: E_F17 = EllipticCurve(GF(17), [0,0,0,1,0])
sage: phi_p = EllipticCurveIsogeny(E_F17, [0,1], algorithm="kohel")
sage: phi_p == phi_v
False
Function that implements the hashability of Isogeny objects.
This hashes the underlying kernel polynomial so that equal isogeny objects have the same hash value. Also, this hashes the base field, and domain and codomain curves as well, so that isogenies with the same kernel polynomial (over different base fields / curves) hash to different values.
EXAMPLES:
sage: E = EllipticCurve(QQ, [0,0,0,1,0])
sage: phi_v = EllipticCurveIsogeny(E, [E(0), E((0,0))], algorithm="velu")
sage: phi_k = EllipticCurveIsogeny(E, [0,1], algorithm="kohel")
sage: phi_k.__hash__() == phi_v.__hash__()
True
sage: E_F17 = EllipticCurve(GF(17), [0,0,0,1,1])
sage: phi_p = EllipticCurveIsogeny(E_F17, [0,1], algorithm="kohel")
sage: phi_p.__hash__() == phi_v.__hash__()
False
Constructor for EllipticCurveIsogeny class.
INPUT:
E - an elliptic curve, the domain of the isogeny to initialize.
If initiating from a domain/codomain, this must be set to None.
then this must be the codomain of a separable normalized isogeny, furthermore, degree must be the degree of the isogeny from E to codomain. If kernel is not None, then this must be isomorphic to the codomain of the normalized separable isogeny defined by kernel, in this case, the isogeny is post composed with an isomorphism so that this parameter is the codomain.
If kernel is None, then this is the degree of the isogeny from E to codomain. If kernel is not None, then this is used to determine whether or not to skip a gcd of the kernel polynomial with the two torsion polynomial of E.
E is a curve over the rationals, then the codomain is set to be the unique global minimum model.
The valid values are “velu” and “kohel”. If “velu” is set, then kernel must be a list of points in E that define a kernel of an isogeny. If “kohel” is set, then the kernel must be either a kernel polynomial or a list of coefficients of a kernel polynomial.
EXAMPLES:
sage: E = EllipticCurve(GF(2), [0,0,1,0,1])
sage: phi = EllipticCurveIsogeny(E, [1,1], algorithm="kohel")
sage: phi
Isogeny of degree 3 from Elliptic Curve defined by y^2 + y = x^3 + 1 over Finite Field of size 2 to Elliptic Curve defined by y^2 + y = x^3 over Finite Field of size 2
sage: E = EllipticCurve(GF(31), [0,0,0,1,0])
sage: P = E((2,17))
sage: ker_list = [j*P for j in xrange(8)]
sage: phi = EllipticCurveIsogeny(E, ker_list, algorithm="velu")
sage: phi
Isogeny of degree 8 from Elliptic Curve defined by y^2 = x^3 + x over Finite Field of size 31 to Elliptic Curve defined by y^2 = x^3 + 10*x + 28 over Finite Field of size 31
sage: E = EllipticCurve('17a1')
sage: phi = EllipticCurveIsogeny(E, [41/3, -55, -1, -1, 1], algorithm="kohel")
sage: phi
Isogeny of degree 9 from Elliptic Curve defined by y^2 + x*y + y = x^3 - x^2 - x - 14 over Rational Field to Elliptic Curve defined by y^2 + x*y + y = x^3 - x^2 - 56*x - 10124 over Rational Field
sage: E = EllipticCurve('37a1')
sage: triv = EllipticCurveIsogeny(E, [E(0)], algorithm="velu")
sage: triv
Isogeny of degree 1 from Elliptic Curve defined by y^2 + y = x^3 - x over Rational Field to Elliptic Curve defined by y^2 + y = x^3 - x over Rational Field
sage: triv.rational_maps()
(x, y)
Function to implement unary negation (-) operator on isogenies. Returns a copy of this isogeny that has been negated.
EXAMPLES:
The following examples inherently exercise this function:
sage: E = EllipticCurve(j=GF(17)(0))
sage: phi = EllipticCurveIsogeny(E, [E(0), E((-1,0))])
sage: negphi = -phi
sage: phi(E((0,1))) + negphi(E((0,1)))
(0 : 1 : 0)
sage: E = EllipticCurve(j=GF(19)(1728))
sage: R.<x> = GF(19)[]
sage: phi = EllipticCurveIsogeny(E, x)
sage: negphi = -phi
sage: phi(E((3,7))) + negphi(E((3,12))) == phi(2*E((3,7)))
True
sage: negphi(E((18,6)))
(17 : 0 : 1)
sage: R.<x> = QQ[]
sage: E = EllipticCurve('17a1')
sage: R.<x> = QQ[]
sage: f = x - 11/4
sage: phi = EllipticCurveIsogeny(E, f)
sage: negphi = -phi
sage: phi.rational_maps()[0] == negphi.rational_maps()[0]
True
sage: P = E((7,13))
sage: phi(P) + negphi(P)
(0 : 1 : 0)
Composition operator function inherited from morphism class.
EXAMPLES:
sage: E = EllipticCurve(j=GF(7)(0))
sage: phi = EllipticCurveIsogeny(E, [E(0), E((0,1)), E((0,-1))])
sage: phi*phi
...
NotImplementedError
sage: phi._composition_(phi, phi.parent())
...
NotImplementedError
Special sage specific function that implements functionality to display an isogeny object as a latex string.
This function returns a latex string representing the isogeny self as the and
coordinate rational functions.
EXAMPLES:
sage: E = EllipticCurve(QQ, [0,0,0,1,-1])
sage: phi = EllipticCurveIsogeny(E, [E(0)], algorithm="velu")
sage: phi._latex_()
'\\left( x , y \\right)'
sage: E = EllipticCurve(GF(17), [0,0,0,1,-1])
sage: phi = EllipticCurveIsogeny(E, [0,1], algorithm="kohel")
sage: phi._latex_()
'\\left( \\frac{x^{3} + 2 x + 13}{x^{2}} , \\frac{x^{3} y + 15 x y + 8 y}{x^{3}} \\right)'
Special sage specific function that implement the functionality to display the isogeny self as a string.
EXAMPLES:
sage: E = EllipticCurve(GF(31), [1,0,1,1,0])
sage: phi = EllipticCurveIsogeny(E, [j*E((0,0)) for j in xrange(17)], algorithm="velu")
sage: phi._repr_()
'Isogeny of degree 17 from Elliptic Curve defined by y^2 + x*y + y = x^3 + x over Finite Field of size 31 to Elliptic Curve defined by y^2 + x*y + y = x^3 + 14*x + 9 over Finite Field of size 31'
sage: E = EllipticCurve(QQ, [1,0,0,1,9])
sage: phi = EllipticCurveIsogeny(E, [2,1], algorithm="kohel")
sage: phi._repr_()
'Isogeny of degree 2 from Elliptic Curve defined by y^2 + x*y = x^3 + x + 9 over Rational Field to Elliptic Curve defined by y^2 + x*y = x^3 - 59*x + 165 over Rational Field'
Returns the codomain (range) curve of this isogeny.
EXAMPLES:
sage: E = EllipticCurve(QQ, [0,0,0,1,0])
sage: phi = EllipticCurveIsogeny(E, [E(0), E((0,0))], algorithm="velu")
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, [14, 27, 4, 1], algorithm="kohel")
sage: phi.codomain()
Elliptic Curve defined by y^2 + x*y = x^3 + 15*x + 8 over Finite Field of size 31
Returns the degree of this isogeny.
EXAMPLES:
sage: E = EllipticCurve(QQ, [0,0,0,1,0])
sage: phi = EllipticCurveIsogeny(E, [E(0), E((0,0))], algorithm="velu")
sage: phi.degree()
2
sage: phi = EllipticCurveIsogeny(E, [0,1,0,1], algorithm="kohel")
sage: phi.degree()
4
sage: E = EllipticCurve(GF(31), [1,0,0,1,2])
sage: phi = EllipticCurveIsogeny(E, [14, 27, 4, 1], algorithm="kohel")
sage: phi.degree()
7
Returns the domain curve of this isogeny.
EXAMPLES:
sage: E = EllipticCurve(QQ, [0,0,0,1,0])
sage: phi = EllipticCurveIsogeny(E, [E(0), E((0,0))], algorithm="velu")
sage: phi.domain() == E
True
sage: E = EllipticCurve(GF(31), [1,0,0,1,2])
sage: phi = EllipticCurveIsogeny(E, [14, 27, 4, 1], algorithm="kohel")
sage: phi.domain()
Elliptic Curve defined by y^2 + x*y = x^3 + x + 2 over Finite Field of size 31
Computes and returns the dual isogeny of this isogeny.
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)
False
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)
False
sage: (Xm, -Ym) == E.multiplication_by_m(5)
True
Returns the post-isomorphism of this isogeny. If there has been no post-isomorphism set, this returns 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 + 77*y = x^3 + 49*x + 28 over Finite Field of size 83
Via: (u,r,s,t) = (82, 7, 41, 3)
Returns the pre-isomorphism of this isogeny. If there has been no pre-isomorphism set, this returns 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)
Method inherited from the morphism class. Returns True iff 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), E((0,-1)), E((0,1))])
sage: phi.is_injective()
False
sage: phi = EllipticCurveIsogeny(E, [E(0)])
sage: phi.is_injective()
True
Returns true if this isogeny is normalized.
This code does two things. If the check_by_pullback flag is set, then the code checks that the coefficient on the pullback of the invariant differential is 1. However, in some cases (after a translation has been applied) the underlying polynomial algebra code can not sufficiently simplify the pullback expression. As such, we also cheat a little by falling back and seeing if the post isomorphism on this isogeny is a translation with no rescaling.
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
This function returns a bool indicating whether or not this isogeny is separable.
This function always returns true. This class only implements separable isogenies.
EXAMPLES:
sage: E = EllipticCurve(GF(17), [0,0,0,3,0])
sage: phi = EllipticCurveIsogeny(E, [E(0), E((0,0))], algorithm="velu")
sage: phi.is_separable()
True
sage: E = EllipticCurve('11a1')
sage: phi = EllipticCurveIsogeny(E, E.torsion_points(), algorithm="velu")
sage: phi.is_separable()
True
For elliptic curve isogenies, always returns True (as a nonconstant map of algebraic curves must be surjective).
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), 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
Member function inherited from morphism class.
EXAMPLES:
sage: E = EllipticCurve(j=GF(7)(0))
sage: phi = EllipticCurveIsogeny(E, [E(0), E((0,1)), E((0,-1))])
sage: phi.is_zero()
...
NotImplementedError
Returns the kernel polynomial of this isogeny.
EXAMPLES:
sage: E = EllipticCurve(QQ, [0,0,0,2,0])
sage: phi = EllipticCurveIsogeny(E, [E(0), E((0,0))], algorithm="velu")
sage: phi.kernel_polynomial()
x
sage: E = EllipticCurve('11a1')
sage: phi = EllipticCurveIsogeny(E, E.torsion_points(), algorithm="velu")
sage: phi.kernel_polynomial()
x^2 - 21*x + 80
sage: E = EllipticCurve(GF(17), [1,-1,1,-1,1])
sage: phi = EllipticCurveIsogeny(E, [1], algorithm="kohel")
sage: phi.kernel_polynomial()
1
sage: E = EllipticCurve(GF(31), [0,0,0,3,0])
sage: phi = EllipticCurveIsogeny(E, [0,3,0,1], algorithm="kohel")
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), E((0,1)), E((0,-1))])
sage: phi.n()
...
NotImplementedError: Numerical approximations do not make sense for Elliptic Curve Isogenies
Member function inherited from morphism class.
EXAMPLES:
sage: E = EllipticCurve(j=GF(7)(0))
sage: phi = EllipticCurveIsogeny(E, [E(0), E((0,1)), E((0,-1))])
sage: phi.post_compose(phi)
...
NotImplementedError
Member function inherited from morphism class.
EXAMPLES:
sage: E = EllipticCurve(j=GF(7)(0))
sage: phi = EllipticCurveIsogeny(E, [E(0), E((0,1)), E((0,-1))])
sage: phi.pre_compose(phi)
...
NotImplementedError
This function returns this isogeny as a pair of rational maps.
EXAMPLES:
sage: E = EllipticCurve(QQ, [0,2,0,1,-1])
sage: phi = EllipticCurveIsogeny(E, [1], algorithm="kohel")
sage: phi.rational_maps()
(x, y)
sage: E = EllipticCurve(GF(17), [0,0,0,3,0])
sage: phi = EllipticCurveIsogeny(E, [E(0), E((0,0))], algorithm="velu")
sage: phi.rational_maps()
((x^2 + 3)/x, (x^2*y - 3*y)/x^2)
Modifies this isogeny object to post compose with the given 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
Modifies this isogeny object to pre compose with the given 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)
This function composes the isogeny with (flipping the
coefficient between +/-1 on the
coordinate rational map).
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))
Given parameters and
(as in velu / kohel / etc formulas)
computes the codomain curve.
EXAMPLES:
This formula is used by every Isogeny Instantiation:
sage: E = EllipticCurve(GF(19), [1,2,3,4,5])
sage: phi = EllipticCurveIsogeny(E, [j*E((1,2)) for j in xrange(8)])
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)
Elliptic Curve defined by y^2 + x*y + 3*y = x^3 + 2*x^2 + 9*x + 13 over Finite Field of size 19
This function computes the codomain from the kernel polynomial as per Kohel’s formulas.
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], algorithm="kohel")
sage: phi.codomain() == isogeny_codomain_from_kernel(E, [9,1], algorithm="kohel")
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, algorithm="kohel")
sage: phi.codomain() == isogeny_codomain_from_kernel(E, x^3 + 14*x^2 + 3*x + 11, algorithm="kohel")
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,algorithm="kohel")
sage: isogeny_codomain_from_kernel(E, x^3 + 7*x^2 + 15*x + 12,algorithm="kohel") == 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
NOTES:
This function uses the formulas of Section 2.4 of [K96].
REFERENCES:
Computes isomorphism from E1 to an intermediate domain and an isomorphism from an intermediate codomain to E2.
Intermediate domain and intermediate codomain, are in short weierstrass form.
This is used so we can compute functions from the short weierstrass model more easily.
The underlying 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 + 77*y = x^3 + 49*x + 28 over Finite Field of size 83
Via: (u,r,s,t) = (1, 7, 42, 80))
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 morphism:
From: 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
To: 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 morphism:
From: 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
To: 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))
Computes the degree ell isogeny between E1 and E2.
There must be a degree ell, separable, normalized isogeny from E1 to E2. If no algorithm is specified, this function determines the best algorithm to use.
INPUT:
E1 - an elliptic curve in short weierstrass form.
E2 - an elliptic curve in short weierstrass form.
ell - the degree of the isogeny from E1 to E2.
Otherwise uses the algorithm specified by the string. Valid values are “starks” or “fastElkies”
OUTPUT:
polynomial – over the field of definition of E1, E2, that is the kernel polynomial of the isogeny from E1 to E2.
Note
When using Stark’s algorithm, occasionally the fast pe computation fails, so we retry with the quadratic algorithm, which works in all cases (for valid inputs.)
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)
z^2 + 7*z + 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)
z^3 + z
Computes the degree ell isogeny between E1 and E2 via Stark’s algorithm. There must be a degree ell, separable, normalized isogeny from E1 to E2.
INPUT:
OUTPUT:
polynomial – over the field of definition of E1, E2, that is the kernel polynomial of the isogeny from E1 to E2.
ALGORITHM:
This function uses Starks Algorithm as presented in section 6.2 of [BMSS].
Note
As published there, the algorithm is incorrect, and a correct version (with slightly different notation) can be found in [M09]. The algorithm originates in [S72]
REFERENCES:
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, pe_algorithm="quadratic")
z^10 + 37*z^9 + 53*z^8 + 66*z^7 + 66*z^6 + 17*z^5 + 57*z^4 + 6*z^3 + 89*z^2 + 53*z + 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)
z^4 + 14*z^3 + z^2 + 34*z + 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, pe_algorithm="fast")
z
Computes the truncated weierstrass function on an elliptic curve defined by short weierstrass model: .
Uses the algorithm specified by the algorithm parameter.
INPUT:
OUTPUT:
polynomial - a polynomial corresponding to the truncated weierstrass function in R.
EXAMPLES:
sage: from sage.schemes.elliptic_curves.ell_curve_isogeny import compute_pe
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: pe = compute_pe(R, E, 7, algorithm="quadratic")
sage: pe = pe(x^2)
sage: pe
(7*x^14 + 9*x^10 + x^8 + 20*x^6 + 22*x^4 + 1)/x^2
Computes the truncated weierstrass function of an elliptic curve defined by short weierstrass model: .
It does this with time complexity
.
Let be the characteristic of the underlying field: Then we must have either
, or
.
INPUT:
- poly_ring - polynomial ring, to compute the function in (assumes that the generator is
for efficiency of storage/operations.)
- A - field element corresponding to the
coefficient in the weierstrass equation of an elliptic curve
- B - field element corresponding to the constant coefficient in the weierstrass equation of an elliptic curve
- ell - degree of
to compute the truncated function to. If
is the characteristic of the underlying field and
, then we must have
.
OUTPUT:
polynomial – the element in poly_ring that corresponds to the truncated function to precision
.
ALGORITHM:
This function uses the algorithm described in section 3.3 of [BMSS].
Note
Occasionally this function will fail to give the right answer, it faithfully implements the above published algorithm, so compute_pe_quadratic should be used as a fallback.
REFERENCES:
EXAMPLES:
sage: from sage.schemes.elliptic_curves.ell_curve_isogeny import compute_pe_fast
sage: R.<x> = QQ[]
sage: compute_pe_fast(R, 1, 8, 7)
(-16/5775*x^7 + 23902/238875*x^6 + 24/385*x^5 + 1/75*x^4 - 8/7*x^3 - 1/5*x^2 + 1)/x
sage: R.<x> = GF(37)[]
sage: compute_pe_fast(R, GF(37)(1), GF(37)(8), 5)
(9*x^5 + x^4 + 20*x^3 + 22*x^2 + 1)/x
Computes the truncated weierstrass function of an elliptic curve defined by short weierstrass model: .
Uses an algorithm that is of complexity
.
Let be the characteristic of the underlying field: Then we must have either
, or
.
INPUT:
- poly_ring - polynomial ring, to compute the
function in (assumes that the generator is
for efficiency of storage/operations.)
- A - field element corresponding to the
coefficient in the weierstrass equation of an elliptic curve
- B - field element corresponding to the constant coefficient in the weierstrass equation of an elliptic curve
- ell - degree of
to compute the truncated function to. If
is the characteristic of the underlying field. If
then we must have
.
OUTPUT:
polynomial – the element in poly_ring that corresponds to the truncated function to precision .
ALGORITHM:
This function uses the algorithm described in section 3.2 of [BMSS].
REFERENCES:
EXAMPLES:
sage: from sage.schemes.elliptic_curves.ell_curve_isogeny import compute_pe_quadratic
sage: R.<x> = GF(37)[]
sage: compute_pe_quadratic(R, GF(37)(1), GF(37)(8), 7)
(7*x^7 + 9*x^5 + x^4 + 20*x^3 + 22*x^2 + 1)/x
sage: R.<x> = QQ[]
sage: compute_pe_quadratic(R, 1, 8, 7)
(-16/5775*x^7 + 23902/238875*x^6 + 24/385*x^5 + 1/75*x^4 - 8/7*x^3 - 1/5*x^2 + 1)/x
Given domain E1 and codomain E2 such that there is a degree ell separable normalized isogeny from E1 to E2, returns pre/post isomorphism, as well as intermediate domain and codomain, and kernel polynomial.
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,
z^2 - 61/3*z + 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 morphism:
From: 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
To: 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 morphism:
From: 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
To: 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,
z^3 + z)
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 + 9*y = x^3 + 83*x + 6 over Finite Field of size 97
Via: (u,r,s,t) = (1, 89, 49, 53),
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,
z^5 + 67*z^4 + 13*z^3 + 35*z^2 + 77*z + 69)
The formula for computing and
using Kohel’s formulas for isogenies of degree 2.
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], algorithm="kohel")
sage: 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)
The formula for computing and
using Kohel’s formulas for isogenies of degree 3.
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,algorithm="kohel")
sage: 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)
This function computes the and
according to Kohel’s formulas.
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, algorithm="kohel")
sage: 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)
This function computes the isogeny codomain given a kernel.
INPUT:
OUTPUT:
(elliptic curve) the codomain of the separable normalized isogeny from this kernel
EXAMPLES:
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.
If kernel is a list of points on the EllipticCurve , then we assume the algorithm to use is Velu.
If kernel is a list of coefficients or a univariate polynomial we 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])
sage: phi = EllipticCurveIsogeny(E, x+3)
sage: phi2 = EllipticCurveIsogeny(E, [GF(5)(3),GF(5)(1)])
sage: phi == phi2
True
sage: phi3 = EllipticCurveIsogeny(E, [E(0), E((2,0))])
sage: phi3 == phi2
True
sage: from sage.schemes.elliptic_curves.ell_curve_isogeny import isogeny_determine_algorithm
sage: isogeny_determine_algorithm(E, x+3, None, None, None)
'kohel'
sage: isogeny_determine_algorithm(E, [3, 1], None, None, None)
'kohel'
sage: isogeny_determine_algorithm(E, [E(0), E((2,0))], None, None, None)
'velu'
Solves a system of linear differential equations:
,
where
,
,
,
,
are polynomials in variable
, and
,
are computed to precision
.
EXAMPLES:
The following examples inherently exercises this function:
sage: from sage.schemes.elliptic_curves.ell_curve_isogeny import compute_isogeny_kernel_polynomial, solve_linear_differential_system
sage: R.<x> = GF(47)[]
sage: E = EllipticCurve(GF(47), [0,0,0,1,0])
sage: E2 = EllipticCurve(GF(47), [0,0,0,16,0])
sage: compute_isogeny_kernel_polynomial(E, E2, 4, algorithm="starks")
z^3 + z
Internal helper function for compute_isogeny_kernel_polynomial.
Given a full kernel polynomial (where two torsion -coordinates
are roots of multiplicity 1, and all other roots have multiplicity
2.) of degree
, returns the maximum separable divisor.
(i.e. the kernel polynomial with roots of multiplicity at most 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, pe_algorithm="quadratic"); ker_poly
z^6 + 2*z^5 + 20*z^4 + 11*z^3 + 36*z^2 + 35*z + 16
sage: split_kernel_polynomial(E, ker_poly, 7)
z^3 + z^2 + 28*z + 33
Helper function for starks algorithm.
INPUT:
OUTPUT:
tuple – where
is the largest exponent such that the
coefficient
of
is nonzero, where
is regarded as
a polynomial in
.
EXAMPLES:
sage: from sage.schemes.elliptic_curves.ell_curve_isogeny import starks_find_r_and_t
sage: from sage.schemes.elliptic_curves.ell_curve_isogeny import compute_pe
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: pe = compute_pe(R, E, 7, algorithm="quadratic")
sage: pe = pe(x^2)
sage: starks_find_r_and_t(pe, x)
(2, 1)
Computes the truncated exponential of , to precision
, using an efficient newton iteration.
This algorithm has complexity
, where
is the cost of a multiplication.
INPUT:
OUTPUT:
polynomial – a polynomial , such that
.
ALGORITHM:
algorithm “fast” uses newton iteration and has complexity O(M(n)), algorithm “quadratic” has complexity O(n*M(n))
EXAMPLES:
sage: from sage.schemes.elliptic_curves.ell_curve_isogeny import truncated_exp
sage: R.<x> = GF(29)[]
sage: f = 1+x+3*x^3
sage: truncated_exp(f, 3)
x + 2
Computes the truncated exponential of , to precision
, using an efficient newton iteration.
This algorithm has complexity
, where
is the cost of a multiplication.
INPUT:
OUTPUT:
polynomial – a polynomial , such that
.
ALGORITHM:
This function uses the algorithm described in section 2.2 of [BMSS]
REFERENCES:
EXAMPLES:
sage: from sage.schemes.elliptic_curves.ell_curve_isogeny import truncated_exp_fast
sage: R.<x> = GF(27, 'a')[]
sage: f = 1 + x^2 + 2*x^3 + x^4
sage: truncated_exp_fast(f, 4)
x^3 + 2
sage: truncated_exp_fast(f, 3)
2*x^2 + 2
Computes the truncated exponential of , to precision
, using a straight forward power series approximation.
This algorithm has complexity , where
is the
cost of a multiplication.
INPUT:
OUTPUT:
polynomial – a polynomial , such that
.
ALGORITHM:
This function uses the algorithm described in section 2.2 of [BMSS]
REFERENCES:
EXAMPLES:
sage: from sage.schemes.elliptic_curves.ell_curve_isogeny import truncated_exp_quadratic
sage: R.<x> = GF(127)[]
sage: f = 1 + x^10
sage: truncated_exp_quadratic(f, 10)
31
sage: truncated_exp_quadratic(f, 11)
31*x^10 + 123
Computes the truncated logarithm of polynomial to precision
.
The complexity of this function is
.
INPUT:
OUTPUT:
polynomial – a polynomial , such that
ALGORITHM:
Uses the algorithm from section 2.2 of [BMSS].
REFERENCES:
EXAMPLES:
sage: from sage.schemes.elliptic_curves.ell_curve_isogeny import truncated_log, truncated_exp
sage: R.<x> = GF(31)[]
sage: f = 1 + x + 2*x^2 + 3*x^3
sage: truncated_log(f, 4)
3*x^3 + 2*x^2 + x
sage: g = truncated_log(f, 4)
sage: truncated_exp(g, 4)
3*x^3 + 2*x^2 + x + 1
Computes the truncated reciprocal of , to precision
.
INPUT:
OUTPUT:
polynomial – a polynomial , such that
ALGORITHM:
“newtoniteration” algorithm has complexity and “quadratic” has algorithm complexity
.
EXAMPLES:
sage: from sage.schemes.elliptic_curves.ell_curve_isogeny import truncated_reciprocal
sage: R.<x> = GF(101)[]
sage: f = 4 + x - x^2 + 50*x^3 + 91*x^5
sage: truncated_reciprocal(f, 5)
37*x^4 + 43*x^3 + 49*x^2 + 82*x + 76
sage: (f*truncated_reciprocal(f, 5)).quo_rem(x**5)[1]
1
Computes the truncated reciprocal of , to precision
by newton iteration.
This algorithm has complexity
, where
is the cost of a multiplication.
INPUT:
OUTPUT:
polynomial – polynomial , such that
ALGORITHM:
This function uses the algorithm described in section 2.1 of [BMSS]
REFERENCES:
EXAMPLES:
sage: from sage.schemes.elliptic_curves.ell_curve_isogeny import truncated_reciprocal_newton
sage: R.<x> = GF(53)[]
sage: f = 4 + 3*x^3 + 47*x^7
sage: truncated_reciprocal_newton(f, 6)
23*x^3 + 40
sage: (f*truncated_reciprocal_newton(f, 6)).quo_rem(x**6)[1]
1
Computes the truncated reciprocal of , to precision
.
This algorithm has complexity
, where
is the cost of a multiplication and
is the degree of
.
INPUT:
OUTPUT:
polynomial – a polynomial , such that
ALGORITHM:
This function uses the algorithm described in section 2.1 of [BMSS]
REFERENCES:
EXAMPLES:
sage: f = 1 + x + 7*x^2 + 9*x^5
sage: from sage.schemes.elliptic_curves.ell_curve_isogeny import truncated_reciprocal_quadratic
sage: R.<x> = QQ[]
sage: f = 1 + x + 7*x^2 + 9*x^5
sage: truncated_reciprocal_quadratic(f, 5)
29*x^4 + 13*x^3 - 6*x^2 - x + 1
sage: (f*truncated_reciprocal_quadratic(f, 5)).quo_rem(x**5)[1]
1
Returns the greatest common divisor 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,algorithm="kohel")
sage: isogeny_codomain_from_kernel(E, x + 13,algorithm="kohel") == phi.codomain()
True
sage: from sage.schemes.elliptic_curves.ell_curve_isogeny import two_torsion_part
sage: two_torsion_part(E, R, x+13, 2)
x + 13