Skip to content

Commit

Permalink
gh-38330: Lattes to curve function code commits for enhancement #38329
Browse files Browse the repository at this point in the history
    
- **Fixing issue #38329**
- **Adding Function with documentation to be reviewed for the above
ticket**
- **Lattes to curve function**

<!-- ^ Please provide a concise and informative title. -->
<!-- ^ Don't put issue numbers in the title, do this in the PR
description below. -->
<!-- ^ For example, instead of "Fixes #12345" use "Introduce new method
to calculate 1 + 2". -->
<!-- v Describe your changes below in detail. -->
<!-- v Why is this change required? What problem does it solve? -->
<!-- v If this PR resolves an open issue, please link to it here. For
example, "Fixes #12345". -->



### 📝 Checklist

<!-- Put an `x` in all the boxes that apply. -->

- [ ] The title is concise and informative.
- [ ] The description explains in detail what this PR is about.
- [ ] I have linked a relevant issue or discussion.
- [ ] I have created tests covering the changes.
- [ ] I have updated the documentation and checked the documentation
preview.

### ⌛ Dependencies

<!-- List all open PRs that this PR logically depends on. For example,
-->
<!-- - #12345: short description why this is a dependency -->
<!-- - #34567: ... -->
    
URL: #38330
Reported by: Nathabolin
Reviewer(s):
  • Loading branch information
Release Manager committed Aug 9, 2024
2 parents 8d9da51 + ff83217 commit 66f679d
Showing 1 changed file with 200 additions and 2 deletions.
202 changes: 200 additions & 2 deletions src/sage/dynamics/arithmetic_dynamics/projective_ds.py
Original file line number Diff line number Diff line change
Expand Up @@ -108,6 +108,8 @@ class initialization directly.
from sage.schemes.projective.projective_space import ProjectiveSpace, ProjectiveSpace_ring
from sage.schemes.projective.projective_subscheme import AlgebraicScheme_subscheme_projective
from sage.structure.element import get_coercion_model
from sage.rings.qqbar import QQbar, number_field_elements_from_algebraics
from sage.schemes.elliptic_curves.constructor import EllipticCurve

