Period lattices of elliptic curves and related functions¶
Let \(E\) be an elliptic curve defined over a number field \(K\) (including \(\QQ\)). We attach a period lattice (a discrete rank 2 subgroup of \(\CC\)) to each embedding of \(K\) into \(\CC\).
In the case of real embeddings, the lattice is stable under complex conjugation and is called a real lattice. These have two types: rectangular, (the real curve has two connected components and positive discriminant) or non-rectangular (one connected component, negative discriminant).
The periods are computed to arbitrary precision using the AGM (Gauss’s Arithmetic-Geometric Mean).
EXAMPLES:
sage: x = polygen(ZZ, 'x')
sage: K.<a> = NumberField(x^3 - 2) # needs sage.rings.number_field
sage: E = EllipticCurve([0,1,0,a,a]) # needs sage.rings.number_field
>>> from sage.all import *
>>> x = polygen(ZZ, 'x')
>>> K = NumberField(x**Integer(3) - Integer(2), names=('a',)); (a,) = K._first_ngens(1)# needs sage.rings.number_field
>>> E = EllipticCurve([Integer(0),Integer(1),Integer(0),a,a]) # needs sage.rings.number_field
x = polygen(ZZ, 'x') K.<a> = NumberField(x^3 - 2) # needs sage.rings.number_field E = EllipticCurve([0,1,0,a,a]) # needs sage.rings.number_field
First we try a real embedding:
sage: emb = K.embeddings(RealField())[0] # needs sage.rings.number_field
sage: L = E.period_lattice(emb); L # needs sage.rings.number_field
Period lattice associated to Elliptic Curve defined by y^2 = x^3 + x^2 + a*x + a
over Number Field in a with defining polynomial x^3 - 2
with respect to the embedding Ring morphism:
From: Number Field in a with defining polynomial x^3 - 2
To: Algebraic Real Field
Defn: a |--> 1.259921049894873?
>>> from sage.all import *
>>> emb = K.embeddings(RealField())[Integer(0)] # needs sage.rings.number_field
>>> L = E.period_lattice(emb); L # needs sage.rings.number_field
Period lattice associated to Elliptic Curve defined by y^2 = x^3 + x^2 + a*x + a
over Number Field in a with defining polynomial x^3 - 2
with respect to the embedding Ring morphism:
From: Number Field in a with defining polynomial x^3 - 2
To: Algebraic Real Field
Defn: a |--> 1.259921049894873?
emb = K.embeddings(RealField())[0] # needs sage.rings.number_field L = E.period_lattice(emb); L # needs sage.rings.number_field
The first basis period is real:
sage: L.basis() # needs sage.rings.number_field
(3.81452977217855, 1.90726488608927 + 1.34047785962440*I)
sage: L.is_real() # needs sage.rings.number_field
True
>>> from sage.all import *
>>> L.basis() # needs sage.rings.number_field
(3.81452977217855, 1.90726488608927 + 1.34047785962440*I)
>>> L.is_real() # needs sage.rings.number_field
True
L.basis() # needs sage.rings.number_field L.is_real() # needs sage.rings.number_field
For a basis \(\omega_1,\omega_2\) normalised so that \(\omega_1/\omega_2\)
is in the fundamental region of the upper half-plane, use the method
normalised_basis()
instead:
sage: L.normalised_basis() # needs sage.rings.number_field
(1.90726488608927 - 1.34047785962440*I, -1.90726488608927 - 1.34047785962440*I)
>>> from sage.all import *
>>> L.normalised_basis() # needs sage.rings.number_field
(1.90726488608927 - 1.34047785962440*I, -1.90726488608927 - 1.34047785962440*I)
L.normalised_basis() # needs sage.rings.number_field
Next a complex embedding:
sage: emb = K.embeddings(ComplexField())[0] # needs sage.rings.number_field
sage: L = E.period_lattice(emb); L # needs sage.rings.number_field
Period lattice associated to Elliptic Curve defined by y^2 = x^3 + x^2 + a*x + a
over Number Field in a with defining polynomial x^3 - 2
with respect to the embedding Ring morphism:
From: Number Field in a with defining polynomial x^3 - 2
To: Algebraic Field
Defn: a |--> -0.6299605249474365? - 1.091123635971722?*I
>>> from sage.all import *
>>> emb = K.embeddings(ComplexField())[Integer(0)] # needs sage.rings.number_field
>>> L = E.period_lattice(emb); L # needs sage.rings.number_field
Period lattice associated to Elliptic Curve defined by y^2 = x^3 + x^2 + a*x + a
over Number Field in a with defining polynomial x^3 - 2
with respect to the embedding Ring morphism:
From: Number Field in a with defining polynomial x^3 - 2
To: Algebraic Field
Defn: a |--> -0.6299605249474365? - 1.091123635971722?*I
emb = K.embeddings(ComplexField())[0] # needs sage.rings.number_field L = E.period_lattice(emb); L # needs sage.rings.number_field
In this case, the basis \(\omega_1\), \(\omega_2\) is always normalised so that \(\tau = \omega_1/\omega_2\) is in the fundamental region in the upper half plane:
sage: # needs sage.rings.number_field
sage: w1, w2 = L.basis(); w1, w2
(-1.37588604166076 - 2.58560946624443*I, -2.10339907847356 + 0.428378776460622*I)
sage: L.is_real()
False
sage: tau = w1/w2; tau
0.387694505032876 + 1.30821088214407*I
sage: L.normalised_basis()
(-1.37588604166076 - 2.58560946624443*I, -2.10339907847356 + 0.428378776460622*I)
>>> from sage.all import *
>>> # needs sage.rings.number_field
>>> w1, w2 = L.basis(); w1, w2
(-1.37588604166076 - 2.58560946624443*I, -2.10339907847356 + 0.428378776460622*I)
>>> L.is_real()
False
>>> tau = w1/w2; tau
0.387694505032876 + 1.30821088214407*I
>>> L.normalised_basis()
(-1.37588604166076 - 2.58560946624443*I, -2.10339907847356 + 0.428378776460622*I)
# needs sage.rings.number_field w1, w2 = L.basis(); w1, w2 L.is_real() tau = w1/w2; tau L.normalised_basis()
We test that bug Issue #8415 (caused by a PARI bug fixed in v2.3.5) is OK:
sage: # needs sage.rings.number_field
sage: E = EllipticCurve('37a')
sage: K.<a> = QuadraticField(-7)
sage: EK = E.change_ring(K)
sage: EK.period_lattice(K.complex_embeddings()[0])
Period lattice associated to Elliptic Curve defined by y^2 + y = x^3 + (-1)*x
over Number Field in a with defining polynomial x^2 + 7
with a = 2.645751311064591?*I
with respect to the embedding Ring morphism:
From: Number Field in a with defining polynomial x^2 + 7
with a = 2.645751311064591?*I
To: Algebraic Field
Defn: a |--> -2.645751311064591?*I
>>> from sage.all import *
>>> # needs sage.rings.number_field
>>> E = EllipticCurve('37a')
>>> K = QuadraticField(-Integer(7), names=('a',)); (a,) = K._first_ngens(1)
>>> EK = E.change_ring(K)
>>> EK.period_lattice(K.complex_embeddings()[Integer(0)])
Period lattice associated to Elliptic Curve defined by y^2 + y = x^3 + (-1)*x
over Number Field in a with defining polynomial x^2 + 7
with a = 2.645751311064591?*I
with respect to the embedding Ring morphism:
From: Number Field in a with defining polynomial x^2 + 7
with a = 2.645751311064591?*I
To: Algebraic Field
Defn: a |--> -2.645751311064591?*I
# needs sage.rings.number_field E = EllipticCurve('37a') K.<a> = QuadraticField(-7) EK = E.change_ring(K) EK.period_lattice(K.complex_embeddings()[0])
REFERENCES:
AUTHORS:
?: initial version.
John Cremona:
Adapted to handle real embeddings of number fields, September 2008.
Added basis_matrix function, November 2008
Added support for complex embeddings, May 2009.
Added complex elliptic logs, March 2010; enhanced, October 2010.
- class sage.schemes.elliptic_curves.period_lattice.PeriodLattice(base_ring, rank, degree, sparse=False, coordinate_ring=None, category=None)[source]¶
Bases:
FreeModule_generic_pid
The class for the period lattice of an algebraic variety.
- class sage.schemes.elliptic_curves.period_lattice.PeriodLattice_ell(E, embedding=None)[source]¶
Bases:
PeriodLattice
The class for the period lattice of an elliptic curve.
Currently supported are elliptic curves defined over \(\QQ\), and elliptic curves defined over a number field with a real or complex embedding, where the lattice constructed depends on that embedding.
- basis(prec=None, algorithm='sage')[source]¶
Return a basis for this period lattice as a 2-tuple.
INPUT:
prec
– (default:None
) precision in bits (default precision ifNone
)algorithm
– string (default:'sage'
); choice of implementation (for real embeddings only) between'sage'
(native Sage implementation) or'pari'
(use the PARI library: only available for real embeddings)
OUTPUT:
(tuple of Complex) \((\omega_1,\omega_2)\) where the lattice is \(\ZZ\omega_1 + \ZZ\omega_2\). If the lattice is real then \(\omega_1\) is real and positive, \(\Im(\omega_2)>0\) and \(\Re(\omega_2/\omega_1)\) is either \(0\) (for rectangular lattices) or \(\frac{1}{2}\) (for non-rectangular lattices). Otherwise, \(\omega_1/\omega_2\) is in the fundamental region of the upper half-plane. If the latter normalisation is required for real lattices, use the method
normalised_basis()
instead.EXAMPLES:
sage: E = EllipticCurve('37a') sage: E.period_lattice().basis() (2.99345864623196, 2.45138938198679*I)
>>> from sage.all import * >>> E = EllipticCurve('37a') >>> E.period_lattice().basis() (2.99345864623196, 2.45138938198679*I)
E = EllipticCurve('37a') E.period_lattice().basis()
This shows that the issue reported at Issue #3954 is fixed:
sage: E = EllipticCurve('37a') sage: b1 = E.period_lattice().basis(prec=30) sage: b2 = E.period_lattice().basis(prec=30) sage: b1 == b2 True
>>> from sage.all import * >>> E = EllipticCurve('37a') >>> b1 = E.period_lattice().basis(prec=Integer(30)) >>> b2 = E.period_lattice().basis(prec=Integer(30)) >>> b1 == b2 True
E = EllipticCurve('37a') b1 = E.period_lattice().basis(prec=30) b2 = E.period_lattice().basis(prec=30) b1 == b2
This shows that the issue reported at Issue #4064 is fixed:
sage: E = EllipticCurve('37a') sage: E.period_lattice().basis(prec=30)[0].parent() Real Field with 30 bits of precision sage: E.period_lattice().basis(prec=100)[0].parent() Real Field with 100 bits of precision
>>> from sage.all import * >>> E = EllipticCurve('37a') >>> E.period_lattice().basis(prec=Integer(30))[Integer(0)].parent() Real Field with 30 bits of precision >>> E.period_lattice().basis(prec=Integer(100))[Integer(0)].parent() Real Field with 100 bits of precision
E = EllipticCurve('37a') E.period_lattice().basis(prec=30)[0].parent() E.period_lattice().basis(prec=100)[0].parent()
sage: # needs sage.rings.number_field sage: x = polygen(ZZ, 'x') sage: K.<a> = NumberField(x^3 - 2) sage: emb = K.embeddings(RealField())[0] sage: E = EllipticCurve([0,1,0,a,a]) sage: L = E.period_lattice(emb) sage: L.basis(64) (3.81452977217854509, 1.90726488608927255 + 1.34047785962440202*I) sage: # needs sage.rings.number_field sage: emb = K.embeddings(ComplexField())[0] sage: L = E.period_lattice(emb) sage: w1, w2 = L.basis(); w1, w2 (-1.37588604166076 - 2.58560946624443*I, -2.10339907847356 + 0.428378776460622*I) sage: L.is_real() False sage: tau = w1/w2; tau 0.387694505032876 + 1.30821088214407*I
>>> from sage.all import * >>> # needs sage.rings.number_field >>> x = polygen(ZZ, 'x') >>> K = NumberField(x**Integer(3) - Integer(2), names=('a',)); (a,) = K._first_ngens(1) >>> emb = K.embeddings(RealField())[Integer(0)] >>> E = EllipticCurve([Integer(0),Integer(1),Integer(0),a,a]) >>> L = E.period_lattice(emb) >>> L.basis(Integer(64)) (3.81452977217854509, 1.90726488608927255 + 1.34047785962440202*I) >>> # needs sage.rings.number_field >>> emb = K.embeddings(ComplexField())[Integer(0)] >>> L = E.period_lattice(emb) >>> w1, w2 = L.basis(); w1, w2 (-1.37588604166076 - 2.58560946624443*I, -2.10339907847356 + 0.428378776460622*I) >>> L.is_real() False >>> tau = w1/w2; tau 0.387694505032876 + 1.30821088214407*I
# needs sage.rings.number_field x = polygen(ZZ, 'x') K.<a> = NumberField(x^3 - 2) emb = K.embeddings(RealField())[0] E = EllipticCurve([0,1,0,a,a]) L = E.period_lattice(emb) L.basis(64) # needs sage.rings.number_field emb = K.embeddings(ComplexField())[0] L = E.period_lattice(emb) w1, w2 = L.basis(); w1, w2 L.is_real() tau = w1/w2; tau
>>> from sage.all import * >>> # needs sage.rings.number_field >>> x = polygen(ZZ, 'x') >>> K = NumberField(x**Integer(3) - Integer(2), names=('a',)); (a,) = K._first_ngens(1) >>> emb = K.embeddings(RealField())[Integer(0)] >>> E = EllipticCurve([Integer(0),Integer(1),Integer(0),a,a]) >>> L = E.period_lattice(emb) >>> L.basis(Integer(64)) (3.81452977217854509, 1.90726488608927255 + 1.34047785962440202*I) >>> # needs sage.rings.number_field >>> emb = K.embeddings(ComplexField())[Integer(0)] >>> L = E.period_lattice(emb) >>> w1, w2 = L.basis(); w1, w2 (-1.37588604166076 - 2.58560946624443*I, -2.10339907847356 + 0.428378776460622*I) >>> L.is_real() False >>> tau = w1/w2; tau 0.387694505032876 + 1.30821088214407*I
# needs sage.rings.number_field x = polygen(ZZ, 'x') K.<a> = NumberField(x^3 - 2) emb = K.embeddings(RealField())[0] E = EllipticCurve([0,1,0,a,a]) L = E.period_lattice(emb) L.basis(64) # needs sage.rings.number_field emb = K.embeddings(ComplexField())[0] L = E.period_lattice(emb) w1, w2 = L.basis(); w1, w2 L.is_real() tau = w1/w2; tau
- basis_matrix(prec=None, normalised=False)[source]¶
Return the basis matrix of this period lattice.
INPUT:
prec
– integer orNone
(default); real precision in bits (default real precision ifNone
)normalised
– boolean (default:False
); ifTrue
and the embedding is real, use the normalised basis (seenormalised_basis()
) instead of the default
OUTPUT:
A \(2\times 2\) real matrix whose rows are the lattice basis vectors, after identifying \(\CC\) with \(\RR^2\).
EXAMPLES:
sage: E = EllipticCurve('37a') sage: E.period_lattice().basis_matrix() [ 2.99345864623196 0.000000000000000] [0.000000000000000 2.45138938198679]
>>> from sage.all import * >>> E = EllipticCurve('37a') >>> E.period_lattice().basis_matrix() [ 2.99345864623196 0.000000000000000] [0.000000000000000 2.45138938198679]
E = EllipticCurve('37a') E.period_lattice().basis_matrix()
sage: # needs sage.rings.number_field sage: x = polygen(ZZ, 'x') sage: K.<a> = NumberField(x^3 - 2) sage: emb = K.embeddings(RealField())[0] sage: E = EllipticCurve([0,1,0,a,a]) sage: L = E.period_lattice(emb) sage: L.basis_matrix(64) [ 3.81452977217854509 0.000000000000000000] [ 1.90726488608927255 1.34047785962440202]
>>> from sage.all import * >>> # needs sage.rings.number_field >>> x = polygen(ZZ, 'x') >>> K = NumberField(x**Integer(3) - Integer(2), names=('a',)); (a,) = K._first_ngens(1) >>> emb = K.embeddings(RealField())[Integer(0)] >>> E = EllipticCurve([Integer(0),Integer(1),Integer(0),a,a]) >>> L = E.period_lattice(emb) >>> L.basis_matrix(Integer(64)) [ 3.81452977217854509 0.000000000000000000] [ 1.90726488608927255 1.34047785962440202]
# needs sage.rings.number_field x = polygen(ZZ, 'x') K.<a> = NumberField(x^3 - 2) emb = K.embeddings(RealField())[0] E = EllipticCurve([0,1,0,a,a]) L = E.period_lattice(emb) L.basis_matrix(64)
>>> from sage.all import * >>> # needs sage.rings.number_field >>> x = polygen(ZZ, 'x') >>> K = NumberField(x**Integer(3) - Integer(2), names=('a',)); (a,) = K._first_ngens(1) >>> emb = K.embeddings(RealField())[Integer(0)] >>> E = EllipticCurve([Integer(0),Integer(1),Integer(0),a,a]) >>> L = E.period_lattice(emb) >>> L.basis_matrix(Integer(64)) [ 3.81452977217854509 0.000000000000000000] [ 1.90726488608927255 1.34047785962440202]
# needs sage.rings.number_field x = polygen(ZZ, 'x') K.<a> = NumberField(x^3 - 2) emb = K.embeddings(RealField())[0] E = EllipticCurve([0,1,0,a,a]) L = E.period_lattice(emb) L.basis_matrix(64)
See Issue #4388:
sage: L = EllipticCurve('11a1').period_lattice() sage: L.basis_matrix() [ 1.26920930427955 0.000000000000000] [0.634604652139777 1.45881661693850] sage: L.basis_matrix(normalised=True) [0.634604652139777 -1.45881661693850] [-1.26920930427955 0.000000000000000]
>>> from sage.all import * >>> L = EllipticCurve('11a1').period_lattice() >>> L.basis_matrix() [ 1.26920930427955 0.000000000000000] [0.634604652139777 1.45881661693850] >>> L.basis_matrix(normalised=True) [0.634604652139777 -1.45881661693850] [-1.26920930427955 0.000000000000000]
L = EllipticCurve('11a1').period_lattice() L.basis_matrix() L.basis_matrix(normalised=True)
sage: L = EllipticCurve('389a1').period_lattice() sage: L.basis_matrix() [ 2.49021256085505 0.000000000000000] [0.000000000000000 1.97173770155165] sage: L.basis_matrix(normalised=True) [ 2.49021256085505 0.000000000000000] [0.000000000000000 -1.97173770155165]
>>> from sage.all import * >>> L = EllipticCurve('389a1').period_lattice() >>> L.basis_matrix() [ 2.49021256085505 0.000000000000000] [0.000000000000000 1.97173770155165] >>> L.basis_matrix(normalised=True) [ 2.49021256085505 0.000000000000000] [0.000000000000000 -1.97173770155165]
L = EllipticCurve('389a1').period_lattice() L.basis_matrix() L.basis_matrix(normalised=True)
>>> from sage.all import * >>> L = EllipticCurve('389a1').period_lattice() >>> L.basis_matrix() [ 2.49021256085505 0.000000000000000] [0.000000000000000 1.97173770155165] >>> L.basis_matrix(normalised=True) [ 2.49021256085505 0.000000000000000] [0.000000000000000 -1.97173770155165]
L = EllipticCurve('389a1').period_lattice() L.basis_matrix() L.basis_matrix(normalised=True)
- complex_area(prec=None)[source]¶
Return the area of a fundamental domain for the period lattice of the elliptic curve.
INPUT:
prec
– integer orNone
(default); real precision in bits (default real precision ifNone
)
EXAMPLES:
sage: E = EllipticCurve('37a') sage: E.period_lattice().complex_area() 7.33813274078958
>>> from sage.all import * >>> E = EllipticCurve('37a') >>> E.period_lattice().complex_area() 7.33813274078958
E = EllipticCurve('37a') E.period_lattice().complex_area()
sage: # needs sage.rings.number_field sage: x = polygen(ZZ, 'x') sage: K.<a> = NumberField(x^3 - 2) sage: embs = K.embeddings(ComplexField()) sage: E = EllipticCurve([0,1,0,a,a]) sage: [E.period_lattice(emb).is_real() for emb in K.embeddings(CC)] [False, False, True] sage: [E.period_lattice(emb).complex_area() for emb in embs] [6.02796894766694, 6.02796894766694, 5.11329270448345]
>>> from sage.all import * >>> # needs sage.rings.number_field >>> x = polygen(ZZ, 'x') >>> K = NumberField(x**Integer(3) - Integer(2), names=('a',)); (a,) = K._first_ngens(1) >>> embs = K.embeddings(ComplexField()) >>> E = EllipticCurve([Integer(0),Integer(1),Integer(0),a,a]) >>> [E.period_lattice(emb).is_real() for emb in K.embeddings(CC)] [False, False, True] >>> [E.period_lattice(emb).complex_area() for emb in embs] [6.02796894766694, 6.02796894766694, 5.11329270448345]
# needs sage.rings.number_field x = polygen(ZZ, 'x') K.<a> = NumberField(x^3 - 2) embs = K.embeddings(ComplexField()) E = EllipticCurve([0,1,0,a,a]) [E.period_lattice(emb).is_real() for emb in K.embeddings(CC)] [E.period_lattice(emb).complex_area() for emb in embs]
>>> from sage.all import * >>> # needs sage.rings.number_field >>> x = polygen(ZZ, 'x') >>> K = NumberField(x**Integer(3) - Integer(2), names=('a',)); (a,) = K._first_ngens(1) >>> embs = K.embeddings(ComplexField()) >>> E = EllipticCurve([Integer(0),Integer(1),Integer(0),a,a]) >>> [E.period_lattice(emb).is_real() for emb in K.embeddings(CC)] [False, False, True] >>> [E.period_lattice(emb).complex_area() for emb in embs] [6.02796894766694, 6.02796894766694, 5.11329270448345]
# needs sage.rings.number_field x = polygen(ZZ, 'x') K.<a> = NumberField(x^3 - 2) embs = K.embeddings(ComplexField()) E = EllipticCurve([0,1,0,a,a]) [E.period_lattice(emb).is_real() for emb in K.embeddings(CC)] [E.period_lattice(emb).complex_area() for emb in embs]
- coordinates(z, rounding=None)[source]¶
Return the coordinates of a complex number w.r.t. the lattice basis.
INPUT:
z
– complex numberrounding
– (default:None
) whether and how to round the output (see below)
OUTPUT:
When
rounding
isNone
(the default), returns a tuple of reals \(x\), \(y\) such that \(z=xw_1+yw_2\) where \(w_1\), \(w_2\) are a basis for the lattice (normalised in the case of complex embeddings).When
rounding
is'round'
, returns a tuple of integers \(n_1\), \(n_2\) which are the closest integers to the \(x\), \(y\) defined above. If \(z\) is in the lattice these are the coordinates of \(z\) with respect to the lattice basis.When
rounding
is'floor'
, returns a tuple of integers \(n_1\), \(n_2\) which are the integer parts to the \(x\), \(y\) defined above. These are used inreduce()
EXAMPLES:
sage: E = EllipticCurve('389a') sage: L = E.period_lattice() sage: w1, w2 = L.basis(prec=100) sage: P = E([-1,1]) sage: zP = P.elliptic_logarithm(precision=100); zP 0.47934825019021931612953301006 + 0.98586885077582410221120384908*I sage: L.coordinates(zP) (0.19249290511394227352563996419, 0.50000000000000000000000000000) sage: sum([x*w for x, w in zip(L.coordinates(zP), L.basis(prec=100))]) 0.47934825019021931612953301006 + 0.98586885077582410221120384908*I sage: L.coordinates(12*w1 + 23*w2) (12.000000000000000000000000000, 23.000000000000000000000000000) sage: L.coordinates(12*w1 + 23*w2, rounding='floor') (11, 22) sage: L.coordinates(12*w1 + 23*w2, rounding='round') (12, 23)
>>> from sage.all import * >>> E = EllipticCurve('389a') >>> L = E.period_lattice() >>> w1, w2 = L.basis(prec=Integer(100)) >>> P = E([-Integer(1),Integer(1)]) >>> zP = P.elliptic_logarithm(precision=Integer(100)); zP 0.47934825019021931612953301006 + 0.98586885077582410221120384908*I >>> L.coordinates(zP) (0.19249290511394227352563996419, 0.50000000000000000000000000000) >>> sum([x*w for x, w in zip(L.coordinates(zP), L.basis(prec=Integer(100)))]) 0.47934825019021931612953301006 + 0.98586885077582410221120384908*I >>> L.coordinates(Integer(12)*w1 + Integer(23)*w2) (12.000000000000000000000000000, 23.000000000000000000000000000) >>> L.coordinates(Integer(12)*w1 + Integer(23)*w2, rounding='floor') (11, 22) >>> L.coordinates(Integer(12)*w1 + Integer(23)*w2, rounding='round') (12, 23)
E = EllipticCurve('389a') L = E.period_lattice() w1, w2 = L.basis(prec=100) P = E([-1,1]) zP = P.elliptic_logarithm(precision=100); zP L.coordinates(zP) sum([x*w for x, w in zip(L.coordinates(zP), L.basis(prec=100))]) L.coordinates(12*w1 + 23*w2) L.coordinates(12*w1 + 23*w2, rounding='floor') L.coordinates(12*w1 + 23*w2, rounding='round')
- curve()[source]¶
Return the elliptic curve associated with this period lattice.
EXAMPLES:
sage: E = EllipticCurve('37a') sage: L = E.period_lattice() sage: L.curve() is E True
>>> from sage.all import * >>> E = EllipticCurve('37a') >>> L = E.period_lattice() >>> L.curve() is E True
E = EllipticCurve('37a') L = E.period_lattice() L.curve() is E
sage: # needs sage.rings.number_field sage: x = polygen(ZZ, 'x') sage: K.<a> = NumberField(x^3 - 2) sage: E = EllipticCurve([0,1,0,a,a]) sage: L = E.period_lattice(K.embeddings(RealField())[0]) sage: L.curve() is E True sage: L = E.period_lattice(K.embeddings(ComplexField())[0]) # needs sage.rings.number_field sage: L.curve() is E # needs sage.rings.number_field True
>>> from sage.all import * >>> # needs sage.rings.number_field >>> x = polygen(ZZ, 'x') >>> K = NumberField(x**Integer(3) - Integer(2), names=('a',)); (a,) = K._first_ngens(1) >>> E = EllipticCurve([Integer(0),Integer(1),Integer(0),a,a]) >>> L = E.period_lattice(K.embeddings(RealField())[Integer(0)]) >>> L.curve() is E True >>> L = E.period_lattice(K.embeddings(ComplexField())[Integer(0)]) # needs sage.rings.number_field >>> L.curve() is E # needs sage.rings.number_field True
# needs sage.rings.number_field x = polygen(ZZ, 'x') K.<a> = NumberField(x^3 - 2) E = EllipticCurve([0,1,0,a,a]) L = E.period_lattice(K.embeddings(RealField())[0]) L.curve() is E L = E.period_lattice(K.embeddings(ComplexField())[0]) # needs sage.rings.number_field L.curve() is E # needs sage.rings.number_field
>>> from sage.all import * >>> # needs sage.rings.number_field >>> x = polygen(ZZ, 'x') >>> K = NumberField(x**Integer(3) - Integer(2), names=('a',)); (a,) = K._first_ngens(1) >>> E = EllipticCurve([Integer(0),Integer(1),Integer(0),a,a]) >>> L = E.period_lattice(K.embeddings(RealField())[Integer(0)]) >>> L.curve() is E True >>> L = E.period_lattice(K.embeddings(ComplexField())[Integer(0)]) # needs sage.rings.number_field >>> L.curve() is E # needs sage.rings.number_field True
# needs sage.rings.number_field x = polygen(ZZ, 'x') K.<a> = NumberField(x^3 - 2) E = EllipticCurve([0,1,0,a,a]) L = E.period_lattice(K.embeddings(RealField())[0]) L.curve() is E L = E.period_lattice(K.embeddings(ComplexField())[0]) # needs sage.rings.number_field L.curve() is E # needs sage.rings.number_field
- e_log_RC(xP, yP, prec=None, reduce=True)[source]¶
Return the elliptic logarithm of a real or complex point.
INPUT:
xP
,yP
(real or complex) – Coordinates of a point on the embedded elliptic curve associated with this period lattice.prec
– (default:None
) real precision in bits (default real precision ifNone
)reduce
– boolean (default:True
); ifTrue
, the result is reduced with respect to the period lattice basis
OUTPUT:
(complex number) The elliptic logarithm of the point \((xP,yP)\) with respect to this period lattice. If \(E\) is the elliptic curve and \(\sigma:K\to\CC\) the embedding, the returned value \(z\) is such that \(z\pmod{L}\) maps to \((xP,yP)=\sigma(P)\) under the standard Weierstrass isomorphism from \(\CC/L\) to \(\sigma(E)\). If
reduce
isTrue
, the output is reduced so that it is in the fundamental period parallelogram with respect to the normalised lattice basis.ALGORITHM:
Uses the complex AGM. See [CT2013] for details.
EXAMPLES:
sage: E = EllipticCurve('389a') sage: L = E.period_lattice() sage: P = E([-1,1]) sage: xP, yP = [RR(c) for c in P.xy()]
>>> from sage.all import * >>> E = EllipticCurve('389a') >>> L = E.period_lattice() >>> P = E([-Integer(1),Integer(1)]) >>> xP, yP = [RR(c) for c in P.xy()]
E = EllipticCurve('389a') L = E.period_lattice() P = E([-1,1]) xP, yP = [RR(c) for c in P.xy()]
The elliptic log from the real coordinates:
sage: L.e_log_RC(xP, yP) 0.479348250190219 + 0.985868850775824*I
>>> from sage.all import * >>> L.e_log_RC(xP, yP) 0.479348250190219 + 0.985868850775824*I
L.e_log_RC(xP, yP)
The same elliptic log from the algebraic point:
sage: L(P) 0.479348250190219 + 0.985868850775824*I
>>> from sage.all import * >>> L(P) 0.479348250190219 + 0.985868850775824*I
L(P)
A number field example:
sage: # needs sage.rings.number_field sage: x = polygen(ZZ, 'x') sage: K.<a> = NumberField(x^3 - 2) sage: E = EllipticCurve([0,0,0,0,a]) sage: v = K.real_places()[0] sage: L = E.period_lattice(v) sage: P = E.lift_x(1/3*a^2 + a + 5/3) sage: L(P) 3.51086196882538 sage: xP, yP = [v(c) for c in P.xy()] sage: L.e_log_RC(xP, yP) 3.51086196882538
>>> from sage.all import * >>> # needs sage.rings.number_field >>> x = polygen(ZZ, 'x') >>> K = NumberField(x**Integer(3) - Integer(2), names=('a',)); (a,) = K._first_ngens(1) >>> E = EllipticCurve([Integer(0),Integer(0),Integer(0),Integer(0),a]) >>> v = K.real_places()[Integer(0)] >>> L = E.period_lattice(v) >>> P = E.lift_x(Integer(1)/Integer(3)*a**Integer(2) + a + Integer(5)/Integer(3)) >>> L(P) 3.51086196882538 >>> xP, yP = [v(c) for c in P.xy()] >>> L.e_log_RC(xP, yP) 3.51086196882538
# needs sage.rings.number_field x = polygen(ZZ, 'x') K.<a> = NumberField(x^3 - 2) E = EllipticCurve([0,0,0,0,a]) v = K.real_places()[0] L = E.period_lattice(v) P = E.lift_x(1/3*a^2 + a + 5/3) L(P) xP, yP = [v(c) for c in P.xy()] L.e_log_RC(xP, yP)
Elliptic logs of real points which do not come from algebraic points:
sage: # needs sage.rings.number_field sage: ER = EllipticCurve([v(ai) for ai in E.a_invariants()]) sage: P = ER.lift_x(12.34) sage: xP, yP = P.xy() sage: xP, yP (12.3400000000000, -43.3628968710567) sage: L.e_log_RC(xP, yP) 0.284656841192041 sage: xP, yP = ER.lift_x(0).xy() sage: L.e_log_RC(xP, yP) 1.34921304541057
>>> from sage.all import * >>> # needs sage.rings.number_field >>> ER = EllipticCurve([v(ai) for ai in E.a_invariants()]) >>> P = ER.lift_x(RealNumber('12.34')) >>> xP, yP = P.xy() >>> xP, yP (12.3400000000000, -43.3628968710567) >>> L.e_log_RC(xP, yP) 0.284656841192041 >>> xP, yP = ER.lift_x(Integer(0)).xy() >>> L.e_log_RC(xP, yP) 1.34921304541057
# needs sage.rings.number_field ER = EllipticCurve([v(ai) for ai in E.a_invariants()]) P = ER.lift_x(12.34) xP, yP = P.xy() xP, yP L.e_log_RC(xP, yP) xP, yP = ER.lift_x(0).xy() L.e_log_RC(xP, yP)
Elliptic logs of complex points:
sage: # needs sage.rings.number_field sage: v = K.complex_embeddings()[0] sage: L = E.period_lattice(v) sage: P = E.lift_x(1/3*a^2 + a + 5/3) sage: L(P) 1.68207104397706 - 1.87873661686704*I sage: xP, yP = [v(c) for c in P.xy()] sage: L.e_log_RC(xP, yP) 1.68207104397706 - 1.87873661686704*I sage: EC = EllipticCurve([v(ai) for ai in E.a_invariants()]) sage: xP, yP = EC.lift_x(0).xy() sage: L.e_log_RC(xP, yP) 2.06711431204080 - 1.73451485683471*I
>>> from sage.all import * >>> # needs sage.rings.number_field >>> v = K.complex_embeddings()[Integer(0)] >>> L = E.period_lattice(v) >>> P = E.lift_x(Integer(1)/Integer(3)*a**Integer(2) + a + Integer(5)/Integer(3)) >>> L(P) 1.68207104397706 - 1.87873661686704*I >>> xP, yP = [v(c) for c in P.xy()] >>> L.e_log_RC(xP, yP) 1.68207104397706 - 1.87873661686704*I >>> EC = EllipticCurve([v(ai) for ai in E.a_invariants()]) >>> xP, yP = EC.lift_x(Integer(0)).xy() >>> L.e_log_RC(xP, yP) 2.06711431204080 - 1.73451485683471*I
# needs sage.rings.number_field v = K.complex_embeddings()[0] L = E.period_lattice(v) P = E.lift_x(1/3*a^2 + a + 5/3) L(P) xP, yP = [v(c) for c in P.xy()] L.e_log_RC(xP, yP) EC = EllipticCurve([v(ai) for ai in E.a_invariants()]) xP, yP = EC.lift_x(0).xy() L.e_log_RC(xP, yP)
- ei()[source]¶
Return the x-coordinates of the 2-division points of the elliptic curve associated with this period lattice, as elements of
QQbar
.EXAMPLES:
sage: E = EllipticCurve('37a') sage: L = E.period_lattice() sage: L.ei() [-1.107159871688768?, 0.2695944364054446?, 0.8375654352833230?]
>>> from sage.all import * >>> E = EllipticCurve('37a') >>> L = E.period_lattice() >>> L.ei() [-1.107159871688768?, 0.2695944364054446?, 0.8375654352833230?]
E = EllipticCurve('37a') L = E.period_lattice() L.ei()
In the following example, we should have one purely real 2-division point coordinate, and two conjugate purely imaginary coordinates.
sage: # needs sage.rings.number_field sage: x = polygen(ZZ, 'x') sage: K.<a> = NumberField(x^3 - 2) sage: E = EllipticCurve([0,1,0,a,a]) sage: L = E.period_lattice(K.embeddings(RealField())[0]) sage: x1,x2,x3 = L.ei() sage: abs(x1.real()) + abs(x2.real()) < 1e-14 True sage: x1.imag(), x2.imag(), x3 (-1.122462048309373?, 1.122462048309373?, -1.000000000000000?)
>>> from sage.all import * >>> # needs sage.rings.number_field >>> x = polygen(ZZ, 'x') >>> K = NumberField(x**Integer(3) - Integer(2), names=('a',)); (a,) = K._first_ngens(1) >>> E = EllipticCurve([Integer(0),Integer(1),Integer(0),a,a]) >>> L = E.period_lattice(K.embeddings(RealField())[Integer(0)]) >>> x1,x2,x3 = L.ei() >>> abs(x1.real()) + abs(x2.real()) < RealNumber('1e-14') True >>> x1.imag(), x2.imag(), x3 (-1.122462048309373?, 1.122462048309373?, -1.000000000000000?)
# needs sage.rings.number_field x = polygen(ZZ, 'x') K.<a> = NumberField(x^3 - 2) E = EllipticCurve([0,1,0,a,a]) L = E.period_lattice(K.embeddings(RealField())[0]) x1,x2,x3 = L.ei() abs(x1.real()) + abs(x2.real()) < 1e-14 x1.imag(), x2.imag(), x3
sage: L = E.period_lattice(K.embeddings(ComplexField())[0]) # needs sage.rings.number_field sage: L.ei() # needs sage.rings.number_field [-1.000000000000000? + 0.?e-1...*I, -0.9720806486198328? - 0.561231024154687?*I, 0.9720806486198328? + 0.561231024154687?*I]
>>> from sage.all import * >>> L = E.period_lattice(K.embeddings(ComplexField())[Integer(0)]) # needs sage.rings.number_field >>> L.ei() # needs sage.rings.number_field [-1.000000000000000? + 0.?e-1...*I, -0.9720806486198328? - 0.561231024154687?*I, 0.9720806486198328? + 0.561231024154687?*I]
L = E.period_lattice(K.embeddings(ComplexField())[0]) # needs sage.rings.number_field L.ei() # needs sage.rings.number_field
>>> from sage.all import * >>> L = E.period_lattice(K.embeddings(ComplexField())[Integer(0)]) # needs sage.rings.number_field >>> L.ei() # needs sage.rings.number_field [-1.000000000000000? + 0.?e-1...*I, -0.9720806486198328? - 0.561231024154687?*I, 0.9720806486198328? + 0.561231024154687?*I]
L = E.period_lattice(K.embeddings(ComplexField())[0]) # needs sage.rings.number_field L.ei() # needs sage.rings.number_field
- elliptic_exponential(z, to_curve=True)[source]¶
Return the elliptic exponential of a complex number.
INPUT:
z
– complex number (viewed modulo this period lattice)to_curve
– boolean (default:True
); see below
OUTPUT:
If
to_curve
is False, a 2-tuple of real or complex numbers representing the point \((x,y) = (\wp(z),\wp'(z))\) where \(\wp\) denotes the Weierstrass \(\wp\)-function with respect to this lattice.If
to_curve
is True, the point \((X,Y) = (x-b_2/12,y-(a_1(x-b_2/12)-a_3)/2)\) as a point in \(E(\RR)\) or \(E(\CC)\), with \((x,y) = (\wp(z),\wp'(z))\) as above, where \(E\) is the elliptic curve over \(\RR\) or \(\CC\) whose period lattice this is.If the lattice is real and \(z\) is also real then the output is a pair of real numbers if
to_curve
isTrue
, or a point in \(E(\RR)\) ifto_curve
isFalse
.
Note
The precision is taken from that of the input
z
.EXAMPLES:
sage: E = EllipticCurve([1,1,1,-8,6]) sage: P = E(1, -2) sage: L = E.period_lattice() sage: z = L(P); z 1.17044757240090 sage: L.elliptic_exponential(z) (0.999999999999999 : -2.00000000000000 : 1.00000000000000) sage: _.curve() Elliptic Curve defined by y^2 + 1.00000000000000*x*y + 1.00000000000000*y = x^3 + 1.00000000000000*x^2 - 8.00000000000000*x + 6.00000000000000 over Real Field with 53 bits of precision sage: L.elliptic_exponential(z,to_curve=False) (1.41666666666667, -2.00000000000000) sage: z = L(P, prec=201); z 1.17044757240089592298992188482371493504472561677451007994189 sage: L.elliptic_exponential(z) (1.00000000000000000000000000000000000000000000000000000000000 : -2.00000000000000000000000000000000000000000000000000000000000 : 1.00000000000000000000000000000000000000000000000000000000000)
>>> from sage.all import * >>> E = EllipticCurve([Integer(1),Integer(1),Integer(1),-Integer(8),Integer(6)]) >>> P = E(Integer(1), -Integer(2)) >>> L = E.period_lattice() >>> z = L(P); z 1.17044757240090 >>> L.elliptic_exponential(z) (0.999999999999999 : -2.00000000000000 : 1.00000000000000) >>> _.curve() Elliptic Curve defined by y^2 + 1.00000000000000*x*y + 1.00000000000000*y = x^3 + 1.00000000000000*x^2 - 8.00000000000000*x + 6.00000000000000 over Real Field with 53 bits of precision >>> L.elliptic_exponential(z,to_curve=False) (1.41666666666667, -2.00000000000000) >>> z = L(P, prec=Integer(201)); z 1.17044757240089592298992188482371493504472561677451007994189 >>> L.elliptic_exponential(z) (1.00000000000000000000000000000000000000000000000000000000000 : -2.00000000000000000000000000000000000000000000000000000000000 : 1.00000000000000000000000000000000000000000000000000000000000)
E = EllipticCurve([1,1,1,-8,6]) P = E(1, -2) L = E.period_lattice() z = L(P); z L.elliptic_exponential(z) _.curve() L.elliptic_exponential(z,to_curve=False) z = L(P, prec=201); z L.elliptic_exponential(z)
Examples over number fields:
sage: # needs sage.rings.number_field sage: x = polygen(QQ) sage: K.<a> = NumberField(x^3 - 2) sage: embs = K.embeddings(CC) sage: E = EllipticCurve('37a') sage: EK = E.change_ring(K) sage: Li = [EK.period_lattice(e) for e in embs] sage: P = EK(-1, -1) sage: Q = EK(a - 1, 1 - a^2) sage: zi = [L.elliptic_logarithm(P) for L in Li] sage: [c.real() for c in Li[0].elliptic_exponential(zi[0])] [-1.00000000000000, -1.00000000000000, 1.00000000000000] sage: [c.real() for c in Li[0].elliptic_exponential(zi[1])] [-1.00000000000000, -1.00000000000000, 1.00000000000000] sage: [c.real() for c in Li[0].elliptic_exponential(zi[2])] [-1.00000000000000, -1.00000000000000, 1.00000000000000] sage: # needs sage.rings.number_field sage: zi = [L.elliptic_logarithm(Q) for L in Li] sage: Li[0].elliptic_exponential(zi[0]) (-1.62996052494744 - 1.09112363597172*I : 1.79370052598410 - 1.37472963699860*I : 1.00000000000000) sage: [embs[0](c) for c in Q] [-1.62996052494744 - 1.09112363597172*I, 1.79370052598410 - 1.37472963699860*I, 1.00000000000000] sage: Li[1].elliptic_exponential(zi[1]) (-1.62996052494744 + 1.09112363597172*I : 1.79370052598410 + 1.37472963699860*I : 1.00000000000000) sage: [embs[1](c) for c in Q] [-1.62996052494744 + 1.09112363597172*I, 1.79370052598410 + 1.37472963699860*I, 1.00000000000000] sage: [c.real() for c in Li[2].elliptic_exponential(zi[2])] [0.259921049894873, -0.587401051968199, 1.00000000000000] sage: [embs[2](c) for c in Q] [0.259921049894873, -0.587401051968200, 1.00000000000000]
>>> from sage.all import * >>> # needs sage.rings.number_field >>> x = polygen(QQ) >>> K = NumberField(x**Integer(3) - Integer(2), names=('a',)); (a,) = K._first_ngens(1) >>> embs = K.embeddings(CC) >>> E = EllipticCurve('37a') >>> EK = E.change_ring(K) >>> Li = [EK.period_lattice(e) for e in embs] >>> P = EK(-Integer(1), -Integer(1)) >>> Q = EK(a - Integer(1), Integer(1) - a**Integer(2)) >>> zi = [L.elliptic_logarithm(P) for L in Li] >>> [c.real() for c in Li[Integer(0)].elliptic_exponential(zi[Integer(0)])] [-1.00000000000000, -1.00000000000000, 1.00000000000000] >>> [c.real() for c in Li[Integer(0)].elliptic_exponential(zi[Integer(1)])] [-1.00000000000000, -1.00000000000000, 1.00000000000000] >>> [c.real() for c in Li[Integer(0)].elliptic_exponential(zi[Integer(2)])] [-1.00000000000000, -1.00000000000000, 1.00000000000000] >>> # needs sage.rings.number_field >>> zi = [L.elliptic_logarithm(Q) for L in Li] >>> Li[Integer(0)].elliptic_exponential(zi[Integer(0)]) (-1.62996052494744 - 1.09112363597172*I : 1.79370052598410 - 1.37472963699860*I : 1.00000000000000) >>> [embs[Integer(0)](c) for c in Q] [-1.62996052494744 - 1.09112363597172*I, 1.79370052598410 - 1.37472963699860*I, 1.00000000000000] >>> Li[Integer(1)].elliptic_exponential(zi[Integer(1)]) (-1.62996052494744 + 1.09112363597172*I : 1.79370052598410 + 1.37472963699860*I : 1.00000000000000) >>> [embs[Integer(1)](c) for c in Q] [-1.62996052494744 + 1.09112363597172*I, 1.79370052598410 + 1.37472963699860*I, 1.00000000000000] >>> [c.real() for c in Li[Integer(2)].elliptic_exponential(zi[Integer(2)])] [0.259921049894873, -0.587401051968199, 1.00000000000000] >>> [embs[Integer(2)](c) for c in Q] [0.259921049894873, -0.587401051968200, 1.00000000000000]
# needs sage.rings.number_field x = polygen(QQ) K.<a> = NumberField(x^3 - 2) embs = K.embeddings(CC) E = EllipticCurve('37a') EK = E.change_ring(K) Li = [EK.period_lattice(e) for e in embs] P = EK(-1, -1) Q = EK(a - 1, 1 - a^2) zi = [L.elliptic_logarithm(P) for L in Li] [c.real() for c in Li[0].elliptic_exponential(zi[0])] [c.real() for c in Li[0].elliptic_exponential(zi[1])] [c.real() for c in Li[0].elliptic_exponential(zi[2])] # needs sage.rings.number_field zi = [L.elliptic_logarithm(Q) for L in Li] Li[0].elliptic_exponential(zi[0]) [embs[0](c) for c in Q] Li[1].elliptic_exponential(zi[1]) [embs[1](c) for c in Q] [c.real() for c in Li[2].elliptic_exponential(zi[2])] [embs[2](c) for c in Q]
Test to show that Issue #8820 is fixed:
sage: # needs sage.rings.number_field sage: E = EllipticCurve('37a') sage: K.<a> = QuadraticField(-5) sage: L = E.change_ring(K).period_lattice(K.places()[0]) sage: L.elliptic_exponential(CDF(.1,.1)) (0.0000142854026029... - 49.9960001066650*I : 249.520141250950 + 250.019855549131*I : 1.00000000000000) sage: L.elliptic_exponential(CDF(.1,.1), to_curve=False) (0.0000142854026029447 - 49.9960001066650*I, 500.040282501900 + 500.039711098263*I)
>>> from sage.all import * >>> # needs sage.rings.number_field >>> E = EllipticCurve('37a') >>> K = QuadraticField(-Integer(5), names=('a',)); (a,) = K._first_ngens(1) >>> L = E.change_ring(K).period_lattice(K.places()[Integer(0)]) >>> L.elliptic_exponential(CDF(RealNumber('.1'),RealNumber('.1'))) (0.0000142854026029... - 49.9960001066650*I : 249.520141250950 + 250.019855549131*I : 1.00000000000000) >>> L.elliptic_exponential(CDF(RealNumber('.1'),RealNumber('.1')), to_curve=False) (0.0000142854026029447 - 49.9960001066650*I, 500.040282501900 + 500.039711098263*I)
# needs sage.rings.number_field E = EllipticCurve('37a') K.<a> = QuadraticField(-5) L = E.change_ring(K).period_lattice(K.places()[0]) L.elliptic_exponential(CDF(.1,.1)) L.elliptic_exponential(CDF(.1,.1), to_curve=False)
\(z=0\) is treated as a special case:
sage: E = EllipticCurve([1,1,1,-8,6]) sage: L = E.period_lattice() sage: L.elliptic_exponential(0) (0.000000000000000 : 1.00000000000000 : 0.000000000000000) sage: L.elliptic_exponential(0, to_curve=False) (+infinity, +infinity)
>>> from sage.all import * >>> E = EllipticCurve([Integer(1),Integer(1),Integer(1),-Integer(8),Integer(6)]) >>> L = E.period_lattice() >>> L.elliptic_exponential(Integer(0)) (0.000000000000000 : 1.00000000000000 : 0.000000000000000) >>> L.elliptic_exponential(Integer(0), to_curve=False) (+infinity, +infinity)
E = EllipticCurve([1,1,1,-8,6]) L = E.period_lattice() L.elliptic_exponential(0) L.elliptic_exponential(0, to_curve=False)
sage: # needs sage.rings.number_field sage: E = EllipticCurve('37a') sage: K.<a> = QuadraticField(-5) sage: L = E.change_ring(K).period_lattice(K.places()[0]) sage: P = L.elliptic_exponential(0); P (0.000000000000000 : 1.00000000000000 : 0.000000000000000) sage: P.parent() Abelian group of points on Elliptic Curve defined by y^2 + 1.00000000000000*y = x^3 + (-1.00000000000000)*x over Complex Field with 53 bits of precision
>>> from sage.all import * >>> # needs sage.rings.number_field >>> E = EllipticCurve('37a') >>> K = QuadraticField(-Integer(5), names=('a',)); (a,) = K._first_ngens(1) >>> L = E.change_ring(K).period_lattice(K.places()[Integer(0)]) >>> P = L.elliptic_exponential(Integer(0)); P (0.000000000000000 : 1.00000000000000 : 0.000000000000000) >>> P.parent() Abelian group of points on Elliptic Curve defined by y^2 + 1.00000000000000*y = x^3 + (-1.00000000000000)*x over Complex Field with 53 bits of precision
# needs sage.rings.number_field E = EllipticCurve('37a') K.<a> = QuadraticField(-5) L = E.change_ring(K).period_lattice(K.places()[0]) P = L.elliptic_exponential(0); P P.parent()
>>> from sage.all import * >>> # needs sage.rings.number_field >>> E = EllipticCurve('37a') >>> K = QuadraticField(-Integer(5), names=('a',)); (a,) = K._first_ngens(1) >>> L = E.change_ring(K).period_lattice(K.places()[Integer(0)]) >>> P = L.elliptic_exponential(Integer(0)); P (0.000000000000000 : 1.00000000000000 : 0.000000000000000) >>> P.parent() Abelian group of points on Elliptic Curve defined by y^2 + 1.00000000000000*y = x^3 + (-1.00000000000000)*x over Complex Field with 53 bits of precision
# needs sage.rings.number_field E = EllipticCurve('37a') K.<a> = QuadraticField(-5) L = E.change_ring(K).period_lattice(K.places()[0]) P = L.elliptic_exponential(0); P P.parent()
Very small \(z\) are handled properly (see Issue #8820):
sage: # needs sage.rings.number_field sage: K.<a> = QuadraticField(-1) sage: E = EllipticCurve([0,0,0,a,0]) sage: L = E.period_lattice(K.complex_embeddings()[0]) sage: L.elliptic_exponential(1e-100) (0.000000000000000 : 1.00000000000000 : 0.000000000000000)
>>> from sage.all import * >>> # needs sage.rings.number_field >>> K = QuadraticField(-Integer(1), names=('a',)); (a,) = K._first_ngens(1) >>> E = EllipticCurve([Integer(0),Integer(0),Integer(0),a,Integer(0)]) >>> L = E.period_lattice(K.complex_embeddings()[Integer(0)]) >>> L.elliptic_exponential(RealNumber('1e-100')) (0.000000000000000 : 1.00000000000000 : 0.000000000000000)
# needs sage.rings.number_field K.<a> = QuadraticField(-1) E = EllipticCurve([0,0,0,a,0]) L = E.period_lattice(K.complex_embeddings()[0]) L.elliptic_exponential(1e-100)
The elliptic exponential of \(z\) is returned as \((0 : 1 : 0)\) if the coordinates of \(z\) with respect to the period lattice are approximately integral:
sage: (100/log(2.0,10))/0.8 415.241011860920 sage: L.elliptic_exponential((RealField(415)(1e-100))).is_zero() # needs sage.rings.number_field True sage: L.elliptic_exponential((RealField(420)(1e-100))).is_zero() # needs sage.rings.number_field False
>>> from sage.all import * >>> (Integer(100)/log(RealNumber('2.0'),Integer(10)))/RealNumber('0.8') 415.241011860920 >>> L.elliptic_exponential((RealField(Integer(415))(RealNumber('1e-100')))).is_zero() # needs sage.rings.number_field True >>> L.elliptic_exponential((RealField(Integer(420))(RealNumber('1e-100')))).is_zero() # needs sage.rings.number_field False
(100/log(2.0,10))/0.8 L.elliptic_exponential((RealField(415)(1e-100))).is_zero() # needs sage.rings.number_field L.elliptic_exponential((RealField(420)(1e-100))).is_zero() # needs sage.rings.number_field
- elliptic_logarithm(P, prec=None, reduce=True)[source]¶
Return the elliptic logarithm of a point.
INPUT:
P
– point on the elliptic curve associated with this period latticeprec
– (default:None
) real precision in bits (default real precision ifNone
)reduce
– boolean (default:True
); ifTrue
, the result is reduced with respect to the period lattice basis
OUTPUT:
(complex number) The elliptic logarithm of the point \(P\) with respect to this period lattice. If \(E\) is the elliptic curve and \(\sigma:K\to\CC\) the embedding, the returned value \(z\) is such that \(z\pmod{L}\) maps to \(\sigma(P)\) under the standard Weierstrass isomorphism from \(\CC/L\) to \(\sigma(E)\). If
reduce
isTrue
, the output is reduced so that it is in the fundamental period parallelogram with respect to the normalised lattice basis.ALGORITHM:
Uses the complex AGM. See [CT2013] for details.
EXAMPLES:
sage: E = EllipticCurve('389a') sage: L = E.period_lattice() sage: E.discriminant() > 0 True sage: L.real_flag 1 sage: P = E([-1,1]) sage: P.is_on_identity_component () False sage: L.elliptic_logarithm(P, prec=96) 0.4793482501902193161295330101 + 0.9858688507758241022112038491*I sage: Q=E([3,5]) sage: Q.is_on_identity_component() True sage: L.elliptic_logarithm(Q, prec=96) 1.931128271542559442488585220
>>> from sage.all import * >>> E = EllipticCurve('389a') >>> L = E.period_lattice() >>> E.discriminant() > Integer(0) True >>> L.real_flag 1 >>> P = E([-Integer(1),Integer(1)]) >>> P.is_on_identity_component () False >>> L.elliptic_logarithm(P, prec=Integer(96)) 0.4793482501902193161295330101 + 0.9858688507758241022112038491*I >>> Q=E([Integer(3),Integer(5)]) >>> Q.is_on_identity_component() True >>> L.elliptic_logarithm(Q, prec=Integer(96)) 1.931128271542559442488585220
E = EllipticCurve('389a') L = E.period_lattice() E.discriminant() > 0 L.real_flag P = E([-1,1]) P.is_on_identity_component () L.elliptic_logarithm(P, prec=96) Q=E([3,5]) Q.is_on_identity_component() L.elliptic_logarithm(Q, prec=96)
Note that this is actually the inverse of the Weierstrass isomorphism:
sage: L.elliptic_exponential(_) # abs tol 1e-26 (3.000000000000000000000000000 : 5.000000000000000000000000000 : 1.000000000000000000000000000)
>>> from sage.all import * >>> L.elliptic_exponential(_) # abs tol 1e-26 (3.000000000000000000000000000 : 5.000000000000000000000000000 : 1.000000000000000000000000000)
L.elliptic_exponential(_) # abs tol 1e-26
An example with negative discriminant, and a torsion point:
sage: E = EllipticCurve('11a1') sage: L = E.period_lattice() sage: E.discriminant() < 0 True sage: L.real_flag -1 sage: P = E([16,-61]) sage: L.elliptic_logarithm(P) 0.253841860855911 sage: L.real_period() / L.elliptic_logarithm(P) 5.00000000000000
>>> from sage.all import * >>> E = EllipticCurve('11a1') >>> L = E.period_lattice() >>> E.discriminant() < Integer(0) True >>> L.real_flag -1 >>> P = E([Integer(16),-Integer(61)]) >>> L.elliptic_logarithm(P) 0.253841860855911 >>> L.real_period() / L.elliptic_logarithm(P) 5.00000000000000
E = EllipticCurve('11a1') L = E.period_lattice() E.discriminant() < 0 L.real_flag P = E([16,-61]) L.elliptic_logarithm(P) L.real_period() / L.elliptic_logarithm(P)
An example where precision is problematic:
sage: E = EllipticCurve([1, 0, 1, -85357462, 303528987048]) #18074g1 sage: P = E([4458713781401/835903744, -64466909836503771/24167649046528, 1]) sage: L = E.period_lattice() sage: L.ei() [5334.003952567705? - 1.964393150436?e-6*I, 5334.003952567705? + 1.964393150436?e-6*I, -10668.25790513541?] sage: L.elliptic_logarithm(P,prec=100) 0.27656204014107061464076203097
>>> from sage.all import * >>> E = EllipticCurve([Integer(1), Integer(0), Integer(1), -Integer(85357462), Integer(303528987048)]) #18074g1 >>> P = E([Integer(4458713781401)/Integer(835903744), -Integer(64466909836503771)/Integer(24167649046528), Integer(1)]) >>> L = E.period_lattice() >>> L.ei() [5334.003952567705? - 1.964393150436?e-6*I, 5334.003952567705? + 1.964393150436?e-6*I, -10668.25790513541?] >>> L.elliptic_logarithm(P,prec=Integer(100)) 0.27656204014107061464076203097
E = EllipticCurve([1, 0, 1, -85357462, 303528987048]) #18074g1 P = E([4458713781401/835903744, -64466909836503771/24167649046528, 1]) L = E.period_lattice() L.ei() L.elliptic_logarithm(P,prec=100)
Some complex examples, taken from the paper by Cremona and Thongjunthug:
sage: # needs sage.rings.number_field sage: K.<i> = QuadraticField(-1) sage: a4 = 9*i - 10 sage: a6 = 21 - i sage: E = EllipticCurve([0,0,0,a4,a6]) sage: e1 = 3 - 2*i; e2 = 1 + i; e3 = -4 + i sage: emb = K.embeddings(CC)[1] sage: L = E.period_lattice(emb) sage: P = E(2 - i, 4 + 2*i)
>>> from sage.all import * >>> # needs sage.rings.number_field >>> K = QuadraticField(-Integer(1), names=('i',)); (i,) = K._first_ngens(1) >>> a4 = Integer(9)*i - Integer(10) >>> a6 = Integer(21) - i >>> E = EllipticCurve([Integer(0),Integer(0),Integer(0),a4,a6]) >>> e1 = Integer(3) - Integer(2)*i; e2 = Integer(1) + i; e3 = -Integer(4) + i >>> emb = K.embeddings(CC)[Integer(1)] >>> L = E.period_lattice(emb) >>> P = E(Integer(2) - i, Integer(4) + Integer(2)*i)
# needs sage.rings.number_field K.<i> = QuadraticField(-1) a4 = 9*i - 10 a6 = 21 - i E = EllipticCurve([0,0,0,a4,a6]) e1 = 3 - 2*i; e2 = 1 + i; e3 = -4 + i emb = K.embeddings(CC)[1] L = E.period_lattice(emb) P = E(2 - i, 4 + 2*i)
By default, the output is reduced with respect to the normalised lattice basis, so that its coordinates with respect to that basis lie in the interval [0,1):
sage: z = L.elliptic_logarithm(P, prec=100); z # needs sage.rings.number_field 0.70448375537782208460499649302 - 0.79246725643650979858266018068*I sage: L.coordinates(z) # needs sage.rings.number_field (0.46247636364807931766105406092, 0.79497588726808704200760395829)
>>> from sage.all import * >>> z = L.elliptic_logarithm(P, prec=Integer(100)); z # needs sage.rings.number_field 0.70448375537782208460499649302 - 0.79246725643650979858266018068*I >>> L.coordinates(z) # needs sage.rings.number_field (0.46247636364807931766105406092, 0.79497588726808704200760395829)
z = L.elliptic_logarithm(P, prec=100); z # needs sage.rings.number_field L.coordinates(z) # needs sage.rings.number_field
Using
reduce=False
this step can be omitted. In this case the coordinates are usually in the interval [-0.5,0.5), but this is not guaranteed. This option is mainly for testing purposes:sage: z = L.elliptic_logarithm(P, prec=100, reduce=False); z # needs sage.rings.number_field 0.57002153834710752778063503023 + 0.46476340520469798857457031393*I sage: L.coordinates(z) # needs sage.rings.number_field (0.46247636364807931766105406092, -0.20502411273191295799239604171)
>>> from sage.all import * >>> z = L.elliptic_logarithm(P, prec=Integer(100), reduce=False); z # needs sage.rings.number_field 0.57002153834710752778063503023 + 0.46476340520469798857457031393*I >>> L.coordinates(z) # needs sage.rings.number_field (0.46247636364807931766105406092, -0.20502411273191295799239604171)
z = L.elliptic_logarithm(P, prec=100, reduce=False); z # needs sage.rings.number_field L.coordinates(z) # needs sage.rings.number_field
The elliptic logs of the 2-torsion points are half-periods:
sage: L.elliptic_logarithm(E(e1, 0), prec=100) # needs sage.rings.number_field 0.64607575874356525952487867052 + 0.22379609053909448304176885364*I sage: L.elliptic_logarithm(E(e2, 0), prec=100) # needs sage.rings.number_field 0.71330686725892253793705940192 - 0.40481924028150941053684639367*I sage: L.elliptic_logarithm(E(e3, 0), prec=100) # needs sage.rings.number_field 0.067231108515357278412180731396 - 0.62861533082060389357861524731*I
>>> from sage.all import * >>> L.elliptic_logarithm(E(e1, Integer(0)), prec=Integer(100)) # needs sage.rings.number_field 0.64607575874356525952487867052 + 0.22379609053909448304176885364*I >>> L.elliptic_logarithm(E(e2, Integer(0)), prec=Integer(100)) # needs sage.rings.number_field 0.71330686725892253793705940192 - 0.40481924028150941053684639367*I >>> L.elliptic_logarithm(E(e3, Integer(0)), prec=Integer(100)) # needs sage.rings.number_field 0.067231108515357278412180731396 - 0.62861533082060389357861524731*I
L.elliptic_logarithm(E(e1, 0), prec=100) # needs sage.rings.number_field L.elliptic_logarithm(E(e2, 0), prec=100) # needs sage.rings.number_field L.elliptic_logarithm(E(e3, 0), prec=100) # needs sage.rings.number_field
We check this by doubling and seeing that the resulting coordinates are integers:
sage: L.coordinates(2*L.elliptic_logarithm(E(e1, 0), prec=100)) # needs sage.rings.number_field (1.0000000000000000000000000000, 0.00000000000000000000000000000) sage: L.coordinates(2*L.elliptic_logarithm(E(e2, 0), prec=100)) # needs sage.rings.number_field (1.0000000000000000000000000000, 1.0000000000000000000000000000) sage: L.coordinates(2*L.elliptic_logarithm(E(e3, 0), prec=100)) # needs sage.rings.number_field (0.00000000000000000000000000000, 1.0000000000000000000000000000)
>>> from sage.all import * >>> L.coordinates(Integer(2)*L.elliptic_logarithm(E(e1, Integer(0)), prec=Integer(100))) # needs sage.rings.number_field (1.0000000000000000000000000000, 0.00000000000000000000000000000) >>> L.coordinates(Integer(2)*L.elliptic_logarithm(E(e2, Integer(0)), prec=Integer(100))) # needs sage.rings.number_field (1.0000000000000000000000000000, 1.0000000000000000000000000000) >>> L.coordinates(Integer(2)*L.elliptic_logarithm(E(e3, Integer(0)), prec=Integer(100))) # needs sage.rings.number_field (0.00000000000000000000000000000, 1.0000000000000000000000000000)
L.coordinates(2*L.elliptic_logarithm(E(e1, 0), prec=100)) # needs sage.rings.number_field L.coordinates(2*L.elliptic_logarithm(E(e2, 0), prec=100)) # needs sage.rings.number_field L.coordinates(2*L.elliptic_logarithm(E(e3, 0), prec=100)) # needs sage.rings.number_field
sage: # needs sage.rings.number_field sage: a4 = -78*i + 104 sage: a6 = -216*i - 312 sage: E = EllipticCurve([0,0,0,a4,a6]) sage: emb = K.embeddings(CC)[1] sage: L = E.period_lattice(emb) sage: P = E(3 + 2*i, 14 - 7*i) sage: L.elliptic_logarithm(P) 0.297147783912228 - 0.546125549639461*I sage: L.coordinates(L.elliptic_logarithm(P)) (0.628653378040238, 0.371417754610223) sage: e1 = 1 + 3*i; e2 = -4 - 12*i; e3 = -e1 - e2 sage: L.coordinates(L.elliptic_logarithm(E(e1, 0))) (0.500000000000000, 0.500000000000000) sage: L.coordinates(L.elliptic_logarithm(E(e2, 0))) (1.00000000000000, 0.500000000000000) sage: L.coordinates(L.elliptic_logarithm(E(e3, 0))) (0.500000000000000, 0.000000000000000)
>>> from sage.all import * >>> # needs sage.rings.number_field >>> a4 = -Integer(78)*i + Integer(104) >>> a6 = -Integer(216)*i - Integer(312) >>> E = EllipticCurve([Integer(0),Integer(0),Integer(0),a4,a6]) >>> emb = K.embeddings(CC)[Integer(1)] >>> L = E.period_lattice(emb) >>> P = E(Integer(3) + Integer(2)*i, Integer(14) - Integer(7)*i) >>> L.elliptic_logarithm(P) 0.297147783912228 - 0.546125549639461*I >>> L.coordinates(L.elliptic_logarithm(P)) (0.628653378040238, 0.371417754610223) >>> e1 = Integer(1) + Integer(3)*i; e2 = -Integer(4) - Integer(12)*i; e3 = -e1 - e2 >>> L.coordinates(L.elliptic_logarithm(E(e1, Integer(0)))) (0.500000000000000, 0.500000000000000) >>> L.coordinates(L.elliptic_logarithm(E(e2, Integer(0)))) (1.00000000000000, 0.500000000000000) >>> L.coordinates(L.elliptic_logarithm(E(e3, Integer(0)))) (0.500000000000000, 0.000000000000000)
# needs sage.rings.number_field a4 = -78*i + 104 a6 = -216*i - 312 E = EllipticCurve([0,0,0,a4,a6]) emb = K.embeddings(CC)[1] L = E.period_lattice(emb) P = E(3 + 2*i, 14 - 7*i) L.elliptic_logarithm(P) L.coordinates(L.elliptic_logarithm(P)) e1 = 1 + 3*i; e2 = -4 - 12*i; e3 = -e1 - e2 L.coordinates(L.elliptic_logarithm(E(e1, 0))) L.coordinates(L.elliptic_logarithm(E(e2, 0))) L.coordinates(L.elliptic_logarithm(E(e3, 0)))
>>> from sage.all import * >>> # needs sage.rings.number_field >>> a4 = -Integer(78)*i + Integer(104) >>> a6 = -Integer(216)*i - Integer(312) >>> E = EllipticCurve([Integer(0),Integer(0),Integer(0),a4,a6]) >>> emb = K.embeddings(CC)[Integer(1)] >>> L = E.period_lattice(emb) >>> P = E(Integer(3) + Integer(2)*i, Integer(14) - Integer(7)*i) >>> L.elliptic_logarithm(P) 0.297147783912228 - 0.546125549639461*I >>> L.coordinates(L.elliptic_logarithm(P)) (0.628653378040238, 0.371417754610223) >>> e1 = Integer(1) + Integer(3)*i; e2 = -Integer(4) - Integer(12)*i; e3 = -e1 - e2 >>> L.coordinates(L.elliptic_logarithm(E(e1, Integer(0)))) (0.500000000000000, 0.500000000000000) >>> L.coordinates(L.elliptic_logarithm(E(e2, Integer(0)))) (1.00000000000000, 0.500000000000000) >>> L.coordinates(L.elliptic_logarithm(E(e3, Integer(0)))) (0.500000000000000, 0.000000000000000)
# needs sage.rings.number_field a4 = -78*i + 104 a6 = -216*i - 312 E = EllipticCurve([0,0,0,a4,a6]) emb = K.embeddings(CC)[1] L = E.period_lattice(emb) P = E(3 + 2*i, 14 - 7*i) L.elliptic_logarithm(P) L.coordinates(L.elliptic_logarithm(P)) e1 = 1 + 3*i; e2 = -4 - 12*i; e3 = -e1 - e2 L.coordinates(L.elliptic_logarithm(E(e1, 0))) L.coordinates(L.elliptic_logarithm(E(e2, 0))) L.coordinates(L.elliptic_logarithm(E(e3, 0)))
- gens(prec=None, algorithm='sage')[source]¶
Return a basis for this period lattice as a 2-tuple.
This is an alias for
basis()
. See the docstring there for a more in-depth explanation and further examples.INPUT:
prec
– (default:None
) precision in bits (default precision ifNone
)algorithm
– string (default:'sage'
); choice of implementation (for real embeddings only) between'sage'
(native Sage implementation) or'pari'
(use the PARI library: only available for real embeddings)
OUTPUT:
(tuple of Complex) \((\omega_1,\omega_2)\) where the lattice is \(\ZZ\omega_1 + \ZZ\omega_2\). If the lattice is real then \(\omega_1\) is real and positive, \(\Im(\omega_2)>0\) and \(\Re(\omega_2/\omega_1)\) is either \(0\) (for rectangular lattices) or \(\frac{1}{2}\) (for non-rectangular lattices). Otherwise, \(\omega_1/\omega_2\) is in the fundamental region of the upper half-plane. If the latter normalisation is required for real lattices, use the method
normalised_basis()
instead.EXAMPLES:
sage: E = EllipticCurve('37a') sage: E.period_lattice().gens() (2.99345864623196, 2.45138938198679*I) sage: E.period_lattice().gens(prec=100) (2.9934586462319596298320099794, 2.4513893819867900608542248319*I)
>>> from sage.all import * >>> E = EllipticCurve('37a') >>> E.period_lattice().gens() (2.99345864623196, 2.45138938198679*I) >>> E.period_lattice().gens(prec=Integer(100)) (2.9934586462319596298320099794, 2.4513893819867900608542248319*I)
E = EllipticCurve('37a') E.period_lattice().gens() E.period_lattice().gens(prec=100)
- is_real()[source]¶
Return
True
if this period lattice is real.EXAMPLES:
sage: f = EllipticCurve('11a') sage: f.period_lattice().is_real() True
>>> from sage.all import * >>> f = EllipticCurve('11a') >>> f.period_lattice().is_real() True
f = EllipticCurve('11a') f.period_lattice().is_real()
sage: # needs sage.rings.number_field sage: K.<i> = QuadraticField(-1) sage: E = EllipticCurve(K, [0,0,0,i,2*i]) sage: emb = K.embeddings(ComplexField())[0] sage: L = E.period_lattice(emb) sage: L.is_real() False
>>> from sage.all import * >>> # needs sage.rings.number_field >>> K = QuadraticField(-Integer(1), names=('i',)); (i,) = K._first_ngens(1) >>> E = EllipticCurve(K, [Integer(0),Integer(0),Integer(0),i,Integer(2)*i]) >>> emb = K.embeddings(ComplexField())[Integer(0)] >>> L = E.period_lattice(emb) >>> L.is_real() False
# needs sage.rings.number_field K.<i> = QuadraticField(-1) E = EllipticCurve(K, [0,0,0,i,2*i]) emb = K.embeddings(ComplexField())[0] L = E.period_lattice(emb) L.is_real()
>>> from sage.all import * >>> # needs sage.rings.number_field >>> K = QuadraticField(-Integer(1), names=('i',)); (i,) = K._first_ngens(1) >>> E = EllipticCurve(K, [Integer(0),Integer(0),Integer(0),i,Integer(2)*i]) >>> emb = K.embeddings(ComplexField())[Integer(0)] >>> L = E.period_lattice(emb) >>> L.is_real() False
# needs sage.rings.number_field K.<i> = QuadraticField(-1) E = EllipticCurve(K, [0,0,0,i,2*i]) emb = K.embeddings(ComplexField())[0] L = E.period_lattice(emb) L.is_real()
sage: x = polygen(ZZ, 'x') sage: K.<a> = NumberField(x^3 - 2) # needs sage.rings.number_field sage: E = EllipticCurve([0,1,0,a,a]) # needs sage.rings.number_field sage: [E.period_lattice(emb).is_real() for emb in K.embeddings(CC)] # needs sage.rings.number_field [False, False, True]
>>> from sage.all import * >>> x = polygen(ZZ, 'x') >>> K = NumberField(x**Integer(3) - Integer(2), names=('a',)); (a,) = K._first_ngens(1)# needs sage.rings.number_field >>> E = EllipticCurve([Integer(0),Integer(1),Integer(0),a,a]) # needs sage.rings.number_field >>> [E.period_lattice(emb).is_real() for emb in K.embeddings(CC)] # needs sage.rings.number_field [False, False, True]
x = polygen(ZZ, 'x') K.<a> = NumberField(x^3 - 2) # needs sage.rings.number_field E = EllipticCurve([0,1,0,a,a]) # needs sage.rings.number_field [E.period_lattice(emb).is_real() for emb in K.embeddings(CC)] # needs sage.rings.number_field
>>> from sage.all import * >>> x = polygen(ZZ, 'x') >>> K = NumberField(x**Integer(3) - Integer(2), names=('a',)); (a,) = K._first_ngens(1)# needs sage.rings.number_field >>> E = EllipticCurve([Integer(0),Integer(1),Integer(0),a,a]) # needs sage.rings.number_field >>> [E.period_lattice(emb).is_real() for emb in K.embeddings(CC)] # needs sage.rings.number_field [False, False, True]
x = polygen(ZZ, 'x') K.<a> = NumberField(x^3 - 2) # needs sage.rings.number_field E = EllipticCurve([0,1,0,a,a]) # needs sage.rings.number_field [E.period_lattice(emb).is_real() for emb in K.embeddings(CC)] # needs sage.rings.number_field
ALGORITHM:
The lattice is real if it is associated to a real embedding; such lattices are stable under conjugation.
- is_rectangular()[source]¶
Return
True
if this period lattice is rectangular.Note
Only defined for real lattices; a
RuntimeError
is raised for non-real lattices.EXAMPLES:
sage: f = EllipticCurve('11a') sage: f.period_lattice().basis() (1.26920930427955, 0.634604652139777 + 1.45881661693850*I) sage: f.period_lattice().is_rectangular() False
>>> from sage.all import * >>> f = EllipticCurve('11a') >>> f.period_lattice().basis() (1.26920930427955, 0.634604652139777 + 1.45881661693850*I) >>> f.period_lattice().is_rectangular() False
f = EllipticCurve('11a') f.period_lattice().basis() f.period_lattice().is_rectangular()
sage: f = EllipticCurve('37b') sage: f.period_lattice().basis() (1.08852159290423, 1.76761067023379*I) sage: f.period_lattice().is_rectangular() True
>>> from sage.all import * >>> f = EllipticCurve('37b') >>> f.period_lattice().basis() (1.08852159290423, 1.76761067023379*I) >>> f.period_lattice().is_rectangular() True
f = EllipticCurve('37b') f.period_lattice().basis() f.period_lattice().is_rectangular()
>>> from sage.all import * >>> f = EllipticCurve('37b') >>> f.period_lattice().basis() (1.08852159290423, 1.76761067023379*I) >>> f.period_lattice().is_rectangular() True
f = EllipticCurve('37b') f.period_lattice().basis() f.period_lattice().is_rectangular()
ALGORITHM:
The period lattice is rectangular precisely if the discriminant of the Weierstrass equation is positive, or equivalently if the number of real components is 2.
- normalised_basis(prec=None, algorithm='sage')[source]¶
Return a normalised basis for this period lattice as a 2-tuple.
INPUT:
prec
– (default:None
) precision in bits (default precision ifNone
)algorithm
– string (default:'sage'
); choice of implementation (for real embeddings only) between'sage'
(native Sage implementation) or'pari'
(use the PARI library: only available for real embeddings)
OUTPUT:
(tuple of Complex) \((\omega_1,\omega_2)\) where the lattice has the form \(\ZZ\omega_1 + \ZZ\omega_2\). The basis is normalised so that \(\omega_1/\omega_2\) is in the fundamental region of the upper half-plane. For an alternative normalisation for real lattices (with the first period real), use the method
basis()
instead.EXAMPLES:
sage: E = EllipticCurve('37a') sage: E.period_lattice().normalised_basis() (2.99345864623196, -2.45138938198679*I)
>>> from sage.all import * >>> E = EllipticCurve('37a') >>> E.period_lattice().normalised_basis() (2.99345864623196, -2.45138938198679*I)
E = EllipticCurve('37a') E.period_lattice().normalised_basis()
sage: # needs sage.rings.number_field sage: x = polygen(ZZ, 'x') sage: K.<a> = NumberField(x^3 - 2) sage: emb = K.embeddings(RealField())[0] sage: E = EllipticCurve([0,1,0,a,a]) sage: L = E.period_lattice(emb) sage: L.normalised_basis(64) (1.90726488608927255 - 1.34047785962440202*I, -1.90726488608927255 - 1.34047785962440202*I) sage: # needs sage.rings.number_field sage: emb = K.embeddings(ComplexField())[0] sage: L = E.period_lattice(emb) sage: w1, w2 = L.normalised_basis(); w1, w2 (-1.37588604166076 - 2.58560946624443*I, -2.10339907847356 + 0.428378776460622*I) sage: L.is_real() False sage: tau = w1/w2; tau 0.387694505032876 + 1.30821088214407*I
>>> from sage.all import * >>> # needs sage.rings.number_field >>> x = polygen(ZZ, 'x') >>> K = NumberField(x**Integer(3) - Integer(2), names=('a',)); (a,) = K._first_ngens(1) >>> emb = K.embeddings(RealField())[Integer(0)] >>> E = EllipticCurve([Integer(0),Integer(1),Integer(0),a,a]) >>> L = E.period_lattice(emb) >>> L.normalised_basis(Integer(64)) (1.90726488608927255 - 1.34047785962440202*I, -1.90726488608927255 - 1.34047785962440202*I) >>> # needs sage.rings.number_field >>> emb = K.embeddings(ComplexField())[Integer(0)] >>> L = E.period_lattice(emb) >>> w1, w2 = L.normalised_basis(); w1, w2 (-1.37588604166076 - 2.58560946624443*I, -2.10339907847356 + 0.428378776460622*I) >>> L.is_real() False >>> tau = w1/w2; tau 0.387694505032876 + 1.30821088214407*I
# needs sage.rings.number_field x = polygen(ZZ, 'x') K.<a> = NumberField(x^3 - 2) emb = K.embeddings(RealField())[0] E = EllipticCurve([0,1,0,a,a]) L = E.period_lattice(emb) L.normalised_basis(64) # needs sage.rings.number_field emb = K.embeddings(ComplexField())[0] L = E.period_lattice(emb) w1, w2 = L.normalised_basis(); w1, w2 L.is_real() tau = w1/w2; tau
>>> from sage.all import * >>> # needs sage.rings.number_field >>> x = polygen(ZZ, 'x') >>> K = NumberField(x**Integer(3) - Integer(2), names=('a',)); (a,) = K._first_ngens(1) >>> emb = K.embeddings(RealField())[Integer(0)] >>> E = EllipticCurve([Integer(0),Integer(1),Integer(0),a,a]) >>> L = E.period_lattice(emb) >>> L.normalised_basis(Integer(64)) (1.90726488608927255 - 1.34047785962440202*I, -1.90726488608927255 - 1.34047785962440202*I) >>> # needs sage.rings.number_field >>> emb = K.embeddings(ComplexField())[Integer(0)] >>> L = E.period_lattice(emb) >>> w1, w2 = L.normalised_basis(); w1, w2 (-1.37588604166076 - 2.58560946624443*I, -2.10339907847356 + 0.428378776460622*I) >>> L.is_real() False >>> tau = w1/w2; tau 0.387694505032876 + 1.30821088214407*I
# needs sage.rings.number_field x = polygen(ZZ, 'x') K.<a> = NumberField(x^3 - 2) emb = K.embeddings(RealField())[0] E = EllipticCurve([0,1,0,a,a]) L = E.period_lattice(emb) L.normalised_basis(64) # needs sage.rings.number_field emb = K.embeddings(ComplexField())[0] L = E.period_lattice(emb) w1, w2 = L.normalised_basis(); w1, w2 L.is_real() tau = w1/w2; tau
- omega(prec=None, bsd_normalise=False)[source]¶
Return the real or complex volume of this period lattice.
INPUT:
prec
– integer orNone
(default); real precision in bits (default real precision ifNone
)bsd_normalise
– boolean (default:False
); flag to use BSD normalisation in the complex case
OUTPUT:
(real) For real lattices, this is the real period times the number of connected components. For non-real lattices it is the complex area, or double the area if
bsd_normalise
isTrue
.Note
If the curve is given by a global minimal Weierstrass equation, then with
bsd_normalise
=True
, this gives the correct period in the BSD conjecture: the product of this quantity over all embeddings appears in the BSD formula. In general a correction factor is required to make allowance for the model.EXAMPLES:
sage: E = EllipticCurve('37a') sage: E.period_lattice().omega() 5.98691729246392
>>> from sage.all import * >>> E = EllipticCurve('37a') >>> E.period_lattice().omega() 5.98691729246392
E = EllipticCurve('37a') E.period_lattice().omega()
This is not a minimal model:
sage: E = EllipticCurve([0, -432*6^2]) sage: E.period_lattice().omega() 0.486109385710056
>>> from sage.all import * >>> E = EllipticCurve([Integer(0), -Integer(432)*Integer(6)**Integer(2)]) >>> E.period_lattice().omega() 0.486109385710056
E = EllipticCurve([0, -432*6^2]) E.period_lattice().omega()
If you were to plug the above omega into the BSD conjecture, you would get an incorrect value, out by a factor of 2. The following works though:
sage: F = E.minimal_model() sage: F.period_lattice().omega() 0.972218771420113
>>> from sage.all import * >>> F = E.minimal_model() >>> F.period_lattice().omega() 0.972218771420113
F = E.minimal_model() F.period_lattice().omega()
sage: # needs sage.rings.number_field sage: x = polygen(ZZ, 'x') sage: K.<a> = NumberField(x^3 - 2) sage: emb = K.embeddings(RealField())[0] sage: E = EllipticCurve([0,1,0,a,a]) sage: L = E.period_lattice(emb) sage: L.omega(64) 3.81452977217854509
>>> from sage.all import * >>> # needs sage.rings.number_field >>> x = polygen(ZZ, 'x') >>> K = NumberField(x**Integer(3) - Integer(2), names=('a',)); (a,) = K._first_ngens(1) >>> emb = K.embeddings(RealField())[Integer(0)] >>> E = EllipticCurve([Integer(0),Integer(1),Integer(0),a,a]) >>> L = E.period_lattice(emb) >>> L.omega(Integer(64)) 3.81452977217854509
# needs sage.rings.number_field x = polygen(ZZ, 'x') K.<a> = NumberField(x^3 - 2) emb = K.embeddings(RealField())[0] E = EllipticCurve([0,1,0,a,a]) L = E.period_lattice(emb) L.omega(64)
>>> from sage.all import * >>> # needs sage.rings.number_field >>> x = polygen(ZZ, 'x') >>> K = NumberField(x**Integer(3) - Integer(2), names=('a',)); (a,) = K._first_ngens(1) >>> emb = K.embeddings(RealField())[Integer(0)] >>> E = EllipticCurve([Integer(0),Integer(1),Integer(0),a,a]) >>> L = E.period_lattice(emb) >>> L.omega(Integer(64)) 3.81452977217854509
# needs sage.rings.number_field x = polygen(ZZ, 'x') K.<a> = NumberField(x^3 - 2) emb = K.embeddings(RealField())[0] E = EllipticCurve([0,1,0,a,a]) L = E.period_lattice(emb) L.omega(64)
A complex example (taken from J.E.Cremona and E.Whitley, Periods of cusp forms and elliptic curves over imaginary quadratic fields, Mathematics of Computation 62 No. 205 (1994), 407-429). See Issue #29645 and Issue #29782:
sage: # needs sage.rings.number_field sage: K.<i> = QuadraticField(-1) sage: E = EllipticCurve([0,1-i,i,-i,0]) sage: L = E.period_lattice(K.embeddings(CC)[0]) sage: L.omega() 8.80694160502647 sage: L.omega(prec=200) 8.8069416050264741493250743632295462227858630765392114070032 sage: L.omega(bsd_normalise=True) 17.6138832100529
>>> from sage.all import * >>> # needs sage.rings.number_field >>> K = QuadraticField(-Integer(1), names=('i',)); (i,) = K._first_ngens(1) >>> E = EllipticCurve([Integer(0),Integer(1)-i,i,-i,Integer(0)]) >>> L = E.period_lattice(K.embeddings(CC)[Integer(0)]) >>> L.omega() 8.80694160502647 >>> L.omega(prec=Integer(200)) 8.8069416050264741493250743632295462227858630765392114070032 >>> L.omega(bsd_normalise=True) 17.6138832100529
# needs sage.rings.number_field K.<i> = QuadraticField(-1) E = EllipticCurve([0,1-i,i,-i,0]) L = E.period_lattice(K.embeddings(CC)[0]) L.omega() L.omega(prec=200) L.omega(bsd_normalise=True)
- real_period(prec=None, algorithm='sage')[source]¶
Return the real period of this period lattice.
INPUT:
prec
– integer orNone
(default); real precision in bits (default real precision ifNone
)algorithm
– string (default:'sage'
); choice of implementation (for real embeddings only) between'sage'
(native Sage implementation) or'pari'
(use the PARI library: only available for real embeddings)
Note
Only defined for real lattices; a
RuntimeError
is raised for non-real lattices.EXAMPLES:
sage: E = EllipticCurve('37a') sage: E.period_lattice().real_period() 2.99345864623196
>>> from sage.all import * >>> E = EllipticCurve('37a') >>> E.period_lattice().real_period() 2.99345864623196
E = EllipticCurve('37a') E.period_lattice().real_period()
sage: # needs sage.rings.number_field sage: x = polygen(ZZ, 'x') sage: K.<a> = NumberField(x^3 - 2) sage: emb = K.embeddings(RealField())[0] sage: E = EllipticCurve([0,1,0,a,a]) sage: L = E.period_lattice(emb) sage: L.real_period(64) 3.81452977217854509
>>> from sage.all import * >>> # needs sage.rings.number_field >>> x = polygen(ZZ, 'x') >>> K = NumberField(x**Integer(3) - Integer(2), names=('a',)); (a,) = K._first_ngens(1) >>> emb = K.embeddings(RealField())[Integer(0)] >>> E = EllipticCurve([Integer(0),Integer(1),Integer(0),a,a]) >>> L = E.period_lattice(emb) >>> L.real_period(Integer(64)) 3.81452977217854509
# needs sage.rings.number_field x = polygen(ZZ, 'x') K.<a> = NumberField(x^3 - 2) emb = K.embeddings(RealField())[0] E = EllipticCurve([0,1,0,a,a]) L = E.period_lattice(emb) L.real_period(64)
>>> from sage.all import * >>> # needs sage.rings.number_field >>> x = polygen(ZZ, 'x') >>> K = NumberField(x**Integer(3) - Integer(2), names=('a',)); (a,) = K._first_ngens(1) >>> emb = K.embeddings(RealField())[Integer(0)] >>> E = EllipticCurve([Integer(0),Integer(1),Integer(0),a,a]) >>> L = E.period_lattice(emb) >>> L.real_period(Integer(64)) 3.81452977217854509
# needs sage.rings.number_field x = polygen(ZZ, 'x') K.<a> = NumberField(x^3 - 2) emb = K.embeddings(RealField())[0] E = EllipticCurve([0,1,0,a,a]) L = E.period_lattice(emb) L.real_period(64)
- reduce(z)[source]¶
Reduce a complex number modulo the lattice.
INPUT:
z
– complex number
OUTPUT:
(complex) the reduction of \(z\) modulo the lattice, lying in the fundamental period parallelogram with respect to the lattice basis. For curves defined over the reals (i.e. real embeddings) the output will be real when possible.
EXAMPLES:
sage: E = EllipticCurve('389a') sage: L = E.period_lattice() sage: w1, w2 = L.basis(prec=100) sage: P = E([-1,1]) sage: zP = P.elliptic_logarithm(precision=100); zP 0.47934825019021931612953301006 + 0.98586885077582410221120384908*I sage: z = zP + 10*w1 - 20*w2; z 25.381473858740770069343110929 - 38.448885180257139986236950114*I sage: L.reduce(z) 0.47934825019021931612953301006 + 0.98586885077582410221120384908*I sage: L.elliptic_logarithm(2*P) 0.958696500380439 sage: L.reduce(L.elliptic_logarithm(2*P)) 0.958696500380439 sage: L.reduce(L.elliptic_logarithm(2*P) + 10*w1 - 20*w2) 0.958696500380444
>>> from sage.all import * >>> E = EllipticCurve('389a') >>> L = E.period_lattice() >>> w1, w2 = L.basis(prec=Integer(100)) >>> P = E([-Integer(1),Integer(1)]) >>> zP = P.elliptic_logarithm(precision=Integer(100)); zP 0.47934825019021931612953301006 + 0.98586885077582410221120384908*I >>> z = zP + Integer(10)*w1 - Integer(20)*w2; z 25.381473858740770069343110929 - 38.448885180257139986236950114*I >>> L.reduce(z) 0.47934825019021931612953301006 + 0.98586885077582410221120384908*I >>> L.elliptic_logarithm(Integer(2)*P) 0.958696500380439 >>> L.reduce(L.elliptic_logarithm(Integer(2)*P)) 0.958696500380439 >>> L.reduce(L.elliptic_logarithm(Integer(2)*P) + Integer(10)*w1 - Integer(20)*w2) 0.958696500380444
E = EllipticCurve('389a') L = E.period_lattice() w1, w2 = L.basis(prec=100) P = E([-1,1]) zP = P.elliptic_logarithm(precision=100); zP z = zP + 10*w1 - 20*w2; z L.reduce(z) L.elliptic_logarithm(2*P) L.reduce(L.elliptic_logarithm(2*P)) L.reduce(L.elliptic_logarithm(2*P) + 10*w1 - 20*w2)
- sigma(z, prec=None, flag=0)[source]¶
Return the value of the Weierstrass sigma function for this elliptic curve period lattice.
INPUT:
z
– a complex numberprec
– (default:None
) real precision in bits(default real precision if
None
)
flag
–0: (default) ???;
1: computes an arbitrary determination of log(sigma(z))
2, 3: same using the product expansion instead of theta series. ???
Note
The reason for the ???’s above, is that the PARI documentation for ellsigma is very vague. Also this is only implemented for curves defined over \(\QQ\).
Todo
This function does not use any of the PeriodLattice functions and so should be moved to ell_rational_field.
EXAMPLES:
sage: EllipticCurve('389a1').period_lattice().sigma(CC(2,1)) 2.60912163570108 - 0.200865080824587*I
>>> from sage.all import * >>> EllipticCurve('389a1').period_lattice().sigma(CC(Integer(2),Integer(1))) 2.60912163570108 - 0.200865080824587*I
EllipticCurve('389a1').period_lattice().sigma(CC(2,1))
- tau(prec=None, algorithm='sage')[source]¶
Return the upper half-plane parameter in the fundamental region.
INPUT:
prec
– (default:None
) precision in bits (default precision ifNone
)algorithm
– string (default:'sage'
); choice of implementation (for real embeddings only) between ‘sage’ (native Sage implementation) or ‘pari’ (use the PARI library: only available for real embeddings)
OUTPUT:
(Complex) \(\tau = \omega_1/\omega_2\) where the lattice has the form \(\ZZ\omega_1 + \ZZ\omega_2\), normalised so that \(\tau = \omega_1/\omega_2\) is in the fundamental region of the upper half-plane.
EXAMPLES:
sage: E = EllipticCurve('37a') sage: L = E.period_lattice() sage: L.tau() 1.22112736076463*I
>>> from sage.all import * >>> E = EllipticCurve('37a') >>> L = E.period_lattice() >>> L.tau() 1.22112736076463*I
E = EllipticCurve('37a') L = E.period_lattice() L.tau()
sage: # needs sage.rings.number_field sage: x = polygen(ZZ, 'x') sage: K.<a> = NumberField(x^3 - 2) sage: emb = K.embeddings(RealField())[0] sage: E = EllipticCurve([0,1,0,a,a]) sage: L = E.period_lattice(emb) sage: tau = L.tau(); tau -0.338718341018919 + 0.940887817679340*I sage: tau.abs() 1.00000000000000 sage: -0.5 <= tau.real() <= 0.5 True sage: # needs sage.rings.number_field sage: emb = K.embeddings(ComplexField())[0] sage: L = E.period_lattice(emb) sage: tau = L.tau(); tau 0.387694505032876 + 1.30821088214407*I sage: tau.abs() 1.36444961115933 sage: -0.5 <= tau.real() <= 0.5 True
>>> from sage.all import * >>> # needs sage.rings.number_field >>> x = polygen(ZZ, 'x') >>> K = NumberField(x**Integer(3) - Integer(2), names=('a',)); (a,) = K._first_ngens(1) >>> emb = K.embeddings(RealField())[Integer(0)] >>> E = EllipticCurve([Integer(0),Integer(1),Integer(0),a,a]) >>> L = E.period_lattice(emb) >>> tau = L.tau(); tau -0.338718341018919 + 0.940887817679340*I >>> tau.abs() 1.00000000000000 >>> -RealNumber('0.5') <= tau.real() <= RealNumber('0.5') True >>> # needs sage.rings.number_field >>> emb = K.embeddings(ComplexField())[Integer(0)] >>> L = E.period_lattice(emb) >>> tau = L.tau(); tau 0.387694505032876 + 1.30821088214407*I >>> tau.abs() 1.36444961115933 >>> -RealNumber('0.5') <= tau.real() <= RealNumber('0.5') True
# needs sage.rings.number_field x = polygen(ZZ, 'x') K.<a> = NumberField(x^3 - 2) emb = K.embeddings(RealField())[0] E = EllipticCurve([0,1,0,a,a]) L = E.period_lattice(emb) tau = L.tau(); tau tau.abs() -0.5 <= tau.real() <= 0.5 # needs sage.rings.number_field emb = K.embeddings(ComplexField())[0] L = E.period_lattice(emb) tau = L.tau(); tau tau.abs() -0.5 <= tau.real() <= 0.5
>>> from sage.all import * >>> # needs sage.rings.number_field >>> x = polygen(ZZ, 'x') >>> K = NumberField(x**Integer(3) - Integer(2), names=('a',)); (a,) = K._first_ngens(1) >>> emb = K.embeddings(RealField())[Integer(0)] >>> E = EllipticCurve([Integer(0),Integer(1),Integer(0),a,a]) >>> L = E.period_lattice(emb) >>> tau = L.tau(); tau -0.338718341018919 + 0.940887817679340*I >>> tau.abs() 1.00000000000000 >>> -RealNumber('0.5') <= tau.real() <= RealNumber('0.5') True >>> # needs sage.rings.number_field >>> emb = K.embeddings(ComplexField())[Integer(0)] >>> L = E.period_lattice(emb) >>> tau = L.tau(); tau 0.387694505032876 + 1.30821088214407*I >>> tau.abs() 1.36444961115933 >>> -RealNumber('0.5') <= tau.real() <= RealNumber('0.5') True
# needs sage.rings.number_field x = polygen(ZZ, 'x') K.<a> = NumberField(x^3 - 2) emb = K.embeddings(RealField())[0] E = EllipticCurve([0,1,0,a,a]) L = E.period_lattice(emb) tau = L.tau(); tau tau.abs() -0.5 <= tau.real() <= 0.5 # needs sage.rings.number_field emb = K.embeddings(ComplexField())[0] L = E.period_lattice(emb) tau = L.tau(); tau tau.abs() -0.5 <= tau.real() <= 0.5
- sage.schemes.elliptic_curves.period_lattice.extended_agm_iteration(a, b, c)[source]¶
Internal function for the extended AGM used in elliptic logarithm computation.
INPUT:
a
,b
,c
– three real or complex numbers
OUTPUT:
(3-tuple) \((a_0,b_0,c_0)\), the limit of the iteration \((a,b,c) \mapsto ((a+b)/2,\sqrt{ab},(c+\sqrt{c^2+b^2-a^2})/2)\).
EXAMPLES:
sage: # needs sage.rings.real_mpfr sage: from sage.schemes.elliptic_curves.period_lattice import extended_agm_iteration sage: extended_agm_iteration(RR(1), RR(2), RR(3)) (1.45679103104691, 1.45679103104691, 3.21245294970054) sage: extended_agm_iteration(CC(1,2), CC(2,3), CC(3,4)) (1.46242448156430 + 2.47791311676267*I, 1.46242448156430 + 2.47791311676267*I, 3.22202144343535 + 4.28383734262540*I)
>>> from sage.all import * >>> # needs sage.rings.real_mpfr >>> from sage.schemes.elliptic_curves.period_lattice import extended_agm_iteration >>> extended_agm_iteration(RR(Integer(1)), RR(Integer(2)), RR(Integer(3))) (1.45679103104691, 1.45679103104691, 3.21245294970054) >>> extended_agm_iteration(CC(Integer(1),Integer(2)), CC(Integer(2),Integer(3)), CC(Integer(3),Integer(4))) (1.46242448156430 + 2.47791311676267*I, 1.46242448156430 + 2.47791311676267*I, 3.22202144343535 + 4.28383734262540*I)
# needs sage.rings.real_mpfr from sage.schemes.elliptic_curves.period_lattice import extended_agm_iteration extended_agm_iteration(RR(1), RR(2), RR(3)) extended_agm_iteration(CC(1,2), CC(2,3), CC(3,4))
- sage.schemes.elliptic_curves.period_lattice.normalise_periods(w1, w2)[source]¶
Normalise the period basis \((w_1,w_2)\) so that \(w_1/w_2\) is in the fundamental region.
INPUT:
w1
,w2
– two complex numbers with non-real ratio
OUTPUT:
(tuple) \(((\omega_1',\omega_2'),[a,b,c,d])\) where \(a,b,c,d\) are integers such that
\(ad-bc=\pm1\);
\((\omega_1',\omega_2') = (a\omega_1+b\omega_2,c\omega_1+d\omega_2)\);
\(\tau=\omega_1'/\omega_2'\) is in the upper half plane;
\(|\tau|\ge1\) and \(|\Re(\tau)|\le\frac{1}{2}\).
EXAMPLES:
sage: # needs sage.rings.real_mpfr sage.symbolic sage: from sage.schemes.elliptic_curves.period_lattice import reduce_tau, normalise_periods sage: w1 = CC(1.234, 3.456) sage: w2 = CC(1.234, 3.456000001) sage: w1/w2 # in lower half plane! 0.999999999743367 - 9.16334785827644e-11*I sage: w1w2, abcd = normalise_periods(w1, w2) sage: a,b,c,d = abcd sage: w1w2 == (a*w1+b*w2, c*w1+d*w2) True sage: w1w2[0]/w1w2[1] 1.23400010389203e9*I sage: a*d-b*c # note change of orientation -1
>>> from sage.all import * >>> # needs sage.rings.real_mpfr sage.symbolic >>> from sage.schemes.elliptic_curves.period_lattice import reduce_tau, normalise_periods >>> w1 = CC(RealNumber('1.234'), RealNumber('3.456')) >>> w2 = CC(RealNumber('1.234'), RealNumber('3.456000001')) >>> w1/w2 # in lower half plane! 0.999999999743367 - 9.16334785827644e-11*I >>> w1w2, abcd = normalise_periods(w1, w2) >>> a,b,c,d = abcd >>> w1w2 == (a*w1+b*w2, c*w1+d*w2) True >>> w1w2[Integer(0)]/w1w2[Integer(1)] 1.23400010389203e9*I >>> a*d-b*c # note change of orientation -1
# needs sage.rings.real_mpfr sage.symbolic from sage.schemes.elliptic_curves.period_lattice import reduce_tau, normalise_periods w1 = CC(1.234, 3.456) w2 = CC(1.234, 3.456000001) w1/w2 # in lower half plane! w1w2, abcd = normalise_periods(w1, w2) a,b,c,d = abcd w1w2 == (a*w1+b*w2, c*w1+d*w2) w1w2[0]/w1w2[1] a*d-b*c # note change of orientation
- sage.schemes.elliptic_curves.period_lattice.reduce_tau(tau)[source]¶
Transform a point in the upper half plane to the fundamental region.
INPUT:
tau
– complex number with positive imaginary part
OUTPUT:
(tuple) \((\tau',[a,b,c,d])\) where \(a,b,c,d\) are integers such that
\(ad-bc=1\);
\(\tau'=(a\tau+b)/(c\tau+d)\);
\(|\tau'|\ge1\);
\(|\Re(\tau')|\le\frac{1}{2}\).
EXAMPLES:
sage: # needs sage.rings.real_mpfr sage.symbolic sage: from sage.schemes.elliptic_curves.period_lattice import reduce_tau sage: reduce_tau(CC(1.23,3.45)) (0.230000000000000 + 3.45000000000000*I, [1, -1, 0, 1]) sage: reduce_tau(CC(1.23,0.0345)) (-0.463960069171512 + 1.35591888067914*I, [-5, 6, 4, -5]) sage: reduce_tau(CC(1.23,0.0000345)) (0.130000000001761 + 2.89855072463768*I, [13, -16, 100, -123])
>>> from sage.all import * >>> # needs sage.rings.real_mpfr sage.symbolic >>> from sage.schemes.elliptic_curves.period_lattice import reduce_tau >>> reduce_tau(CC(RealNumber('1.23'),RealNumber('3.45'))) (0.230000000000000 + 3.45000000000000*I, [1, -1, 0, 1]) >>> reduce_tau(CC(RealNumber('1.23'),RealNumber('0.0345'))) (-0.463960069171512 + 1.35591888067914*I, [-5, 6, 4, -5]) >>> reduce_tau(CC(RealNumber('1.23'),RealNumber('0.0000345'))) (0.130000000001761 + 2.89855072463768*I, [13, -16, 100, -123])
# needs sage.rings.real_mpfr sage.symbolic from sage.schemes.elliptic_curves.period_lattice import reduce_tau reduce_tau(CC(1.23,3.45)) reduce_tau(CC(1.23,0.0345)) reduce_tau(CC(1.23,0.0000345))