Skip to content
New issue

Have a question about this project? Sign up for a free GitHub account to open an issue and contact its maintainers and the community.

By clicking “Sign up for GitHub”, you agree to our terms of service and privacy statement. We’ll occasionally send you account related emails.

Already on GitHub? Sign in to your account

Compute composite degree (separable) isogenies of EllipticCurves #35949

Open
wants to merge 10 commits into
base: develop
Choose a base branch
from
113 changes: 110 additions & 3 deletions src/sage/schemes/elliptic_curves/ell_field.py
Original file line number Diff line number Diff line change
Expand Up @@ -1598,9 +1598,9 @@ def kernel_polynomial_from_divisor(self, f, l, *, check=True):

def isogenies_prime_degree(self, l=None, max_l=31):
"""
Return a list of all separable isogenies of given prime degree(s)
with domain equal to ``self``, which are defined over the base
field of ``self``.
Return a list of all separable isogenies (up to post-composition with
isomorphisms) of given prime degree(s) with domain equal to ``self``,
which are defined over the base field of ``self``.

INPUT:

Expand Down Expand Up @@ -1882,6 +1882,113 @@ def isogenies_prime_degree(self, l=None, max_l=31):
from .isogeny_small_degree import isogenies_prime_degree
return sum([isogenies_prime_degree(self, d) for d in L], [])

def isogenies_degree(self, n, *, _intermediate=False):
r"""
Return an iterator of all separable isogenies of given degree (up to
post-composition with isomorphisms) with domain equal to ``self``,
which are defined over the base field of ``self``.

ALGORITHM:

The prime factors `p` of `n` are processed one by one, each time
"branching" out by taking isogenies of degree `p`.

INPUT:

- ``n`` -- integer, or its
:class:`~sage.structure.factorization.Factorization`

- ``_intermediate`` -- (bool, default: False): If set, the curves
traversed within the depth-first search are returned. This is for
internal use only.

EXAMPLES::

sage: E = EllipticCurve(GF(11), [1, 1])
sage: list(E.isogenies_degree(23 * 19))
[Composite morphism of degree 437 = 19*23:
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 + 8*x + 7 over Finite Field of size 11,
Composite morphism of degree 437 = 19*23:
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 + 2*x + 6 over Finite Field of size 11,
Composite morphism of degree 437 = 19*23:
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 + 6*x + 2 over Finite Field of size 11,
Composite morphism of degree 437 = 19*23:
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]

grhkm21 marked this conversation as resolved.
Show resolved Hide resolved
::

sage: E = EllipticCurve(GF(next_prime(2^32)), j=1728)
sage: sorted([phi.codomain().j_invariant() for phi in E.isogenies_degree(11 * 17 * 19^2)])
[1348157279, 1348157279, 1713365879, 1713365879, 3153894341, 3153894341, 3153894341,
3153894341, 3225140514, 3225140514, 3673460198, 3673460198, 3994312564, 3994312564,
3994312564, 3994312564]
sage: it = E.isogenies_degree(2^2); it
<generator object EllipticCurve_field.isogenies_degree at 0x...>
sage: all(phi.degree() == 2^2 for phi in it)
True

::

sage: pol = PolynomialRing(QQ, 'x')([1, -3, 5, -5, 5, -3, 1])
sage: L.<a> = NumberField(pol)
sage: js = hilbert_class_polynomial(-23).roots(L, multiplicities=False)
sage: len(js)
3
sage: E = EllipticCurve(j=js[0])
sage: len(list(E.isogenies_degree(2**2)))
7
sage: len(list(E.isogenies_degree(2**5))) # long time (15s)
99

TESTS::

sage: E = EllipticCurve(GF(next_prime(2^32)), j=1728)
sage: list(E.isogenies_degree(2^2, _intermediate=True))
[Elliptic-curve endomorphism of Elliptic Curve defined by y^2 = x^3 + x over Finite Field of size 4294967311
Via: (u,r,s,t) = (1, 0, 0, 0),
Isogeny of degree 2 from Elliptic Curve defined by y^2 = x^3 + x over Finite Field of size 4294967311 to Elliptic Curve defined by
y^2 = x^3 + 4294967307*x over Finite Field of size 4294967311,
Isogeny of degree 2 from Elliptic Curve defined by y^2 = x^3 + x over Finite Field of size 4294967311 to Elliptic Curve defined by
y^2 = x^3 + 4294967307*x over Finite Field of size 4294967311,
Composite morphism of degree 4 = 2^2:
From: Elliptic Curve defined by y^2 = x^3 + x over Finite Field of size 4294967311
To: Elliptic Curve defined by y^2 = x^3 + 16*x over Finite Field of size 4294967311,
Composite morphism of degree 4 = 2^2:
From: Elliptic Curve defined by y^2 = x^3 + x over Finite Field of size 4294967311
To: Elliptic Curve defined by y^2 = x^3 + 4294967267*x + 4294967199 over Finite Field of size 4294967311,
Composite morphism of degree 4 = 2^2:
From: Elliptic Curve defined by y^2 = x^3 + x over Finite Field of size 4294967311
To: Elliptic Curve defined by y^2 = x^3 + 4294967267*x + 112 over Finite Field of size 4294967311]

The following curve has no degree-`5` isogenies, so the code is quick::

sage: E = EllipticCurve(GF(103), [3, 5])
sage: list(E.isogenies_degree(5 * product(prime_range(7, 100))))
[]
"""
from sage.structure.factorization import Factorization
from sage.schemes.elliptic_curves.weierstrass_morphism import identity_morphism

if not isinstance(n, Factorization):
n = Integer(n).factor()

if n.value() == 1:
yield identity_morphism(self)
return

p = n[-1][0]
for iso in self.isogenies_degree(n / p, _intermediate=_intermediate):
if _intermediate:
yield iso

Eiso = iso.codomain()
for next_iso in Eiso.isogenies_prime_degree(p):
yield next_iso * iso

def is_isogenous(self, other, field=None):
"""
Return whether or not self is isogenous to other.
Expand Down
Loading