Book a Demo!
CoCalc Logo Icon
StoreFeaturesDocsShareSupportNewsAboutPoliciesSign UpSign In
sagemath
GitHub Repository: sagemath/sagesmc
Path: blob/master/src/sage/schemes/elliptic_curves/period_lattice.py
8820 views
1
# -*- coding: utf-8 -*-
2
r"""
3
Period lattices of elliptic curves and related functions
4
5
Let `E` be an elliptic curve defined over a number field `K`
6
(including `\QQ`). We attach a period lattice (a discrete rank 2
7
subgroup of `\CC`) to each embedding of `K` into `\CC`.
8
9
In the case of real embeddings, the lattice is stable under complex
10
conjugation and is called a real lattice. These have two types:
11
rectangular, (the real curve has two connected components and positive
12
discriminant) or non-rectangular (one connected component, negative
13
discriminant).
14
15
The periods are computed to arbitrary precision using the AGM (Gauss's
16
Arithmetic-Geometric Mean).
17
18
EXAMPLES::
19
20
sage: K.<a> = NumberField(x^3-2)
21
sage: E = EllipticCurve([0,1,0,a,a])
22
23
First we try a real embedding::
24
25
sage: emb = K.embeddings(RealField())[0]
26
sage: L = E.period_lattice(emb); L
27
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:
28
From: Number Field in a with defining polynomial x^3 - 2
29
To: Algebraic Real Field
30
Defn: a |--> 1.259921049894873?
31
32
The first basis period is real::
33
34
sage: L.basis()
35
(3.81452977217855, 1.90726488608927 + 1.34047785962440*I)
36
sage: L.is_real()
37
True
38
39
For a basis `\omega_1,\omega_2` normalised so that `\omega_1/\omega_2`
40
is in the fundamental region of the upper half-plane, use the function
41
``normalised_basis()`` instead::
42
43
sage: L.normalised_basis()
44
(1.90726488608927 - 1.34047785962440*I, -1.90726488608927 - 1.34047785962440*I)
45
46
Next a complex embedding::
47
48
sage: emb = K.embeddings(ComplexField())[0]
49
sage: L = E.period_lattice(emb); L
50
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:
51
From: Number Field in a with defining polynomial x^3 - 2
52
To: Algebraic Field
53
Defn: a |--> -0.6299605249474365? - 1.091123635971722?*I
54
55
In this case, the basis `\omega_1`, `\omega_2` is always normalised so
56
that `\tau = \omega_1/\omega_2` is in the fundamental region in the
57
upper half plane::
58
59
sage: w1,w2 = L.basis(); w1,w2
60
(-1.37588604166076 - 2.58560946624443*I, -2.10339907847356 + 0.428378776460622*I)
61
sage: L.is_real()
62
False
63
sage: tau = w1/w2; tau
64
0.387694505032876 + 1.30821088214407*I
65
sage: L.normalised_basis()
66
(-1.37588604166076 - 2.58560946624443*I, -2.10339907847356 + 0.428378776460622*I)
67
68
We test that bug #8415 (caused by a PARI bug fixed in v2.3.5) is OK::
69
70
sage: E = EllipticCurve('37a')
71
sage: K.<a> = QuadraticField(-7)
72
sage: EK = E.change_ring(K)
73
sage: EK.period_lattice(K.complex_embeddings()[0])
74
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 respect to the embedding Ring morphism:
75
From: Number Field in a with defining polynomial x^2 + 7
76
To: Algebraic Field
77
Defn: a |--> -2.645751311064591?*I
78
79
80
81
AUTHORS:
82
83
- ?: initial version.
84
85
- John Cremona:
86
87
- Adapted to handle real embeddings of number fields, September 2008.
88
89
- Added basis_matrix function, November 2008
90
91
- Added support for complex embeddings, May 2009.
92
93
- Added complex elliptic logs, March 2010; enhanced, October 2010.
94
95
"""
96
97
from sage.modules.free_module import FreeModule_generic_pid
98
from sage.rings.all import ZZ, QQ, RealField, ComplexField, QQbar, AA
99
from sage.rings.real_mpfr import is_RealField
100
from sage.rings.complex_field import is_ComplexField
101
from sage.rings.real_mpfr import RealNumber as RealNumber
102
from sage.rings.complex_number import ComplexNumber as ComplexNumber
103
from sage.rings.number_field.number_field import refine_embedding
104
from sage.rings.infinity import Infinity
105
from sage.schemes.elliptic_curves.constructor import EllipticCurve
106
from sage.misc.cachefunc import cached_method
107
108
class PeriodLattice(FreeModule_generic_pid):
109
"""
110
The class for the period lattice of an algebraic variety.
111
"""
112
pass
113
114
class PeriodLattice_ell(PeriodLattice):
115
r"""
116
The class for the period lattice of an elliptic curve.
117
118
Currently supported are elliptic curves defined over `\QQ`, and
119
elliptic curves defined over a number field with a real or complex
120
embedding, where the lattice constructed depends on that
121
embedding.
122
"""
123
124
def __init__(self, E, embedding=None):
125
r"""
126
Initialises the period lattice by storing the elliptic curve and the embedding.
127
128
INPUT:
129
130
- ``E`` -- an elliptic curve
131
132
- ``embedding`` (defult: ``None``) -- an embedding of the base
133
field `K` of ``E`` into a real or complex field. If
134
``None``:
135
136
- use the built-in coercion to `\RR` for `K=\QQ`;
137
138
- use the first embedding into `\RR` given by
139
``K.embeddings(RealField())``, if there are any;
140
141
- use the first embedding into `\CC` given by
142
``K.embeddings(ComplexField())``, if `K` is totally complex.
143
144
.. note::
145
146
No periods are computed on creation of the lattice; see the
147
functions ``basis()``, ``normalised_basis()`` and
148
``real_period()`` for precision setting.
149
150
EXAMPLES:
151
152
This function is not normally called directly, but will be
153
called by the period_lattice() function of classes
154
ell_number_field and ell_rational_field::
155
156
sage: from sage.schemes.elliptic_curves.period_lattice import PeriodLattice_ell
157
sage: E = EllipticCurve('37a')
158
sage: PeriodLattice_ell(E)
159
Period lattice associated to Elliptic Curve defined by y^2 + y = x^3 - x over Rational Field
160
161
::
162
163
sage: K.<a> = NumberField(x^3-2)
164
sage: emb = K.embeddings(RealField())[0]
165
sage: E = EllipticCurve([0,1,0,a,a])
166
sage: L = PeriodLattice_ell(E,emb); L
167
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:
168
From: Number Field in a with defining polynomial x^3 - 2
169
To: Algebraic Real Field
170
Defn: a |--> 1.259921049894873?
171
172
sage: emb = K.embeddings(ComplexField())[0]
173
sage: L = PeriodLattice_ell(E,emb); L
174
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:
175
From: Number Field in a with defining polynomial x^3 - 2
176
To: Algebraic Field
177
Defn: a |--> -0.6299605249474365? - 1.091123635971722?*I
178
179
TESTS::
180
181
sage: from sage.schemes.elliptic_curves.period_lattice import PeriodLattice_ell
182
sage: K.<a> = NumberField(x^3-2)
183
sage: emb = K.embeddings(RealField())[0]
184
sage: E = EllipticCurve([0,1,0,a,a])
185
sage: L = PeriodLattice_ell(E,emb)
186
sage: L == loads(dumps(L))
187
True
188
"""
189
# First we cache the elliptic curve with this period lattice:
190
191
self.E = E
192
193
# Next we cache the embedding into QQbar or AA which extends
194
# the given embedding:
195
196
K = E.base_field()
197
if embedding is None:
198
embs = K.embeddings(AA)
199
real = len(embs)>0
200
if not real:
201
embs = K.embeddings(QQbar)
202
embedding = embs[0]
203
else:
204
embedding = refine_embedding(embedding,Infinity)
205
real = embedding(K.gen()).imag().is_zero()
206
207
self.embedding = embedding
208
209
# Next we compute and cache (in self.real_flag) the type of
210
# the lattice: +1 for real rectangular, -1 for real
211
# non-rectangular, 0 for non-real:
212
213
self.real_flag = 0
214
if real:
215
self.real_flag = +1
216
if embedding(E.discriminant())<0:
217
self.real_flag = -1
218
219
# The following algebraic data associated to E and the
220
# embedding is cached:
221
#
222
# Ebar: the curve E base-changed to QQbar (or AA)
223
# f2: the 2-division polynomial of Ebar
224
# ei: the roots e1, e2, e3 of f2, as elements of QQbar (or AA)
225
#
226
# The ei are used both for period computation and elliptic
227
# logarithms.
228
229
self.Ebar = self.E.change_ring(self.embedding)
230
self.f2 = self.Ebar.two_division_polynomial()
231
if self.real_flag == 1: # positive discriminant
232
self._ei = self.f2.roots(AA,multiplicities=False)
233
self._ei.sort() # e1 < e2 < e3
234
e1, e2, e3 = self._ei
235
elif self.real_flag == -1: # negative discriminant
236
self._ei = self.f2.roots(QQbar,multiplicities=False)
237
self._ei = list(sorted(self._ei,key=lambda z: z.imag()))
238
e1, e3, e2 = self._ei # so e3 is real
239
e3 = AA(e3)
240
self._ei = [e1, e2, e3]
241
else:
242
self._ei = self.f2.roots(QQbar,multiplicities=False)
243
e1, e2, e3 = self._ei
244
245
# The quantities sqrt(e_i-e_j) are cached (as elements of
246
# QQbar) to be used in period computations:
247
248
self._abc = (e3-e1).sqrt(), (e3-e2).sqrt(), (e2-e1).sqrt()
249
250
PeriodLattice.__init__(self, base_ring=ZZ, rank=2, degree=1, sparse=False)
251
252
def __cmp__(self, other):
253
r"""
254
Comparison function for period lattices
255
256
TESTS::
257
258
sage: from sage.schemes.elliptic_curves.period_lattice import PeriodLattice_ell
259
sage: K.<a> = NumberField(x^3-2)
260
sage: E = EllipticCurve([0,1,0,a,a])
261
sage: embs = K.embeddings(ComplexField())
262
sage: L1,L2,L3 = [PeriodLattice_ell(E,e) for e in embs]
263
sage: L1 < L2 < L3
264
True
265
266
"""
267
if not isinstance(other, PeriodLattice_ell): return -1
268
t = cmp(self.E, other.E)
269
if t: return t
270
a = self.E.base_field().gen()
271
return cmp(self.embedding(a), other.embedding(a))
272
273
def __repr__(self):
274
"""
275
Returns the string representation of this period lattice.
276
277
EXAMPLES::
278
279
sage: E = EllipticCurve('37a')
280
sage: E.period_lattice()
281
Period lattice associated to Elliptic Curve defined by y^2 + y = x^3 - x over Rational Field
282
283
::
284
285
sage: K.<a> = NumberField(x^3-2)
286
sage: emb = K.embeddings(RealField())[0]
287
sage: E = EllipticCurve([0,1,0,a,a])
288
sage: L = E.period_lattice(emb); L
289
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:
290
From: Number Field in a with defining polynomial x^3 - 2
291
To: Algebraic Real Field
292
Defn: a |--> 1.259921049894873?
293
"""
294
if self.E.base_field() is QQ:
295
return "Period lattice associated to %s"%(self.E)
296
else:
297
return "Period lattice associated to %s with respect to the embedding %s"%(self.E, self.embedding)
298
299
def __call__(self, P, prec=None):
300
r"""
301
Return the elliptic logarithm of a point `P`.
302
303
INPUT:
304
305
- ``P`` (point) -- a point on the elliptic curve associated
306
with this period lattice.
307
308
- ``prec`` (default: ``None``) -- precision in bits (default
309
precision if ``None``).
310
311
OUTPUT:
312
313
(complex number) The elliptic logarithm of the point `P` with
314
respect to this period lattice. If `E` is the elliptic curve
315
and `\sigma:K\to\CC` the embedding, then the returned value `z`
316
is such that `z\pmod{L}` maps to `\sigma(P)` under the
317
standard Weierstrass isomorphism from `\CC/L` to `\sigma(E)`.
318
319
EXAMPLES::
320
321
sage: E = EllipticCurve('389a')
322
sage: L = E.period_lattice()
323
sage: E.discriminant() > 0
324
True
325
sage: L.real_flag
326
1
327
sage: P = E([-1,1])
328
sage: P.is_on_identity_component ()
329
False
330
sage: L(P, prec=96)
331
0.4793482501902193161295330101 + 0.985868850775824102211203849...*I
332
sage: Q=E([3,5])
333
sage: Q.is_on_identity_component()
334
True
335
sage: L(Q, prec=96)
336
1.931128271542559442488585220
337
338
Note that this is actually the inverse of the Weierstrass isomorphism::
339
340
sage: L.elliptic_exponential(L(Q))
341
(3.00000000000000 : 5.00000000000000 : 1.00000000000000)
342
343
An example with negative discriminant, and a torsion point::
344
345
sage: E = EllipticCurve('11a1')
346
sage: L = E.period_lattice()
347
sage: E.discriminant() < 0
348
True
349
sage: L.real_flag
350
-1
351
sage: P = E([16,-61])
352
sage: L(P)
353
0.253841860855911
354
sage: L.real_period() / L(P)
355
5.00000000000000
356
"""
357
return self.elliptic_logarithm(P,prec)
358
359
@cached_method
360
def basis(self, prec=None, algorithm='sage'):
361
r"""
362
Return a basis for this period lattice as a 2-tuple.
363
364
INPUT:
365
366
- ``prec`` (default: ``None``) -- precision in bits (default
367
precision if ``None``).
368
369
- ``algorithm`` (string, default 'sage') -- choice of
370
implementation (for real embeddings only) between 'sage'
371
(native Sage implementation) or 'pari' (use the PARI
372
library: only available for real embeddings).
373
374
OUTPUT:
375
376
(tuple of Complex) `(\omega_1,\omega_2)` where the lattice is
377
`\ZZ\omega_1 + \ZZ\omega_2`. If the lattice is real then
378
`\omega_1` is real and positive, `\Im(\omega_2)>0` and
379
`\Re(\omega_1/\omega_2)` is either `0` (for rectangular
380
lattices) or `\frac{1}{2}` (for non-rectangular lattices).
381
Otherwise, `\omega_1/\omega_2` is in the fundamental region of
382
the upper half-plane. If the latter normalisation is required
383
for real lattices, use the function ``normalised_basis()``
384
instead.
385
386
EXAMPLES::
387
388
sage: E = EllipticCurve('37a')
389
sage: E.period_lattice().basis()
390
(2.99345864623196, 2.45138938198679*I)
391
392
This shows that the issue reported at trac \#3954 is fixed::
393
394
sage: E = EllipticCurve('37a')
395
sage: b1 = E.period_lattice().basis(prec=30)
396
sage: b2 = E.period_lattice().basis(prec=30)
397
sage: b1 == b2
398
True
399
400
This shows that the issue reported at trac \#4064 is fixed::
401
402
sage: E = EllipticCurve('37a')
403
sage: E.period_lattice().basis(prec=30)[0].parent()
404
Real Field with 30 bits of precision
405
sage: E.period_lattice().basis(prec=100)[0].parent()
406
Real Field with 100 bits of precision
407
408
::
409
410
sage: K.<a> = NumberField(x^3-2)
411
sage: emb = K.embeddings(RealField())[0]
412
sage: E = EllipticCurve([0,1,0,a,a])
413
sage: L = E.period_lattice(emb)
414
sage: L.basis(64)
415
(3.81452977217854509, 1.90726488608927255 + 1.34047785962440202*I)
416
417
sage: emb = K.embeddings(ComplexField())[0]
418
sage: L = E.period_lattice(emb)
419
sage: w1,w2 = L.basis(); w1,w2
420
(-1.37588604166076 - 2.58560946624443*I, -2.10339907847356 + 0.428378776460622*I)
421
sage: L.is_real()
422
False
423
sage: tau = w1/w2; tau
424
0.387694505032876 + 1.30821088214407*I
425
"""
426
# We divide into two cases: (1) Q, or a number field with a
427
# real embedding; (2) a number field with a complex embedding.
428
# In each case the periods are computed by a different
429
# internal function.
430
431
if self.is_real():
432
return self._compute_periods_real(prec=prec, algorithm=algorithm)
433
else:
434
return self._compute_periods_complex(prec=prec)
435
436
@cached_method
437
def normalised_basis(self, prec=None, algorithm='sage'):
438
r"""
439
Return a normalised basis for this period lattice as a 2-tuple.
440
441
INPUT:
442
443
- ``prec`` (default: ``None``) -- precision in bits (default
444
precision if ``None``).
445
446
- ``algorithm`` (string, default 'sage') -- choice of
447
implementation (for real embeddings only) between 'sage'
448
(native Sage implementation) or 'pari' (use the PARI
449
library: only available for real embeddings).
450
451
OUTPUT:
452
453
(tuple of Complex) `(\omega_1,\omega_2)` where the lattice has
454
the form `\ZZ\omega_1 + \ZZ\omega_2`. The basis is normalised
455
so that `\omega_1/\omega_2` is in the fundamental region of
456
the upper half-plane. For an alternative normalisation for
457
real lattices (with the first period real), use the function
458
basis() instead.
459
460
EXAMPLES::
461
462
sage: E = EllipticCurve('37a')
463
sage: E.period_lattice().normalised_basis()
464
(2.99345864623196, -2.45138938198679*I)
465
466
::
467
468
sage: K.<a> = NumberField(x^3-2)
469
sage: emb = K.embeddings(RealField())[0]
470
sage: E = EllipticCurve([0,1,0,a,a])
471
sage: L = E.period_lattice(emb)
472
sage: L.normalised_basis(64)
473
(1.90726488608927255 - 1.34047785962440202*I, -1.90726488608927255 - 1.34047785962440202*I)
474
475
sage: emb = K.embeddings(ComplexField())[0]
476
sage: L = E.period_lattice(emb)
477
sage: w1,w2 = L.normalised_basis(); w1,w2
478
(-1.37588604166076 - 2.58560946624443*I, -2.10339907847356 + 0.428378776460622*I)
479
sage: L.is_real()
480
False
481
sage: tau = w1/w2; tau
482
0.387694505032876 + 1.30821088214407*I
483
"""
484
w1, w2 = periods = self.basis(prec=prec, algorithm=algorithm)
485
periods, mat = normalise_periods(w1,w2)
486
return periods
487
488
@cached_method
489
def _compute_periods_real(self, prec=None, algorithm='sage'):
490
r"""
491
Internal function to compute the periods (real embedding case).
492
493
INPUT:
494
495
496
- `prec` (int or ``None`` (default)) -- floating point
497
precision (in bits); if None, use the default precision.
498
499
- `algorithm` (string, default 'sage') -- choice of implementation between
500
- `pari`: use the PARI library
501
502
- `sage`: use a native Sage implementation (with the same underlying algorithm).
503
504
505
OUTPUT:
506
507
(tuple of Complex) `(\omega_1,\omega_2)` where the lattice has
508
the form `\ZZ\omega_1 + \ZZ\omega_2`, `\omega_1` is real and
509
`\omega_1/\omega_2` has real part either `0` or `frac{1}{2}`.
510
511
EXAMPLES::
512
513
sage: K.<a> = NumberField(x^3-2)
514
sage: E = EllipticCurve([0,1,0,a,a])
515
sage: embs = K.embeddings(CC)
516
sage: Ls = [E.period_lattice(e) for e in embs]
517
sage: [L.is_real() for L in Ls]
518
[False, False, True]
519
sage: Ls[2]._compute_periods_real(100)
520
(3.8145297721785450936365098936,
521
1.9072648860892725468182549468 + 1.3404778596244020196600112394*I)
522
sage: Ls[2]._compute_periods_real(100, algorithm='pari')
523
(3.8145297721785450936365098936,
524
1.9072648860892725468182549468 - 1.3404778596244020196600112394*I)
525
"""
526
if prec is None:
527
prec = RealField().precision()
528
R = RealField(prec)
529
C = ComplexField(prec)
530
531
if algorithm=='pari':
532
if self.E.base_field() is QQ:
533
periods = self.E.pari_curve(prec).omega().python()
534
return (R(periods[0]), C(periods[1]))
535
536
from sage.libs.pari.all import pari
537
E_pari = pari([R(self.embedding(ai).real()) for ai in self.E.a_invariants()]).ellinit(precision=prec)
538
periods = E_pari.omega().python()
539
return (R(periods[0]), C(periods[1]))
540
541
if algorithm!='sage':
542
raise ValueError, "invalid value of 'algorithm' parameter"
543
544
pi = R.pi()
545
# Up to now everything has been exact in AA or QQbar, but now
546
# we must go transcendental. Only now is the desired
547
# precision used!
548
if self.real_flag == 1: # positive discriminant
549
a, b, c = (R(x) for x in self._abc)
550
w1 = R(pi/a.agm(b)) # least real period
551
w2 = C(0,pi/a.agm(c)) # least pure imaginary period
552
else:
553
a = C(self._abc[0])
554
x, y, r = a.real().abs(), a.imag().abs(), a.abs()
555
w1 = R(pi/r.agm(x)) # least real period
556
w2 = R(pi/r.agm(y)) # least pure imaginary period /i
557
w2 = C(w1,w2)/2
558
559
return (w1,w2)
560
561
@cached_method
562
def _compute_periods_complex(self, prec=None, normalise=True):
563
r"""
564
Internal function to compute the periods (complex embedding case).
565
566
INPUT:
567
568
- `prec` (int or ``None`` (default)) -- floating point precision (in bits); if None,
569
use the default precision.
570
571
- `normalise` (bool, default True) -- whether to normalise the
572
basis after computation.
573
574
OUTPUT:
575
576
(tuple of Complex) `(\omega_1,\omega_2)` where the lattice has
577
the form `\ZZ\omega_1 + \ZZ\omega_2`. If `normalise` is
578
`True`, the basis is normalised so that `(\omega_1/\omega_2)`
579
is in the fundamental region of the upper half plane.
580
581
EXAMPLES::
582
583
sage: K.<a> = NumberField(x^3-2)
584
sage: E = EllipticCurve([0,1,0,a,a])
585
sage: embs = K.embeddings(CC)
586
sage: Ls = [E.period_lattice(e) for e in embs]
587
sage: [L.is_real() for L in Ls]
588
[False, False, True]
589
sage: L = Ls[0]
590
sage: w1,w2 = L._compute_periods_complex(100); w1,w2
591
(-1.3758860416607626645495991458 - 2.5856094662444337042877901304*I, -2.1033990784735587243397865076 + 0.42837877646062187766760569686*I)
592
sage: tau = w1/w2; tau
593
0.38769450503287609349437509561 + 1.3082108821440725664008561928*I
594
sage: tau.real()
595
0.38769450503287609349437509561
596
sage: tau.abs()
597
1.3644496111593345713923386773
598
599
Without normalisation::
600
601
sage: w1,w2 = L._compute_periods_complex(normalise=False); w1,w2
602
(2.10339907847356 - 0.428378776460622*I, 0.727513036812796 - 3.01398824270506*I)
603
sage: tau = w1/w2; tau
604
0.293483964608883 + 0.627038168678760*I
605
sage: tau.real()
606
0.293483964608883
607
sage: tau.abs() # > 1
608
0.692321964451917
609
"""
610
if prec is None:
611
prec = RealField().precision()
612
C = ComplexField(prec)
613
614
# Up to now everything has been exact in AA, but now we
615
# must go transcendental. Only now is the desired
616
# precision used!
617
pi = C.pi()
618
a, b, c = (C(x) for x in self._abc)
619
if (a+b).abs() < (a-b).abs(): b=-b
620
if (a+c).abs() < (a-c).abs(): c=-c
621
w1 = pi/a.agm(b)
622
w2 = pi*C.gen()/a.agm(c)
623
if (w1/w2).imag()<0: w2=-w2
624
if normalise:
625
w1w2, mat = normalise_periods(w1,w2)
626
return w1w2
627
return (w1,w2)
628
629
def is_real(self):
630
r"""
631
Return True if this period lattice is real.
632
633
EXAMPLES::
634
635
sage: f = EllipticCurve('11a')
636
sage: f.period_lattice().is_real()
637
True
638
639
::
640
641
sage: K.<i> = QuadraticField(-1)
642
sage: E = EllipticCurve(K,[0,0,0,i,2*i])
643
sage: emb = K.embeddings(ComplexField())[0]
644
sage: L = E.period_lattice(emb)
645
sage: L.is_real()
646
False
647
648
::
649
650
sage: K.<a> = NumberField(x^3-2)
651
sage: E = EllipticCurve([0,1,0,a,a])
652
sage: [E.period_lattice(emb).is_real() for emb in K.embeddings(CC)]
653
[False, False, True]
654
655
656
ALGORITHM:
657
658
The lattice is real if it is associated to a real embedding;
659
such lattices are stable under conjugation.
660
"""
661
return self.real_flag!=0
662
663
def is_rectangular(self):
664
r"""
665
Return True if this period lattice is rectangular.
666
667
.. note::
668
669
Only defined for real lattices; a RuntimeError is raised for
670
non-real lattices.
671
672
EXAMPLES::
673
674
sage: f = EllipticCurve('11a')
675
sage: f.period_lattice().basis()
676
(1.26920930427955, 0.634604652139777 + 1.45881661693850*I)
677
sage: f.period_lattice().is_rectangular()
678
False
679
680
::
681
682
sage: f = EllipticCurve('37b')
683
sage: f.period_lattice().basis()
684
(1.08852159290423, 1.76761067023379*I)
685
sage: f.period_lattice().is_rectangular()
686
True
687
688
ALGORITHM:
689
690
The period lattice is rectangular precisely if the
691
discriminant of the Weierstrass equation is positive, or
692
equivalently if the number of real components is 2.
693
"""
694
if self.is_real():
695
return self.real_flag == +1
696
raise RuntimeError, "Not defined for non-real lattices."
697
698
def real_period(self, prec = None, algorithm='sage'):
699
"""
700
Returns the real period of this period lattice.
701
702
INPUT:
703
704
- ``prec`` (int or ``None`` (default)) -- real precision in
705
bits (default real precision if ``None``)
706
707
- ``algorithm`` (string, default 'sage') -- choice of
708
implementation (for real embeddings only) between 'sage'
709
(native Sage implementation) or 'pari' (use the PARI
710
library: only available for real embeddings).
711
712
.. note::
713
714
Only defined for real lattices; a RuntimeError is raised for
715
non-real lattices.
716
717
EXAMPLES::
718
719
sage: E = EllipticCurve('37a')
720
sage: E.period_lattice().real_period()
721
2.99345864623196
722
723
::
724
725
sage: K.<a> = NumberField(x^3-2)
726
sage: emb = K.embeddings(RealField())[0]
727
sage: E = EllipticCurve([0,1,0,a,a])
728
sage: L = E.period_lattice(emb)
729
sage: L.real_period(64)
730
3.81452977217854509
731
"""
732
if self.is_real():
733
return self.basis(prec,algorithm)[0]
734
raise RuntimeError, "Not defined for non-real lattices."
735
736
def omega(self, prec = None):
737
r"""
738
Returns the real or complex volume of this period lattice.
739
740
INPUT:
741
742
- ``prec`` (int or ``None``(default)) -- real precision in
743
bits (default real precision if ``None``)
744
745
OUTPUT:
746
747
(real) For real lattices, this is the real period times the
748
number of connected components. For non-real lattices it is
749
the complex area.
750
751
.. note::
752
753
If the curve is defined over `\QQ` and is given by a
754
*minimal* Weierstrass equation, then this is the correct
755
period in the BSD conjecture, i.e., it is the least real
756
period * 2 when the period lattice is rectangular. More
757
generally the product of this quantity over all embeddings
758
appears in the generalised BSD formula.
759
760
761
EXAMPLES::
762
763
sage: E = EllipticCurve('37a')
764
sage: E.period_lattice().omega()
765
5.98691729246392
766
767
This is not a minimal model::
768
769
sage: E = EllipticCurve([0,-432*6^2])
770
sage: E.period_lattice().omega()
771
0.486109385710056
772
773
If you were to plug the above omega into the BSD conjecture, you
774
would get nonsense. The following works though::
775
776
sage: F = E.minimal_model()
777
sage: F.period_lattice().omega()
778
0.972218771420113
779
780
::
781
782
sage: K.<a> = NumberField(x^3-2)
783
sage: emb = K.embeddings(RealField())[0]
784
sage: E = EllipticCurve([0,1,0,a,a])
785
sage: L = E.period_lattice(emb)
786
sage: L.omega(64)
787
3.81452977217854509
788
789
A complex example (taken from J.E.Cremona and E.Whitley,
790
*Periods of cusp forms and elliptic curves over imaginary
791
quadratic fields*, Mathematics of Computation 62 No. 205
792
(1994), 407-429)::
793
794
sage: K.<i> = QuadraticField(-1)
795
sage: E = EllipticCurve([0,1-i,i,-i,0])
796
sage: L = E.period_lattice(K.embeddings(CC)[0])
797
sage: L.omega()
798
8.80694160502647
799
"""
800
if self.is_real():
801
n_components = (self.real_flag+3)//2
802
return self.real_period(prec) * n_components
803
else:
804
return self.complex_area()
805
806
@cached_method
807
def basis_matrix(self, prec=None, normalised=False):
808
r"""
809
Return the basis matrix of this period lattice.
810
811
INPUT:
812
813
- ``prec`` (int or ``None``(default)) -- real precision in
814
bits (default real precision if ``None``).
815
816
- ``normalised`` (bool, default None) -- if True and the
817
embedding is real, use the normalised basis (see
818
``normalised_basis()``) instead of the default.
819
820
OUTPUT:
821
822
A 2x2 real matrix whose rows are the lattice basis vectors,
823
after identifying `\CC` with `\RR^2`.
824
825
EXAMPLES::
826
827
sage: E = EllipticCurve('37a')
828
sage: E.period_lattice().basis_matrix()
829
[ 2.99345864623196 0.000000000000000]
830
[0.000000000000000 2.45138938198679]
831
832
::
833
834
sage: K.<a> = NumberField(x^3-2)
835
sage: emb = K.embeddings(RealField())[0]
836
sage: E = EllipticCurve([0,1,0,a,a])
837
sage: L = E.period_lattice(emb)
838
sage: L.basis_matrix(64)
839
[ 3.81452977217854509 0.000000000000000000]
840
[ 1.90726488608927255 1.34047785962440202]
841
842
See \#4388::
843
844
sage: L = EllipticCurve('11a1').period_lattice()
845
sage: L.basis_matrix()
846
[ 1.26920930427955 0.000000000000000]
847
[0.634604652139777 1.45881661693850]
848
sage: L.basis_matrix(normalised=True)
849
[0.634604652139777 -1.45881661693850]
850
[-1.26920930427955 0.000000000000000]
851
852
::
853
854
sage: L = EllipticCurve('389a1').period_lattice()
855
sage: L.basis_matrix()
856
[ 2.49021256085505 0.000000000000000]
857
[0.000000000000000 1.97173770155165]
858
sage: L.basis_matrix(normalised=True)
859
[ 2.49021256085505 0.000000000000000]
860
[0.000000000000000 -1.97173770155165]
861
"""
862
from sage.matrix.all import Matrix
863
864
if normalised:
865
return Matrix([list(w) for w in self.normalised_basis(prec)])
866
867
w1,w2 = self.basis(prec)
868
if self.is_real():
869
return Matrix([[w1,0],list(w2)])
870
else:
871
return Matrix([list(w) for w in w1,w2])
872
873
def complex_area(self, prec=None):
874
"""
875
Return the area of a fundamental domain for the period lattice
876
of the elliptic curve.
877
878
INPUT:
879
880
- ``prec`` (int or ``None``(default)) -- real precision in
881
bits (default real precision if ``None``).
882
883
EXAMPLES::
884
885
sage: E = EllipticCurve('37a')
886
sage: E.period_lattice().complex_area()
887
7.33813274078958
888
889
::
890
891
sage: K.<a> = NumberField(x^3-2)
892
sage: embs = K.embeddings(ComplexField())
893
sage: E = EllipticCurve([0,1,0,a,a])
894
sage: [E.period_lattice(emb).is_real() for emb in K.embeddings(CC)]
895
[False, False, True]
896
sage: [E.period_lattice(emb).complex_area() for emb in embs]
897
[6.02796894766694, 6.02796894766694, 5.11329270448345]
898
"""
899
w1,w2 = self.basis(prec)
900
return (w1*w2.conjugate()).imag().abs()
901
902
def sigma(self, z, prec = None, flag=0):
903
r"""
904
Returns the value of the Weierstrass sigma function for this elliptic curve period lattice.
905
906
INPUT:
907
908
- ``z`` -- a complex number
909
910
- ``prec`` (default: ``None``) -- real precision in bits
911
(default real precision if None).
912
913
- ``flag`` --
914
915
0: (default) ???;
916
917
1: computes an arbitrary determination of log(sigma(z))
918
919
2, 3: same using the product expansion instead of theta series. ???
920
921
.. note::
922
923
The reason for the ???'s above, is that the PARI
924
documentation for ellsigma is very vague. Also this is
925
only implemented for curves defined over `\QQ`.
926
927
TODO:
928
929
This function does not use any of the PeriodLattice functions
930
and so should be moved to ell_rational_field.
931
932
EXAMPLES::
933
934
sage: EllipticCurve('389a1').period_lattice().sigma(CC(2,1))
935
2.60912163570108 - 0.200865080824587*I
936
"""
937
if prec is None:
938
prec = RealField().precision()
939
try:
940
return self.E.pari_curve(prec).ellsigma(z, flag)
941
except AttributeError:
942
raise NotImplementedError, "sigma function not yet implemented for period lattices of curves not defined over Q."
943
944
def curve(self):
945
r"""
946
Return the elliptic curve associated with this period lattice.
947
948
EXAMPLES::
949
950
sage: E = EllipticCurve('37a')
951
sage: L = E.period_lattice()
952
sage: L.curve() is E
953
True
954
955
::
956
957
sage: K.<a> = NumberField(x^3-2)
958
sage: E = EllipticCurve([0,1,0,a,a])
959
sage: L = E.period_lattice(K.embeddings(RealField())[0])
960
sage: L.curve() is E
961
True
962
963
sage: L = E.period_lattice(K.embeddings(ComplexField())[0])
964
sage: L.curve() is E
965
True
966
"""
967
return self.E
968
969
def ei(self):
970
r"""
971
Return the x-coordinates of the 2-division points of the elliptic curve associated with this period lattice, as elements of QQbar.
972
973
EXAMPLES::
974
975
sage: E = EllipticCurve('37a')
976
sage: L = E.period_lattice()
977
sage: L.ei()
978
[-1.107159871688768?, 0.2695944364054446?, 0.8375654352833230?]
979
980
::
981
982
sage: K.<a> = NumberField(x^3-2)
983
sage: E = EllipticCurve([0,1,0,a,a])
984
sage: L = E.period_lattice(K.embeddings(RealField())[0])
985
sage: L.ei()
986
[0.?e-19 - 1.122462048309373?*I, 0.?e-19 + 1.122462048309373?*I, -1]
987
988
sage: L = E.period_lattice(K.embeddings(ComplexField())[0])
989
sage: L.ei()
990
[-1.000000000000000? + 0.?e-1...*I,
991
-0.9720806486198328? - 0.561231024154687?*I,
992
0.9720806486198328? + 0.561231024154687?*I]
993
"""
994
return self._ei
995
996
def coordinates(self, z, rounding=None):
997
r"""
998
Returns the coordinates of a complex number w.r.t. the lattice basis
999
1000
INPUT:
1001
1002
- ``z`` (complex) -- A complex number.
1003
1004
- ``rounding`` (default ``None``) -- whether and how to round the
1005
output (see below).
1006
1007
OUTPUT:
1008
1009
When ``rounding`` is ``None`` (the default), returns a tuple
1010
of reals `x`, `y` such that `z=xw_1+yw_2` where `w_1`, `w_2`
1011
are a basis for the lattice (normalised in the case of complex
1012
embeddings).
1013
1014
When ``rounding`` is 'round', returns a tuple of integers `n_1`,
1015
`n_2` which are the closest integers to the `x`, `y` defined
1016
above. If `z` is in the lattice these are the coordinates of
1017
`z` with respect to the lattice basis.
1018
1019
When ``rounding`` is 'floor', returns a tuple of integers
1020
`n_1`, `n_2` which are the integer parts to the `x`, `y`
1021
defined above. These are used in :meth:``.reduce``
1022
1023
EXAMPLES::
1024
1025
sage: E = EllipticCurve('389a')
1026
sage: L = E.period_lattice()
1027
sage: w1, w2 = L.basis(prec=100)
1028
sage: P = E([-1,1])
1029
sage: zP = P.elliptic_logarithm(precision=100); zP
1030
0.47934825019021931612953301006 + 0.98586885077582410221120384908*I
1031
sage: L.coordinates(zP)
1032
(0.19249290511394227352563996419, 0.50000000000000000000000000000)
1033
sage: sum([x*w for x,w in zip(L.coordinates(zP), L.basis(prec=100))])
1034
0.47934825019021931612953301006 + 0.98586885077582410221120384908*I
1035
1036
sage: L.coordinates(12*w1+23*w2)
1037
(12.000000000000000000000000000, 23.000000000000000000000000000)
1038
sage: L.coordinates(12*w1+23*w2, rounding='floor')
1039
(11, 22)
1040
sage: L.coordinates(12*w1+23*w2, rounding='round')
1041
(12, 23)
1042
"""
1043
C = z.parent()
1044
if is_RealField(C):
1045
C = ComplexField(C.precision())
1046
z = C(z)
1047
else:
1048
if is_ComplexField(C):
1049
pass
1050
else:
1051
try:
1052
C = ComplexField()
1053
z = C(z)
1054
except TypeError:
1055
raise TypeError, "%s is not a complex number"%z
1056
prec = C.precision()
1057
from sage.matrix.all import Matrix
1058
from sage.modules.all import vector
1059
if self.real_flag:
1060
w1,w2 = self.basis(prec)
1061
M = Matrix([[w1,0], list(w2)])**(-1)
1062
else:
1063
w1,w2 = self.normalised_basis(prec)
1064
M = Matrix([list(w1), list(w2)])**(-1)
1065
u,v = vector(z)*M
1066
# Now z = u*w1+v*w2
1067
if rounding=='round':
1068
return u.round(), v.round()
1069
if rounding=='floor':
1070
return u.floor(), v.floor()
1071
return u,v
1072
1073
def reduce(self, z):
1074
r"""
1075
Reduce a complex number modulo the lattice
1076
1077
INPUT:
1078
1079
- ``z`` (complex) -- A complex number.
1080
1081
OUTPUT:
1082
1083
(complex) the reduction of `z` modulo the lattice, lying in
1084
the fundamental period parallelogram with respect to the
1085
lattice basis. For curves defined over the reals (i.e. real
1086
embeddings) the output will be real when possible.
1087
1088
EXAMPLES::
1089
1090
sage: E = EllipticCurve('389a')
1091
sage: L = E.period_lattice()
1092
sage: w1, w2 = L.basis(prec=100)
1093
sage: P = E([-1,1])
1094
sage: zP = P.elliptic_logarithm(precision=100); zP
1095
0.47934825019021931612953301006 + 0.98586885077582410221120384908*I
1096
sage: z = zP+10*w1-20*w2; z
1097
25.381473858740770069343110929 - 38.448885180257139986236950114*I
1098
sage: L.reduce(z)
1099
0.47934825019021931612953301006 + 0.98586885077582410221120384908*I
1100
sage: L.elliptic_logarithm(2*P)
1101
0.958696500380439
1102
sage: L.reduce(L.elliptic_logarithm(2*P))
1103
0.958696500380439
1104
sage: L.reduce(L.elliptic_logarithm(2*P)+10*w1-20*w2)
1105
0.958696500380444
1106
"""
1107
C = z.parent()
1108
z_is_real = False
1109
if is_RealField(C):
1110
z_is_real = True
1111
C = ComplexField(C.precision())
1112
z = C(z)
1113
else:
1114
if is_ComplexField(C):
1115
z_is_real = z.is_real()
1116
else:
1117
try:
1118
C = ComplexField()
1119
z = C(z)
1120
z_is_real = z.is_real()
1121
except TypeError:
1122
raise TypeError, "%s is not a complex number"%z
1123
prec = C.precision()
1124
if self.real_flag:
1125
w1,w2 = self.basis(prec) # w1 real
1126
else:
1127
w1,w2 = self.normalised_basis(prec)
1128
# print "z = ",z
1129
# print "w1,w2 = ",w1,w2
1130
u,v = self.coordinates(z, rounding='floor')
1131
# print "u,v = ",u,v
1132
z = z-u*w1-v*w2
1133
1134
# Final adjustments for the real case.
1135
1136
# NB We assume here that when the embedding is real then the
1137
# point is also real!
1138
1139
if self.real_flag == 0: return z
1140
if self.real_flag == -1:
1141
k = (z.imag()/w2.imag()).round()
1142
z = z-k*w2
1143
return C(z.real(),0)
1144
1145
if ((2*z.imag()/w2.imag()).round())%2:
1146
return C(z.real(),w2.imag()/2)
1147
else:
1148
return C(z.real(),0)
1149
1150
def elliptic_logarithm(self, P, prec=None, reduce=True):
1151
r"""
1152
Return the elliptic logarithm of a point.
1153
1154
INPUT:
1155
1156
- ``P`` (point) -- A point on the elliptic curve associated
1157
with this period lattice.
1158
1159
- ``prec`` (default: ``None``) -- real precision in bits
1160
(default real precision if None).
1161
1162
- ``reduce`` (default: ``True``) -- if ``True``, the result
1163
is reduced with respect to the period lattice basis.
1164
1165
OUTPUT:
1166
1167
(complex number) The elliptic logarithm of the point `P` with
1168
respect to this period lattice. If `E` is the elliptic curve
1169
and `\sigma:K\to\CC` the embedding, the the returned value `z`
1170
is such that `z\pmod{L}` maps to `\sigma(P)` under the
1171
standard Weierstrass isomorphism from `\CC/L` to `\sigma(E)`.
1172
If ``reduce`` is ``True``, the output is reduced so that it is
1173
in the fundamental period parallelogram with respect to the
1174
normalised lattice basis.
1175
1176
ALGORITHM: Uses the complex AGM. See [Cremona2010]_ for details.
1177
1178
.. [Cremona2010] J. E. Cremona and T. Thongjunthug, The
1179
Complex AGM, periods of elliptic curves over $\CC$ and
1180
complex elliptic logarithms. Preprint 2010.
1181
1182
EXAMPLES::
1183
1184
sage: E = EllipticCurve('389a')
1185
sage: L = E.period_lattice()
1186
sage: E.discriminant() > 0
1187
True
1188
sage: L.real_flag
1189
1
1190
sage: P = E([-1,1])
1191
sage: P.is_on_identity_component ()
1192
False
1193
sage: L.elliptic_logarithm(P, prec=96)
1194
0.4793482501902193161295330101 + 0.9858688507758241022112038491*I
1195
sage: Q=E([3,5])
1196
sage: Q.is_on_identity_component()
1197
True
1198
sage: L.elliptic_logarithm(Q, prec=96)
1199
1.931128271542559442488585220
1200
1201
Note that this is actually the inverse of the Weierstrass isomorphism::
1202
1203
sage: L.elliptic_exponential(_)
1204
(3.00000000000000000000000000... : 5.00000000000000000000000000... : 1.000000000000000000000000000)
1205
1206
An example with negative discriminant, and a torsion point::
1207
1208
sage: E = EllipticCurve('11a1')
1209
sage: L = E.period_lattice()
1210
sage: E.discriminant() < 0
1211
True
1212
sage: L.real_flag
1213
-1
1214
sage: P = E([16,-61])
1215
sage: L.elliptic_logarithm(P)
1216
0.253841860855911
1217
sage: L.real_period() / L.elliptic_logarithm(P)
1218
5.00000000000000
1219
1220
An example where precision is problematic::
1221
1222
sage: E = EllipticCurve([1, 0, 1, -85357462, 303528987048]) #18074g1
1223
sage: P = E([4458713781401/835903744, -64466909836503771/24167649046528, 1])
1224
sage: L = E.period_lattice()
1225
sage: L.ei()
1226
[5334.003952567705? - 1.964393150436?e-6*I, 5334.003952567705? + 1.964393150436?e-6*I, -10668.25790513541?]
1227
sage: L.elliptic_logarithm(P,prec=100)
1228
0.27656204014107061464076203097
1229
1230
Some complex examples, taken from the paper by Cremona and Thongjunthug::
1231
1232
sage: K.<i> = QuadraticField(-1)
1233
sage: a4 = 9*i-10
1234
sage: a6 = 21-i
1235
sage: E = EllipticCurve([0,0,0,a4,a6])
1236
sage: e1 = 3-2*i; e2 = 1+i; e3 = -4+i
1237
sage: emb = K.embeddings(CC)[1]
1238
sage: L = E.period_lattice(emb)
1239
sage: P = E(2-i,4+2*i)
1240
1241
By default, the output is reduced with respect to the
1242
normalised lattice basis, so that its coordinates with respect
1243
to that basis lie in the interval [0,1)::
1244
1245
sage: z = L.elliptic_logarithm(P,prec=100); z
1246
0.70448375537782208460499649302 - 0.79246725643650979858266018068*I
1247
sage: L.coordinates(z)
1248
(0.46247636364807931766105406092, 0.79497588726808704200760395829)
1249
1250
Using ``reduce=False`` this step can be omitted. In this case
1251
the coordinates are usually in the interval [-0.5,0.5), but
1252
this is not guaranteed. This option is mainly for testing
1253
purposes::
1254
1255
sage: z = L.elliptic_logarithm(P,prec=100, reduce=False); z
1256
0.57002153834710752778063503023 + 0.46476340520469798857457031393*I
1257
sage: L.coordinates(z)
1258
(0.46247636364807931766105406092, -0.20502411273191295799239604171)
1259
1260
The elliptic logs of the 2-torsion points are half-periods::
1261
1262
sage: L.elliptic_logarithm(E(e1,0),prec=100)
1263
0.64607575874356525952487867052 + 0.22379609053909448304176885364*I
1264
sage: L.elliptic_logarithm(E(e2,0),prec=100)
1265
0.71330686725892253793705940192 - 0.40481924028150941053684639367*I
1266
sage: L.elliptic_logarithm(E(e3,0),prec=100)
1267
0.067231108515357278412180731396 - 0.62861533082060389357861524731*I
1268
1269
We check this by doubling and seeing that the resulting
1270
coordinates are integers::
1271
1272
sage: L.coordinates(2*L.elliptic_logarithm(E(e1,0),prec=100))
1273
(1.0000000000000000000000000000, 0.00000000000000000000000000000)
1274
sage: L.coordinates(2*L.elliptic_logarithm(E(e2,0),prec=100))
1275
(1.0000000000000000000000000000, 1.0000000000000000000000000000)
1276
sage: L.coordinates(2*L.elliptic_logarithm(E(e3,0),prec=100))
1277
(0.00000000000000000000000000000, 1.0000000000000000000000000000)
1278
1279
::
1280
1281
sage: a4 = -78*i + 104
1282
sage: a6 = -216*i - 312
1283
sage: E = EllipticCurve([0,0,0,a4,a6])
1284
sage: emb = K.embeddings(CC)[1]
1285
sage: L = E.period_lattice(emb)
1286
sage: P = E(3+2*i,14-7*i)
1287
sage: L.elliptic_logarithm(P)
1288
0.297147783912228 - 0.546125549639461*I
1289
sage: L.coordinates(L.elliptic_logarithm(P))
1290
(0.628653378040238, 0.371417754610223)
1291
sage: e1 = 1+3*i; e2 = -4-12*i; e3=-e1-e2
1292
sage: L.coordinates(L.elliptic_logarithm(E(e1,0)))
1293
(0.500000000000000, 0.500000000000000)
1294
sage: L.coordinates(L.elliptic_logarithm(E(e2,0)))
1295
(1.00000000000000, 0.500000000000000)
1296
sage: L.coordinates(L.elliptic_logarithm(E(e3,0)))
1297
(0.500000000000000, 0.000000000000000)
1298
1299
TESTS (see #10026 and #11767)::
1300
1301
sage: K.<w> = QuadraticField(2)
1302
sage: E = EllipticCurve([ 0, -1, 1, -3*w -4, 3*w + 4 ])
1303
sage: T = E.simon_two_descent()
1304
sage: P,Q = T[2]
1305
sage: embs = K.embeddings(CC)
1306
sage: Lambda = E.period_lattice(embs[0])
1307
sage: Lambda.elliptic_logarithm(P,100)
1308
4.7100131126199672766973600998
1309
sage: R.<x> = QQ[]
1310
sage: K.<a> = NumberField(x^2 + x + 5)
1311
sage: E = EllipticCurve(K, [0,0,1,-3,-5])
1312
sage: P = E([0,a])
1313
sage: Lambda = P.curve().period_lattice(K.embeddings(ComplexField(600))[0])
1314
sage: Lambda.elliptic_logarithm(P, prec=600)
1315
-0.842248166487739393375018008381693990800588864069506187033873183845246233548058477561706400464057832396643843146464236956684557207157300006542470428493573195030603817094900751609464 - 0.571366031453267388121279381354098224265947866751130917440598461117775339240176310729173301979590106474259885638797913383502735083088736326391919063211421189027226502851390118943491*I
1316
sage: K.<a> = QuadraticField(-5)
1317
sage: E = EllipticCurve([1,1,a,a,0])
1318
sage: P = E(0,0)
1319
sage: L = P.curve().period_lattice(K.embeddings(ComplexField())[0])
1320
sage: L.elliptic_logarithm(P, prec=500)
1321
1.17058357737548897849026170185581196033579563441850967539191867385734983296504066660506637438866628981886518901958717288150400849746892393771983141354 - 1.13513899565966043682474529757126359416758251309237866586896869548539516543734207347695898664875799307727928332953834601460994992792519799260968053875*I
1322
sage: L.elliptic_logarithm(P, prec=1000)
1323
1.17058357737548897849026170185581196033579563441850967539191867385734983296504066660506637438866628981886518901958717288150400849746892393771983141354014895386251320571643977497740116710952913769943240797618468987304985625823413440999754037939123032233879499904283600304184828809773650066658885672885 - 1.13513899565966043682474529757126359416758251309237866586896869548539516543734207347695898664875799307727928332953834601460994992792519799260968053875387282656993476491590607092182964878750169490985439873220720963653658829712494879003124071110818175013453207439440032582917366703476398880865439217473*I
1324
"""
1325
if not P.curve() is self.E:
1326
raise ValueError, "Point is on the wrong curve"
1327
if prec is None:
1328
prec = RealField().precision()
1329
if P.is_zero():
1330
return ComplexField(prec)(0)
1331
# Note: using log2(prec) + 3 guard bits is usually enough.
1332
# To avoid computing a logarithm, we use 40 guard bits which
1333
# should be largely enough in practice.
1334
prec2 = prec + 40
1335
R = RealField(prec2)
1336
C = ComplexField(prec2)
1337
pi = R.pi()
1338
e1,e2,e3 = self._ei
1339
a1,a2,a3 = [self.embedding(a) for a in self.E.ainvs()[:3]]
1340
xP, yP = [self.embedding(coord) for coord in P.xy()]
1341
wP = 2*yP+a1*xP+a3
1342
1343
# We treat the case of 2-torsion points separately. (Note
1344
# that Cohen's algorithm does not handle these properly.)
1345
1346
if wP.is_zero(): # 2-torsion treated separately
1347
w1,w2 = self._compute_periods_complex(prec,normalise=False)
1348
if xP==e1:
1349
z = w2/2
1350
else:
1351
if xP==e3:
1352
z = w1/2
1353
else:
1354
z = (w1+w2)/2
1355
if reduce:
1356
z = self.reduce(z)
1357
return z
1358
1359
# NB The first block of code works fine for real embeddings as
1360
# well as complex embeddings. The special code for real
1361
# embeddings uses only real arithmetic in the iteration, and is
1362
# based on Cremona and Thongjunthug.
1363
1364
# An older version, based on Cohen's Algorithm 7.4.8 also uses
1365
# only real arithmetic, and gives different normalisations,
1366
# but also causes problems (see #10026). It is left in but
1367
# commented out below.
1368
1369
if self.real_flag==0: # complex case
1370
1371
a = C((e1-e3).sqrt())
1372
b = C((e1-e2).sqrt())
1373
if (a+b).abs() < (a-b).abs(): b=-b
1374
r = C(((xP-e3)/(xP-e2)).sqrt())
1375
if r.real()<0: r=-r
1376
t = -C(wP)/(2*r*(xP-e2))
1377
# eps controls the end of the loop. Since we aim at a target
1378
# precision of prec bits, eps = 2^(-prec) is enough.
1379
eps = R(1) >> prec
1380
while True:
1381
s = b*r+a
1382
a, b = (a+b)/2, (a*b).sqrt()
1383
if (a+b).abs() < (a-b).abs(): b=-b
1384
r = (a*(r+1)/s).sqrt()
1385
if (r.abs()-1).abs() < eps: break
1386
if r.real()<0: r=-r
1387
t *= r
1388
z = ((a/t).arctan())/a
1389
z = ComplexField(prec)(z)
1390
if reduce:
1391
z = self.reduce(z)
1392
return z
1393
1394
if self.real_flag==-1: # real, connected case
1395
z = C(self._abc[0]) # sqrt(e3-e1)
1396
a, y, b = z.real(), z.imag(), z.abs()
1397
uv = (xP-e1).sqrt()
1398
u, v = uv.real().abs(), uv.imag().abs()
1399
r = (u*a/(u*a+v*y)).sqrt()
1400
t = -r*R(wP)/(2*(u**2+v**2))
1401
on_egg = False
1402
else: # real, disconnected case
1403
a = R(e3-e1).sqrt()
1404
b = R(e3-e2).sqrt()
1405
if (a+b).abs() < (a-b).abs(): b=-b
1406
on_egg = (xP<e3)
1407
if on_egg:
1408
r = a/R(e3-xP).sqrt()
1409
t = r*R(wP)/(2*R(xP-e1))
1410
else:
1411
r = R((xP-e1)/(xP-e2)).sqrt()
1412
t = -R(wP)/(2*r*R(xP-e2))
1413
1414
# eps controls the end of the loop. Since we aim at a target
1415
# precision of prec bits, eps = 2^(-prec) is enough.
1416
eps = R(1) >> prec
1417
while True:
1418
s = b*r+a
1419
a, b = (a+b)/2, (a*b).sqrt()
1420
r = (a*(r+1)/s).sqrt()
1421
if (r-1).abs() < eps: break
1422
t *= r
1423
z = ((a/t).arctan())/a
1424
if on_egg:
1425
w1,w2 = self._compute_periods_real(prec)
1426
z += w2/2
1427
z = ComplexField(prec)(z)
1428
if reduce:
1429
z = self.reduce(z)
1430
return z
1431
1432
1433
# if self.real_flag==-1:
1434
# z = self._abc[1] # sqrt(e3-e2)
1435
# beta = z.norm()
1436
# alpha = 2*(e3-e2).real()
1437
# a = beta.sqrt()*2
1438
# b = (alpha+2*beta).sqrt()
1439
# c = (xP-e3+beta)/(xP-e3).sqrt()
1440
# else:
1441
# on_egg = (xP<e3)
1442
# if on_egg:
1443
# y3 = -(a1*e1+a3)/2
1444
# lam = (yP-y3)/(xP-e1)
1445
# x3 = lam*(lam+a1)-a2-xP-e1
1446
# yP = lam*(xP-x3)-yP-a1*x3-a3
1447
# xP = x3
1448
# wP = 2*yP+a1*xP+a3
1449
# a = self._abc[0] # sqrt(e3-e1)
1450
# b = self._abc[1] # sqrt(e3-e2)
1451
# c = (xP-e1).sqrt()
1452
# # So far the values of a,b,c are algebraic (in AA)
1453
# a = R(a)
1454
# b = R(b)
1455
# c = R(c)
1456
# a,b,c = extended_agm_iteration(a,b,c)
1457
# if self.real_flag==-1:
1458
# z = (a/c).arcsin()
1459
# if wP*((xP-e3)*(xP-e3)-beta*beta) >= 0:
1460
# z = pi - z
1461
# if wP > 0:
1462
# z += pi
1463
# z /= a
1464
# elif self.real_flag==+1:
1465
# z = (a/c).arcsin()/a
1466
# w1 = w2 = None
1467
# if wP > 0:
1468
# if w1 is None:
1469
# w1, w2 = self.basis(prec)
1470
# z = w1 - z
1471
# if on_egg:
1472
# if w2 is None:
1473
# w1, w2 = self.basis(prec)
1474
# z += w2/2
1475
# z = ComplexField(prec)(z)
1476
# if reduce:
1477
# z = self.reduce(z)
1478
# return z
1479
1480
def elliptic_exponential(self, z, to_curve=True):
1481
r"""
1482
Return the elliptic exponential of a complex number.
1483
1484
INPUT:
1485
1486
- ``z`` (complex) -- A complex number (viewed modulo this period lattice).
1487
1488
- ``to_curve`` (bool, default True): see below.
1489
1490
OUTPUT:
1491
1492
- If ``to_curve`` is False, a 2-tuple of real or complex
1493
numbers representing the point `(x,y) = (\wp(z),\wp'(z))`
1494
where `\wp` denotes the Weierstrass `\wp`-function with
1495
respect to this lattice.
1496
1497
- If ``to_curve`` is True, the point `(X,Y) =
1498
(x-b_2/12,y-(a_1(x-b_2/12)-a_3)/2)` as a point in `E(\RR)`
1499
or `E(\CC)`, with `(x,y) = (\wp(z),\wp'(z))` as above, where
1500
`E` is the elliptic curve over `\RR` or `\CC` whose period
1501
lattice this is.
1502
1503
- If the lattice is real and `z` is also real then the output
1504
is a pair of real numbers if ``to_curve`` is True, or a
1505
point in `E(\RR)` if ``to_curve`` is False.
1506
1507
.. note::
1508
1509
The precision is taken from that of the input ``z``.
1510
1511
EXAMPLES::
1512
1513
sage: E = EllipticCurve([1,1,1,-8,6])
1514
sage: P = E(1,-2)
1515
sage: L = E.period_lattice()
1516
sage: z = L(P); z
1517
1.17044757240090
1518
sage: L.elliptic_exponential(z)
1519
(0.999999999999999 : -2.00000000000000 : 1.00000000000000)
1520
sage: _.curve()
1521
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
1522
sage: L.elliptic_exponential(z,to_curve=False)
1523
(1.41666666666667, -1.00000000000000)
1524
sage: z = L(P,prec=201); z
1525
1.17044757240089592298992188482371493504472561677451007994189
1526
sage: L.elliptic_exponential(z)
1527
(1.00000000000000000000000000000000000000000000000000000000000 : -2.00000000000000000000000000000000000000000000000000000000000 : 1.00000000000000000000000000000000000000000000000000000000000)
1528
1529
Examples over number fields::
1530
1531
sage: x = polygen(QQ)
1532
sage: K.<a> = NumberField(x^3-2)
1533
sage: embs = K.embeddings(CC)
1534
sage: E = EllipticCurve('37a')
1535
sage: EK = E.change_ring(K)
1536
sage: Li = [EK.period_lattice(e) for e in embs]
1537
sage: P = EK(-1,-1)
1538
sage: Q = EK(a-1,1-a^2)
1539
sage: zi = [L.elliptic_logarithm(P) for L in Li]
1540
sage: [c.real() for c in Li[0].elliptic_exponential(zi[0])]
1541
[-1.00000000000000, -1.00000000000000, 1.00000000000000]
1542
sage: [c.real() for c in Li[0].elliptic_exponential(zi[1])]
1543
[-1.00000000000000, -1.00000000000000, 1.00000000000000]
1544
sage: [c.real() for c in Li[0].elliptic_exponential(zi[2])]
1545
[-1.00000000000000, -1.00000000000000, 1.00000000000000]
1546
1547
sage: zi = [L.elliptic_logarithm(Q) for L in Li]
1548
sage: Li[0].elliptic_exponential(zi[0])
1549
(-1.62996052494744 - 1.09112363597172*I : 1.79370052598410 - 1.37472963699860*I : 1.00000000000000)
1550
sage: [embs[0](c) for c in Q]
1551
[-1.62996052494744 - 1.09112363597172*I, 1.79370052598410 - 1.37472963699860*I, 1.00000000000000]
1552
sage: Li[1].elliptic_exponential(zi[1])
1553
(-1.62996052494744 + 1.09112363597172*I : 1.79370052598410 + 1.37472963699860*I : 1.00000000000000)
1554
sage: [embs[1](c) for c in Q]
1555
[-1.62996052494744 + 1.09112363597172*I, 1.79370052598410 + 1.37472963699860*I, 1.00000000000000]
1556
sage: [c.real() for c in Li[2].elliptic_exponential(zi[2])]
1557
[0.259921049894873, -0.587401051968199, 1.00000000000000]
1558
sage: [embs[2](c) for c in Q]
1559
[0.259921049894873, -0.587401051968200, 1.00000000000000]
1560
1561
Test to show that #8820 is fixed::
1562
1563
sage: E = EllipticCurve('37a')
1564
sage: K.<a> = QuadraticField(-5)
1565
sage: L = E.change_ring(K).period_lattice(K.places()[0])
1566
sage: L.elliptic_exponential(CDF(.1,.1))
1567
(0.0000142854026029... - 49.9960001066650*I : 249.520141250950 + 250.019855549131*I : 1.00000000000000)
1568
sage: L.elliptic_exponential(CDF(.1,.1), to_curve=False)
1569
(0.0000142854026029... - 49.9960001066650*I, 250.020141250950 + 250.019855549131*I)
1570
1571
`z=0` is treated as a special case::
1572
1573
sage: E = EllipticCurve([1,1,1,-8,6])
1574
sage: L = E.period_lattice()
1575
sage: L.elliptic_exponential(0)
1576
(0.000000000000000 : 1.00000000000000 : 0.000000000000000)
1577
sage: L.elliptic_exponential(0, to_curve=False)
1578
(+infinity, +infinity)
1579
1580
::
1581
1582
sage: E = EllipticCurve('37a')
1583
sage: K.<a> = QuadraticField(-5)
1584
sage: L = E.change_ring(K).period_lattice(K.places()[0])
1585
sage: P = L.elliptic_exponential(0); P
1586
(0.000000000000000 : 1.00000000000000 : 0.000000000000000)
1587
sage: P.parent()
1588
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
1589
1590
Very small `z` are handled properly (see #8820)::
1591
1592
sage: K.<a> = QuadraticField(-1)
1593
sage: E = EllipticCurve([0,0,0,a,0])
1594
sage: L = E.period_lattice(K.complex_embeddings()[0])
1595
sage: L.elliptic_exponential(1e-100)
1596
(0.000000000000000 : 1.00000000000000 : 0.000000000000000)
1597
1598
The elliptic exponential of `z` is returned as (0 : 1 : 0) if
1599
the coordinates of z with respect to the period lattice are
1600
approximately integral::
1601
1602
sage: (100/log(2.0,10))/0.8
1603
415.241011860920
1604
sage: L.elliptic_exponential((RealField(415)(1e-100))).is_zero()
1605
True
1606
sage: L.elliptic_exponential((RealField(420)(1e-100))).is_zero()
1607
False
1608
"""
1609
C = z.parent()
1610
z_is_real = False
1611
if is_RealField(C):
1612
z_is_real = True
1613
C = ComplexField(C.precision())
1614
z = C(z)
1615
else:
1616
if is_ComplexField(C):
1617
z_is_real = z.is_real()
1618
else:
1619
try:
1620
C = ComplexField()
1621
z = C(z)
1622
z_is_real = z.is_real()
1623
except TypeError:
1624
raise TypeError, "%s is not a complex number"%z
1625
prec = C.precision()
1626
1627
# test for the point at infinity:
1628
1629
eps = (C(2)**(-0.8*prec)).real() ## to test integrality w.r.t. lattice within 20%
1630
if all([(t.round()-t).abs() < eps for t in self.coordinates(z)]):
1631
K = z.parent()
1632
if to_curve:
1633
return self.curve().change_ring(K)(0)
1634
else:
1635
return (K('+infinity'),K('+infinity'))
1636
1637
# general number field code (including QQ):
1638
1639
# We do not use PARI's ellztopoint function since it is only
1640
# defined for curves over the reals (note that PARI only
1641
# computes the period lattice basis in that case). But Sage
1642
# can compute the period lattice basis over CC, and then
1643
# PARI's ellwp function works fine.
1644
1645
# NB converting the PARI values to Sage values might land up
1646
# in real/complex fields of spuriously higher precision than
1647
# the input, since PARI's precision is in word-size chunks.
1648
# So we force the results back into the real/complex fields of
1649
# the same precision as the input.
1650
1651
from sage.libs.all import pari
1652
1653
x,y = pari(self.basis(prec=prec)).ellwp(pari(z),flag=1)
1654
x,y = [C(t) for t in (x,y)]
1655
1656
if self.real_flag and z_is_real:
1657
x = x.real()
1658
y = y.real()
1659
1660
if to_curve:
1661
a1,a2,a3,a4,a6 = [self.embedding(a) for a in self.E.ainvs()]
1662
b2 = self.embedding(self.E.b2())
1663
x = x - (b2/12)
1664
y = y - (a1*x+a3)/2
1665
K = x.parent()
1666
EK = EllipticCurve(K,[a1,a2,a3,a4,a6])
1667
return EK.point((x,y,K(1)), check=False)
1668
else:
1669
return (x,y)
1670
1671
def reduce_tau(tau):
1672
r"""
1673
Transform a point in the upper half plane to the fundamental region.
1674
1675
INPUT:
1676
1677
- ``tau`` (complex) -- a complex number with positive imaginary part
1678
1679
OUTPUT:
1680
1681
(tuple) `(\tau',[a,b,c,d])` where `a,b,c,d` are integers such that
1682
1683
- `ad-bc=1`;
1684
- `\tau`=(a\tau+b)/(c\tau+d)`;
1685
- `|\tau'|\ge1`;
1686
- `|\Re(\tau')|\le\frac{1}{2}`.
1687
1688
EXAMPLES::
1689
1690
sage: from sage.schemes.elliptic_curves.period_lattice import reduce_tau
1691
sage: reduce_tau(CC(1.23,3.45))
1692
(0.230000000000000 + 3.45000000000000*I, [1, -1, 0, 1])
1693
sage: reduce_tau(CC(1.23,0.0345))
1694
(-0.463960069171512 + 1.35591888067914*I, [-5, 6, 4, -5])
1695
sage: reduce_tau(CC(1.23,0.0000345))
1696
(0.130000000001761 + 2.89855072463768*I, [13, -16, 100, -123])
1697
"""
1698
assert tau.imag()>0
1699
tau_orig = tau
1700
a, b = ZZ(1), ZZ(0)
1701
c, d = b, a
1702
k = tau.real().round()
1703
tau -= k
1704
a -= k*c
1705
b -= k*d
1706
while tau.abs()<0.999:
1707
tau = -1/tau
1708
a, b, c, d = c, d, -a, -b
1709
k = tau.real().round()
1710
tau -= k
1711
a -= k*c
1712
b -= k*d
1713
assert a*d-b*c==1
1714
assert tau.abs()>=0.999 and tau.real().abs() <= 0.5
1715
return tau, [a,b,c,d]
1716
1717
def normalise_periods(w1,w2):
1718
r"""
1719
Normalise the period basis `(w_1,w_2)` so that `w_1/w_2` is in the fundamental region.
1720
1721
INPUT:
1722
1723
- ``w1,w2`` (complex) -- two complex numbers with non-real ratio
1724
1725
OUTPUT:
1726
1727
(tuple) `((\omega_1',\omega_2'),[a,b,c,d])` where `a,b,c,d` are
1728
integers such that
1729
1730
- `ad-bc=\pm1`;
1731
- `(\omega_1',\omega_2') = (a\omega_1+b\omega_2,c\omega_1+d\omega_2)`;
1732
- `\tau=\omega_1'/\omega_2'` is in the upper half plane;
1733
- `|\tau|\ge1` and `|\Re(\tau)|\le\frac{1}{2}`.
1734
1735
EXAMPLES::
1736
1737
sage: from sage.schemes.elliptic_curves.period_lattice import reduce_tau, normalise_periods
1738
sage: w1 = CC(1.234, 3.456)
1739
sage: w2 = CC(1.234, 3.456000001)
1740
sage: w1/w2 # in lower half plane!
1741
0.999999999743367 - 9.16334785827644e-11*I
1742
sage: w1w2, abcd = normalise_periods(w1,w2)
1743
sage: a,b,c,d = abcd
1744
sage: w1w2 == (a*w1+b*w2, c*w1+d*w2)
1745
True
1746
sage: w1w2[0]/w1w2[1]
1747
1.23400010389203e9*I
1748
sage: a*d-b*c # note change of orientation
1749
-1
1750
1751
"""
1752
tau = w1/w2
1753
s = +1
1754
if tau.imag()<0:
1755
w2 = -w2
1756
tau = -tau
1757
s = -1
1758
tau, abcd = reduce_tau(tau)
1759
a, b, c, d = abcd
1760
if s<0:
1761
abcd = (a,-b,c,-d)
1762
return (a*w1+b*w2,c*w1+d*w2), abcd
1763
1764
1765
def extended_agm_iteration(a,b,c):
1766
r"""
1767
Internal function for the extended AGM used in elliptic logarithm computation.
1768
INPUT:
1769
1770
- ``a``, ``b``, ``c`` (real or complex) -- three real or complex numbers.
1771
1772
OUTPUT:
1773
1774
(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)`.
1775
1776
EXAMPLES::
1777
1778
sage: from sage.schemes.elliptic_curves.period_lattice import extended_agm_iteration
1779
sage: extended_agm_iteration(RR(1),RR(2),RR(3))
1780
(1.45679103104691, 1.45679103104691, 3.21245294970054)
1781
sage: extended_agm_iteration(CC(1,2),CC(2,3),CC(3,4))
1782
(1.46242448156430 + 2.47791311676267*I,
1783
1.46242448156430 + 2.47791311676267*I,
1784
3.22202144343535 + 4.28383734262540*I)
1785
1786
TESTS::
1787
1788
sage: extended_agm_iteration(1,2,3)
1789
Traceback (most recent call last):
1790
...
1791
ValueError: values must be real or complex numbers
1792
1793
"""
1794
if not isinstance(a, (RealNumber,ComplexNumber)):
1795
raise ValueError, "values must be real or complex numbers"
1796
eps = a.parent().one().real()>>(a.parent().precision()-10)
1797
while True:
1798
a1 = (a + b)/2
1799
b1 = (a*b).sqrt()
1800
delta = (b**2 - a**2)/c**2
1801
f = (1 + (1 + delta).sqrt())/2
1802
if (f.abs()-1).abs() < eps:
1803
return a,b,c
1804
c*=f
1805
a,b = a1,b1
1806
1807