lazy_import('sage.rings.algebraic_closure_finite_field', 'AlgebraicClosureFiniteField_generic')
lazy_import('sage.rings.number_field.number_field_ideal', 'NumberFieldFractionalIdeal')
Expand Down Expand Up @@ -6739,8 +6741,9 @@ def is_Lattes(self):
"""
# We need `f` to be defined over a number field for
# the function `is_postcrtically_finite` to work
if self.base_ring() not in NumberFields():
raise NotImplementedError("Base ring must be a number field")
if self.base_ring() is not QQbar:
if self.base_ring() not in NumberFields():
raise NotImplementedError("Base ring must be a number field")

if self.domain().dimension() != 1:
return False
Expand Down Expand Up @@ -6829,6 +6832,201 @@ def is_Lattes(self):
r_vals = sorted([val for val in r.values() if val != 1])
return r_vals in r_lattes_cases

def Lattes_to_curve(self, return_conjugation=False, check_lattes=False):
r"""
Finds a Short Weierstrass Model Elliptic curve of self
self assumed to be Lattes map and not in charateristic 2 or 3
INPUT:
`return_conjugation`` -- (default: ``False``) if ``True``, then
return the conjugation that moves self to a map that comes from a
Short Weierstrass Model Elliptic curve
`check_lattes``.-.(default:.``False``) if ``True``, then will ValueError if not Lattes
OUTPUT: a Short Weierstrass Model Elliptic curve which is isogenous to
the Elliptic curve of 'self',
If ``return_conjugation`` is ``True``
then also returns conjugation of 'self' to short form as a matrix
EXAMPLES::
sage: P.<x,y> = ProjectiveSpace(QQ, 1)
sage: f = P.Lattes_map(EllipticCurve([0, 0, 0, 10, 2]), 2)
sage: f.Lattes_to_curve()
Elliptic Curve defined by y^2 = x^3 + 10*x + 2 over Rational Field
::
sage: P.<x,y> = ProjectiveSpace(QQ, 1)
sage: M = matrix(QQ,2,2,[[1,2],[-1,2]])
sage: f = P.Lattes_map(EllipticCurve([1, 1, 1, 1, 2]), 2)
sage: f = f.conjugate(M)
sage: f.Lattes_to_curve(return_conjugation = True)
(
[ -7/36*a^2 + 7/12*a + 7/3 -17/18*a^2 + 17/6*a + 34/3]
[ -1/8*a^2 + 1/4*a + 3/2 1/4*a^2 - 1/2*a - 3],
Elliptic Curve defined by y^2 = x^3 + (-94/27*a^2+94/9*a+376/9)*x +
12232/243 over Number Field in a with defining polynomial y^3 - 18*y - 30
)
::
sage: P.<x,y> = ProjectiveSpace(QQ,1)
sage: f = P.Lattes_map(EllipticCurve([1, 1, 1, 2, 2]), 2)
sage: L.<i> = CyclotomicField(4)
sage: M = Matrix([[1+i,2*i], [0, -i]])
sage: f = f.conjugate(M)
sage: f.Lattes_to_curve(return_conjugation = True)
(
[ 1 19/24*a + 19/24]
[ 0 1],
Elliptic Curve defined by y^2 = x^3 + 95/96*a*x + (-1169/3456*a+1169/3456)
over Number Field in a with defining polynomial y^2 + 1
)
::
sage: P.<x,y> = ProjectiveSpace(QQ, 1)
sage: M = matrix(QQ,2,2,[[1,3],[2,1]])
sage: E = EllipticCurve([1, 1, 1, 2, 3])
sage: f = P.Lattes_map(E, 2)
sage: f = f.conjugate(M)
sage: f.Lattes_to_curve(return_conjugation = True)
(
[11/1602*a^2 41/3204*a^2]
[ -2/5*a -1/5*a],
Elliptic Curve defined by y^2 = x^3 + 2375/3421872*a^2*x + (-254125/61593696)
over Number Field in a with defining polynomial y^3 - 267
)
::
sage: P.<x,y> = ProjectiveSpace(QQ , 1)
sage: M = matrix(QQ,2,2,[[1 , 3],[2 , 1]])
sage: E = EllipticCurve([1, 1, 1, 2, 3])
sage: f = P.Lattes_map(E , 2)
sage: f = f.conjugate(M)
sage: m,H = f.Lattes_to_curve(true)
sage: J.<x,y> = ProjectiveSpace(H.base_ring(), 1)
sage: K = J.Lattes_map(H,2)
sage: K = K.conjugate(m)
sage: K.scale_by(f[0].lc()/K[0].lc())
sage: K == f.change_ring(K.base_ring())
True
::
sage: P.<x,y> = ProjectiveSpace(RR, 1)
sage: F = DynamicalSystem_projective([x^4, y^4])
sage: F.Lattes_to_curve(check_lattes=True)
Traceback (most recent call last):
...
NotImplementedError: Base ring must be a number field
::
sage: P.<x,y> = ProjectiveSpace(QQ, 1)
sage: F = DynamicalSystem_projective([x^4, y^4])
sage: F.Lattes_to_curve(check_lattes=True)
Traceback (most recent call last):
...
ValueError: Map is not Lattes
::
sage: P.<x,y> = ProjectiveSpace(QQ, 1)
sage: F = DynamicalSystem_projective([x^4, y^4])
sage: F.Lattes_to_curve()
Traceback (most recent call last):
...
ValueError: No Solutions found. Check if map is Lattes
::
sage: P.<x,y> = ProjectiveSpace(QQ, 1)
sage: F = DynamicalSystem_projective([x^3, y^3])
sage: F.Lattes_to_curve(check_lattes=True)
Traceback (most recent call last):
...
NotImplementedError: Map is not Lattes or is Complex Lattes
::
sage: K.<x>=QuadraticField(2)
sage: P.<a,y>=ProjectiveSpace(K, 1)
sage: E=EllipticCurve([1, x])
sage: f=P.Lattes_map(E, 2)
sage: f.Lattes_to_curve()
Elliptic Curve defined by y^2 = x^3 + x + a
over Number Field in a with defining polynomial y^2 - 2
::
sage: P.<x,y>=ProjectiveSpace(QQbar, 1)
sage: E=EllipticCurve([1, 2])
sage: f=P.Lattes_map(E, 2)
sage: f.Lattes_to_curve(check_lattes=true)
Elliptic Curve defined by y^2 = x^3 + x + 2 over Rational Field
"""
if self.base_ring() is not QQbar:
if self.base_ring() not in NumberFields():
raise NotImplementedError("Base ring must be a number field")
#The Complex case is hard to implement and needs to be done later
if sqrt(self.degree()) != int(sqrt(self.degree())):
raise NotImplementedError("Map is not Lattes or is Complex Lattes")
if check_lattes:
V = self.is_Lattes()
if not V:
raise ValueError("Map is not Lattes")
n = int(sqrt(self.degree()))
#Creating a Symbolic Lattes map f_sym from a short Elliptic curve
R = PolynomialRing(self.base_ring(), 6, "avar, bvar, uvar, vvar, wvar, tvar")
a, b, u, v, w, t = R.gens()
P = ProjectiveSpace(R, 1, self.domain().gens())
E_sym = EllipticCurve([a, b])
f_sym = P.Lattes_map(E_sym, n)
# Conjugating f_sym map to have the right form so we can solve for the conjugating matrix later
m = matrix(R, 2, [u, v, t, w])
f_sym = f_sym.conjugate(m)
f_sym.scale_by(u*w - v*t)
F_sym = f_sym.dehomogenize(1)
#extracting the base variables to do term by term matching
self.scale_by(1/self[0].lc())
F = self.dehomogenize(1)
#Creating a set of equations, eq, from term by term matching
eq = [u*w - v*t-1]
for j in range(2):
if j == 0:
g = F[0].numerator()
g_sym = F_sym[0].numerator()
else:
g = F[0].denominator()
g_sym = F_sym[0].denominator()
eq += (g - g_sym).coefficients()
#Solving the equations
phi = QQbar.coerce_map_from(R.base_ring())
if phi is None:
phi = R.base_ring().embeddings(QQbar)[0]
eq = [poly.numerator().change_ring(phi) for poly in eq]
I = eq[0].parent().ideal(eq)
pts = I.variety()
if len(pts) == 0:
raise ValueError("No Solutions found. Check if map is Lattes")
a = pts[0]['avar']
b = pts[0]['bvar']
u = pts[0]['uvar']
v = pts[0]['vvar']
t = pts[0]['tvar']
w = pts[0]['wvar']
K, [a, b, u, v, t, w], phi = number_field_elements_from_algebraics([a, b, u, v, t, w])
#creating our end products
E = EllipticCurve([a, b])
if return_conjugation:
M = matrix(K, 2, 2, [u, v, t, w])
return (M, E)
return E

class DynamicalSystem_projective_field(DynamicalSystem_projective,
SchemeMorphism_polynomial_projective_space_field):
Expand Down

0 comments on commit 66f679d

Please sign in to comment.