Skip to content

Commit cf322e4

Browse files
author
Release Manager
committed
gh-37081: convenience method to construct elliptic-curve isomorphisms from u,r,s,t The `EllipticCurve_generic` class has methods for computing an isomorphism between two given curves as well as computing the codomain of a particular isomorphism (specified by its coefficients $u,r,s,t$). What it does not currently have is a method for computing an isomorphism from its (co)domain and the coefficients $u,r,s,t$. The best way to construct such an isomorphism is currently to `from sage.schemes.elliptic_curves.weierstrass_morphism import WeierstrassIsomorphism`, which is less than convenient. In this patch we add the desired method. #sd123 URL: #37081 Reported by: Lorenz Panny Reviewer(s): Travis Scrimshaw
2 parents e788164 + 89f8531 commit cf322e4

File tree

1 file changed

+88
-0
lines changed

1 file changed

+88
-0
lines changed

src/sage/schemes/elliptic_curves/ell_generic.py

Lines changed: 88 additions & 0 deletions
Original file line numberDiff line numberDiff line change
@@ -1555,6 +1555,94 @@ def scale_curve(self, u):
15551555
u = self.base_ring()(u) # because otherwise 1/u would round!
15561556
return self.change_weierstrass_model(1/u, 0, 0, 0)
15571557

1558+
def isomorphism(self, u, r=0, s=0, t=0, *, is_codomain=False):
1559+
r"""
1560+
Given four values `u,r,s,t` in the base ring of this curve, return
1561+
the :class:`WeierstrassIsomorphism` defined by `u,r,s,t` with this
1562+
curve as its codomain.
1563+
(The value `u` must be a unit; the values `r,s,t` default to zero.)
1564+
1565+
Optionally, if the keyword argument ``is_codomain`` is set to ``True``,
1566+
return the isomorphism defined by `u,r,s,t` with this curve as its
1567+
**co**\domain.
1568+
1569+
EXAMPLES::
1570+
1571+
sage: E = EllipticCurve([1, 2, 3, 4, 5])
1572+
sage: iso = E.isomorphism(6); iso
1573+
Elliptic-curve morphism:
1574+
From: Elliptic Curve defined by y^2 + x*y + 3*y = x^3 + 2*x^2 + 4*x + 5 over Rational Field
1575+
To: Elliptic Curve defined by y^2 + 1/6*x*y + 1/72*y = x^3 + 1/18*x^2 + 1/324*x + 5/46656 over Rational Field
1576+
Via: (u,r,s,t) = (6, 0, 0, 0)
1577+
sage: iso.domain() == E
1578+
True
1579+
sage: iso.codomain() == E.scale_curve(1 / 6)
1580+
True
1581+
1582+
sage: iso = E.isomorphism(1, 7, 8, 9); iso
1583+
Elliptic-curve morphism:
1584+
From: Elliptic Curve defined by y^2 + x*y + 3*y = x^3 + 2*x^2 + 4*x + 5 over Rational Field
1585+
To: Elliptic Curve defined by y^2 + 17*x*y + 28*y = x^3 - 49*x^2 - 54*x + 303 over Rational Field
1586+
Via: (u,r,s,t) = (1, 7, 8, 9)
1587+
sage: iso.domain() == E
1588+
True
1589+
sage: iso.codomain() == E.rst_transform(7, 8, 9)
1590+
True
1591+
1592+
sage: iso = E.isomorphism(6, 7, 8, 9); iso
1593+
Elliptic-curve morphism:
1594+
From: Elliptic Curve defined by y^2 + x*y + 3*y = x^3 + 2*x^2 + 4*x + 5 over Rational Field
1595+
To: Elliptic Curve defined by y^2 + 17/6*x*y + 7/54*y = x^3 - 49/36*x^2 - 1/24*x + 101/15552 over Rational Field
1596+
Via: (u,r,s,t) = (6, 7, 8, 9)
1597+
sage: iso.domain() == E
1598+
True
1599+
sage: iso.codomain() == E.rst_transform(7, 8, 9).scale_curve(1 / 6)
1600+
True
1601+
1602+
The ``is_codomain`` argument reverses the role of domain and codomain::
1603+
1604+
sage: E = EllipticCurve([1, 2, 3, 4, 5])
1605+
sage: iso = E.isomorphism(6, is_codomain=True); iso
1606+
Elliptic-curve morphism:
1607+
From: Elliptic Curve defined by y^2 + 6*x*y + 648*y = x^3 + 72*x^2 + 5184*x + 233280 over Rational Field
1608+
To: Elliptic Curve defined by y^2 + x*y + 3*y = x^3 + 2*x^2 + 4*x + 5 over Rational Field
1609+
Via: (u,r,s,t) = (6, 0, 0, 0)
1610+
sage: iso.domain() == E.scale_curve(6)
1611+
True
1612+
sage: iso.codomain() == E
1613+
True
1614+
1615+
sage: iso = E.isomorphism(1, 7, 8, 9, is_codomain=True); iso
1616+
Elliptic-curve morphism:
1617+
From: Elliptic Curve defined by y^2 - 15*x*y + 90*y = x^3 - 75*x^2 + 796*x - 2289 over Rational Field
1618+
To: Elliptic Curve defined by y^2 + x*y + 3*y = x^3 + 2*x^2 + 4*x + 5 over Rational Field
1619+
Via: (u,r,s,t) = (1, 7, 8, 9)
1620+
sage: iso.domain().rst_transform(7, 8, 9) == E
1621+
True
1622+
sage: iso.codomain() == E
1623+
True
1624+
1625+
sage: iso = E.isomorphism(6, 7, 8, 9, is_codomain=True); iso
1626+
Elliptic-curve morphism:
1627+
From: Elliptic Curve defined by y^2 - 10*x*y + 700*y = x^3 + 35*x^2 + 9641*x + 169486 over Rational Field
1628+
To: Elliptic Curve defined by y^2 + x*y + 3*y = x^3 + 2*x^2 + 4*x + 5 over Rational Field
1629+
Via: (u,r,s,t) = (6, 7, 8, 9)
1630+
sage: iso.domain().rst_transform(7, 8, 9) == E.scale_curve(6)
1631+
True
1632+
sage: iso.codomain() == E
1633+
True
1634+
1635+
.. SEEALSO::
1636+
1637+
- :class:`~sage.schemes.elliptic_curves.weierstrass_morphism.WeierstrassIsomorphism`
1638+
- :meth:`rst_transform`
1639+
- :meth:`scale_curve`
1640+
"""
1641+
from sage.schemes.elliptic_curves.weierstrass_morphism import WeierstrassIsomorphism
1642+
if is_codomain:
1643+
return WeierstrassIsomorphism(None, (u,r,s,t), self)
1644+
return WeierstrassIsomorphism(self, (u,r,s,t))
1645+
15581646
# ###########################################################
15591647
#
15601648
# Explanation of the division (also known as torsion) polynomial

0 commit comments

Comments
 (0)