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
59 changes: 59 additions & 0 deletions src/sage/schemes/elliptic_curves/ell_field.py
Original file line number Diff line number Diff line change
Expand Up @@ -1588,6 +1588,65 @@
from .isogeny_small_degree import isogenies_prime_degree
return sum([isogenies_prime_degree(self, d) for d in L], [])

def isogenies_degree(self, n):
r"""
Return a list of all separable isogenies of given degree with domain
equal to ``self``, which are defined over the base field of ``self``.
grhkm21 marked this conversation as resolved.
Show resolved Hide resolved

INPUT:

- ``n`` -- an integer.

TESTS::

sage: E = EllipticCurve(GF(11), [1, 1])
sage: E.isogenies_degree(23 * 19)
[Composite morphism of degree 437 = 1*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 = 1*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 = 1*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 = 1*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: pol = PolynomialRing(QQ, 'x')([1, -3, 5, -5, 5, -3, 1])

Check warning on line 1619 in src/sage/schemes/elliptic_curves/ell_field.py

View check run for this annotation

Codecov / codecov/patch

src/sage/schemes/elliptic_curves/ell_field.py#L1619

Added line #L1619 was not covered by tests
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(E.isogenies_degree(2**2))
7
sage: len(E.isogenies_degree(2**5)) # long time (15s)
99
"""
n = rings.Integer(n)
grhkm21 marked this conversation as resolved.
Show resolved Hide resolved
if n.is_prime():
grhkm21 marked this conversation as resolved.
Show resolved Hide resolved
return self.isogenies_prime_degree(n)

prime_divisors = sum([[p] * e for p, e in n.factor()], [])

Check warning on line 1634 in src/sage/schemes/elliptic_curves/ell_field.py

View check run for this annotation

Codecov / codecov/patch

src/sage/schemes/elliptic_curves/ell_field.py#L1634

Added line #L1634 was not covered by tests
grhkm21 marked this conversation as resolved.
Show resolved Hide resolved
isos = [self.isogeny(self(0))]
grhkm21 marked this conversation as resolved.
Show resolved Hide resolved

for p in prime_divisors:
if not isos:
break

new_isos = []
for iso in isos:
Eiso = iso.codomain()
for next_iso in Eiso.isogenies_prime_degree(p):
new_isos.append(next_iso * iso)
isos = new_isos

return isos

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