Book a Demo!
CoCalc Logo Icon
StoreFeaturesDocsShareSupportNewsAboutPoliciesSign UpSign In
sagemath
GitHub Repository: sagemath/sagelib
Path: blob/master/sage/schemes/elliptic_curves/ell_local_data.py
4107 views
1
# -*- coding: utf-8 -*-
2
r"""
3
Local data for elliptic curves over number fields
4
5
Let `E` be an elliptic curve over a number field `K` (including `\QQ`).
6
There are several local invariants at a finite place `v` that
7
can be computed via Tate's algorithm (see [Sil2] IV.9.4 or [Ta]).
8
9
These include the type of reduction (good, additive, multiplicative),
10
a minimal equation of `E` over `K_v`,
11
the Tamagawa number `c_v`, defined to be the index `[E(K_v):E^0(K_v)]`
12
of the points with good reduction among the local points, and the
13
exponent of the conductor `f_v`.
14
15
The functions in this file will typically be called by using ``local_data``.
16
17
EXAMPLES::
18
19
sage: K.<i> = NumberField(x^2+1)
20
sage: E = EllipticCurve([(2+i)^2,(2+i)^7])
21
sage: pp = K.fractional_ideal(2+i)
22
sage: da = E.local_data(pp)
23
sage: da.has_bad_reduction()
24
True
25
sage: da.has_multiplicative_reduction()
26
False
27
sage: da.kodaira_symbol()
28
I0*
29
sage: da.tamagawa_number()
30
4
31
sage: da.minimal_model()
32
Elliptic Curve defined by y^2 = x^3 + (4*i+3)*x + (-29*i-278) over Number Field in i with defining polynomial x^2 + 1
33
34
An example to show how the Neron model can change as one extends the field::
35
36
sage: E = EllipticCurve([0,-1])
37
sage: E.local_data(2)
38
Local data at Principal ideal (2) of Integer Ring:
39
Reduction type: bad additive
40
Local minimal model: Elliptic Curve defined by y^2 = x^3 - 1 over Rational Field
41
Minimal discriminant valuation: 4
42
Conductor exponent: 4
43
Kodaira Symbol: II
44
Tamagawa Number: 1
45
46
sage: EK = E.base_extend(K)
47
sage: EK.local_data(1+i)
48
Local data at Fractional ideal (i + 1):
49
Reduction type: bad additive
50
Local minimal model: Elliptic Curve defined by y^2 = x^3 + (-1) over Number Field in i with defining polynomial x^2 + 1
51
Minimal discriminant valuation: 8
52
Conductor exponent: 2
53
Kodaira Symbol: IV*
54
Tamagawa Number: 3
55
56
Or how the minimal equation changes::
57
58
sage: E = EllipticCurve([0,8])
59
sage: E.is_minimal()
60
True
61
sage: EK = E.base_extend(K)
62
sage: da = EK.local_data(1+i)
63
sage: da.minimal_model()
64
Elliptic Curve defined by y^2 = x^3 + i over Number Field in i with defining polynomial x^2 + 1
65
66
REFERENCES:
67
68
- [Sil2] Silverman, Joseph H., Advanced topics in the arithmetic of elliptic curves.
69
Graduate Texts in Mathematics, 151. Springer-Verlag, New York, 1994.
70
71
- [Ta] Tate, John, Algorithm for determining the type of a singular fiber in an elliptic pencil.
72
Modular functions of one variable, IV, pp. 33--52. Lecture Notes in Math., Vol. 476,
73
Springer, Berlin, 1975.
74
75
AUTHORS:
76
77
- John Cremona: First version 2008-09-21 (refactoring code from
78
``ell_number_field.py`` and ``ell_rational_field.py``)
79
80
- Chris Wuthrich: more documentation 2010-01
81
82
"""
83
84
#*****************************************************************************
85
# Copyright (C) 2005 William Stein <[email protected]>
86
#
87
# Distributed under the terms of the GNU General Public License (GPL)
88
#
89
# This code is distributed in the hope that it will be useful,
90
# but WITHOUT ANY WARRANTY; without even the implied warranty of
91
# MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the GNU
92
# General Public License for more details.
93
#
94
# The full text of the GPL is available at:
95
#
96
# http://www.gnu.org/licenses/
97
#*****************************************************************************
98
99
100
from sage.structure.sage_object import SageObject
101
from sage.misc.misc import verbose
102
103
from sage.rings.all import PolynomialRing, QQ, ZZ, Integer, is_Ideal, is_NumberFieldElement, is_NumberFieldFractionalIdeal, is_NumberField
104
from constructor import EllipticCurve
105
from kodaira_symbol import KodairaSymbol
106
107
class EllipticCurveLocalData(SageObject):
108
r"""
109
The class for the local reduction data of an elliptic curve.
110
111
Currently supported are elliptic curves defined over `\QQ`, and
112
elliptic curves defined over a number field, at an arbitrary prime
113
or prime ideal.
114
115
INPUT:
116
117
- ``E`` -- an elliptic curve defined over a number field, or `\QQ`.
118
119
- ``P`` -- a prime ideal of the field, or a prime integer if the field is `\QQ`.
120
121
- ``proof`` (bool)-- if True, only use provably correct
122
methods (default controlled by global proof module). Note
123
that the proof module is number_field, not elliptic_curves,
124
since the functions that actually need the flag are in
125
number fields.
126
127
- ``algorithm`` (string, default: "pari") -- Ignored unless the
128
base field is `\QQ`. If "pari", use the PARI C-library
129
``ellglobalred`` implementation of Tate's algorithm over
130
`\QQ`. If "generic", use the general number field
131
implementation.
132
133
.. note::
134
135
This function is not normally called directly by users, who
136
may access the data via methods of the EllipticCurve
137
classes.
138
139
EXAMPLES::
140
141
sage: from sage.schemes.elliptic_curves.ell_local_data import EllipticCurveLocalData
142
sage: E = EllipticCurve('14a1')
143
sage: EllipticCurveLocalData(E,2)
144
Local data at Principal ideal (2) of Integer Ring:
145
Reduction type: bad non-split multiplicative
146
Local minimal model: Elliptic Curve defined by y^2 + x*y + y = x^3 + 4*x - 6 over Rational Field
147
Minimal discriminant valuation: 6
148
Conductor exponent: 1
149
Kodaira Symbol: I6
150
Tamagawa Number: 2
151
152
"""
153
154
def __init__(self, E, P, proof=None, algorithm="pari"):
155
r"""
156
Initializes the reduction data for the elliptic curve `E` at the prime `P`.
157
158
INPUT:
159
160
- ``E`` -- an elliptic curve defined over a number field, or `\QQ`.
161
162
- ``P`` -- a prime ideal of the field, or a prime integer if the field is `\QQ`.
163
164
- ``proof`` (bool)-- if True, only use provably correct
165
methods (default controlled by global proof module). Note
166
that the proof module is number_field, not elliptic_curves,
167
since the functions that actually need the flag are in
168
number fields.
169
170
- ``algorithm`` (string, default: "pari") -- Ignored unless the
171
base field is `\QQ`. If "pari", use the PARI C-library
172
``ellglobalred`` implementation of Tate's algorithm over
173
`\QQ`. If "generic", use the general number field
174
implementation.
175
176
.. note::
177
178
This function is not normally called directly by users, who
179
may access the data via methods of the EllipticCurve
180
classes.
181
182
EXAMPLES::
183
184
sage: from sage.schemes.elliptic_curves.ell_local_data import EllipticCurveLocalData
185
sage: E = EllipticCurve('14a1')
186
sage: EllipticCurveLocalData(E,2)
187
Local data at Principal ideal (2) of Integer Ring:
188
Reduction type: bad non-split multiplicative
189
Local minimal model: Elliptic Curve defined by y^2 + x*y + y = x^3 + 4*x - 6 over Rational Field
190
Minimal discriminant valuation: 6
191
Conductor exponent: 1
192
Kodaira Symbol: I6
193
Tamagawa Number: 2
194
195
::
196
197
sage: EllipticCurveLocalData(E,2,algorithm="generic")
198
Local data at Principal ideal (2) of Integer Ring:
199
Reduction type: bad non-split multiplicative
200
Local minimal model: Elliptic Curve defined by y^2 + x*y + y = x^3 + 4*x - 6 over Rational Field
201
Minimal discriminant valuation: 6
202
Conductor exponent: 1
203
Kodaira Symbol: I6
204
Tamagawa Number: 2
205
206
::
207
208
sage: EllipticCurveLocalData(E,2,algorithm="pari")
209
Local data at Principal ideal (2) of Integer Ring:
210
Reduction type: bad non-split multiplicative
211
Local minimal model: Elliptic Curve defined by y^2 + x*y + y = x^3 + 4*x - 6 over Rational Field
212
Minimal discriminant valuation: 6
213
Conductor exponent: 1
214
Kodaira Symbol: I6
215
Tamagawa Number: 2
216
217
::
218
219
sage: EllipticCurveLocalData(E,2,algorithm="unknown")
220
Traceback (most recent call last):
221
...
222
ValueError: algorithm must be one of 'pari', 'generic'
223
224
::
225
226
sage: EllipticCurveLocalData(E,3)
227
Local data at Principal ideal (3) of Integer Ring:
228
Reduction type: good
229
Local minimal model: Elliptic Curve defined by y^2 + x*y + y = x^3 + 4*x - 6 over Rational Field
230
Minimal discriminant valuation: 0
231
Conductor exponent: 0
232
Kodaira Symbol: I0
233
Tamagawa Number: 1
234
235
::
236
237
sage: EllipticCurveLocalData(E,7)
238
Local data at Principal ideal (7) of Integer Ring:
239
Reduction type: bad split multiplicative
240
Local minimal model: Elliptic Curve defined by y^2 + x*y + y = x^3 + 4*x - 6 over Rational Field
241
Minimal discriminant valuation: 3
242
Conductor exponent: 1
243
Kodaira Symbol: I3
244
Tamagawa Number: 3
245
"""
246
self._curve = E
247
K = E.base_field()
248
p = check_prime(K,P) # error handling done in that function
249
if algorithm != "pari" and algorithm != "generic":
250
raise ValueError, "algorithm must be one of 'pari', 'generic'"
251
252
self._reduction_type = None
253
if K is QQ:
254
self._prime = ZZ.ideal(p)
255
else:
256
self._prime = p
257
258
if algorithm=="pari" and K is QQ:
259
Eint = E.integral_model()
260
data = Eint.pari_curve().elllocalred(p)
261
self._fp = data[0].python()
262
self._KS = KodairaSymbol(data[1].python())
263
self._cp = data[3].python()
264
# We use a global minimal model since we can:
265
self._Emin_reduced = Eint.minimal_model()
266
self._val_disc = self._Emin_reduced.discriminant().valuation(p)
267
if self._fp>0:
268
self._reduction_type = Eint.ap(p) # = 0,-1 or +1
269
else:
270
self._Emin, ch, self._val_disc, self._fp, self._KS, self._cp, self._split = self._tate(proof)
271
if self._fp>0:
272
if self._Emin.c4().valuation(p)>0:
273
self._reduction_type = 0
274
elif self._split:
275
self._reduction_type = +1
276
else:
277
self._reduction_type = -1
278
279
def __repr__(self):
280
r"""
281
Returns the string representation of this reduction data.
282
283
EXAMPLES::
284
285
sage: from sage.schemes.elliptic_curves.ell_local_data import EllipticCurveLocalData
286
sage: E = EllipticCurve('14a1')
287
sage: EllipticCurveLocalData(E,2).__repr__()
288
'Local data at Principal ideal (2) of Integer Ring:\nReduction type: bad non-split multiplicative\nLocal minimal model: Elliptic Curve defined by y^2 + x*y + y = x^3 + 4*x - 6 over Rational Field\nMinimal discriminant valuation: 6\nConductor exponent: 1\nKodaira Symbol: I6\nTamagawa Number: 2'
289
sage: EllipticCurveLocalData(E,3).__repr__()
290
'Local data at Principal ideal (3) of Integer Ring:\nReduction type: good\nLocal minimal model: Elliptic Curve defined by y^2 + x*y + y = x^3 + 4*x - 6 over Rational Field\nMinimal discriminant valuation: 0\nConductor exponent: 0\nKodaira Symbol: I0\nTamagawa Number: 1'
291
sage: EllipticCurveLocalData(E,7).__repr__()
292
'Local data at Principal ideal (7) of Integer Ring:\nReduction type: bad split multiplicative\nLocal minimal model: Elliptic Curve defined by y^2 + x*y + y = x^3 + 4*x - 6 over Rational Field\nMinimal discriminant valuation: 3\nConductor exponent: 1\nKodaira Symbol: I3\nTamagawa Number: 3'
293
"""
294
red_type = "good"
295
if not self._reduction_type is None:
296
red_type = ["bad non-split multiplicative","bad additive","bad split multiplicative"][1+self._reduction_type]
297
return "Local data at %s:\nReduction type: %s\nLocal minimal model: %s\nMinimal discriminant valuation: %s\nConductor exponent: %s\nKodaira Symbol: %s\nTamagawa Number: %s"%(self._prime,red_type,self.minimal_model(),self._val_disc,self._fp,self._KS,self._cp)
298
299
def minimal_model(self, reduce=True):
300
"""
301
Return the (local) minimal model from this local reduction data.
302
303
INPUT:
304
305
- ``reduce`` -- (default: True) if set to True the EC returned
306
by Tate's algorithm will be
307
"reduced" as specified in _reduce_model() for curves over
308
number fields.
309
310
EXAMPLES::
311
312
sage: from sage.schemes.elliptic_curves.ell_local_data import EllipticCurveLocalData
313
sage: E = EllipticCurve([0,0,0,0,64]); E
314
Elliptic Curve defined by y^2 = x^3 + 64 over Rational Field
315
sage: data = EllipticCurveLocalData(E,2)
316
sage: data.minimal_model()
317
Elliptic Curve defined by y^2 = x^3 + 1 over Rational Field
318
sage: data.minimal_model() == E.local_minimal_model(2)
319
True
320
321
To demonstrate the behaviour of the parameter ``reduce``::
322
323
sage: K.<a> = NumberField(x^3+x+1)
324
sage: E = EllipticCurve(K, [0, 0, a, 0, 1])
325
sage: E.local_data(K.ideal(a-1)).minimal_model()
326
Elliptic Curve defined by y^2 + a*y = x^3 + 1 over Number Field in a with defining polynomial x^3 + x + 1
327
sage: E.local_data(K.ideal(a-1)).minimal_model(reduce=False)
328
Elliptic Curve defined by y^2 + (a+2)*y = x^3 + 3*x^2 + 3*x + (-a+1) over Number Field in a with defining polynomial x^3 + x + 1
329
330
sage: E = EllipticCurve([2, 1, 0, -2, -1])
331
sage: E.local_data(ZZ.ideal(2), algorithm="generic").minimal_model(reduce=False)
332
Elliptic Curve defined by y^2 + 2*x*y + 2*y = x^3 + x^2 - 4*x - 2 over Rational Field
333
sage: E.local_data(ZZ.ideal(2), algorithm="pari").minimal_model(reduce=False)
334
Traceback (most recent call last):
335
...
336
ValueError: the argument reduce must not be False if algorithm=pari is used
337
sage: E.local_data(ZZ.ideal(2), algorithm="generic").minimal_model()
338
Elliptic Curve defined by y^2 = x^3 - x^2 - 3*x + 2 over Rational Field
339
sage: E.local_data(ZZ.ideal(2), algorithm="pari").minimal_model()
340
Elliptic Curve defined by y^2 = x^3 - x^2 - 3*x + 2 over Rational Field
341
"""
342
if reduce:
343
try:
344
return self._Emin_reduced
345
except AttributeError:
346
pass
347
self._Emin_reduced = self._Emin._reduce_model()
348
return self._Emin_reduced
349
else:
350
try:
351
return self._Emin
352
except AttributeError:
353
raise ValueError, "the argument reduce must not be False if algorithm=pari is used"
354
355
def prime(self):
356
"""
357
Return the prime ideal associated with this local reduction data.
358
359
EXAMPLES::
360
361
sage: from sage.schemes.elliptic_curves.ell_local_data import EllipticCurveLocalData
362
sage: E = EllipticCurve([0,0,0,0,64]); E
363
Elliptic Curve defined by y^2 = x^3 + 64 over Rational Field
364
sage: data = EllipticCurveLocalData(E,2)
365
sage: data.prime()
366
Principal ideal (2) of Integer Ring
367
"""
368
return self._prime
369
370
def conductor_valuation(self):
371
"""
372
Return the valuation of the conductor from this local reduction data.
373
374
EXAMPLES::
375
376
sage: from sage.schemes.elliptic_curves.ell_local_data import EllipticCurveLocalData
377
sage: E = EllipticCurve([0,0,0,0,64]); E
378
Elliptic Curve defined by y^2 = x^3 + 64 over Rational Field
379
sage: data = EllipticCurveLocalData(E,2)
380
sage: data.conductor_valuation()
381
2
382
"""
383
return self._fp
384
385
def kodaira_symbol(self):
386
r"""
387
Return the Kodaira symbol from this local reduction data.
388
389
EXAMPLES::
390
391
sage: from sage.schemes.elliptic_curves.ell_local_data import EllipticCurveLocalData
392
sage: E = EllipticCurve([0,0,0,0,64]); E
393
Elliptic Curve defined by y^2 = x^3 + 64 over Rational Field
394
sage: data = EllipticCurveLocalData(E,2)
395
sage: data.kodaira_symbol()
396
IV
397
"""
398
return self._KS
399
400
def tamagawa_number(self):
401
r"""
402
Return the Tamagawa number from this local reduction data.
403
404
This is the index `[E(K_v):E^0(K_v)]`.
405
406
EXAMPLES::
407
408
sage: from sage.schemes.elliptic_curves.ell_local_data import EllipticCurveLocalData
409
sage: E = EllipticCurve([0,0,0,0,64]); E
410
Elliptic Curve defined by y^2 = x^3 + 64 over Rational Field
411
sage: data = EllipticCurveLocalData(E,2)
412
sage: data.tamagawa_number()
413
3
414
"""
415
return self._cp
416
417
def tamagawa_exponent(self):
418
r"""
419
Return the Tamagawa index from this local reduction data.
420
421
This is the exponent of `E(K_v)/E^0(K_v)`; in most cases it is
422
the same as the Tamagawa index.
423
424
EXAMPLES::
425
426
sage: from sage.schemes.elliptic_curves.ell_local_data import EllipticCurveLocalData
427
sage: E = EllipticCurve('816a1')
428
sage: data = EllipticCurveLocalData(E,2)
429
sage: data.kodaira_symbol()
430
I2*
431
sage: data.tamagawa_number()
432
4
433
sage: data.tamagawa_exponent()
434
2
435
436
sage: E = EllipticCurve('200c4')
437
sage: data = EllipticCurveLocalData(E,5)
438
sage: data.kodaira_symbol()
439
I4*
440
sage: data.tamagawa_number()
441
4
442
sage: data.tamagawa_exponent()
443
2
444
"""
445
cp = self._cp
446
if cp!=4:
447
return cp
448
ks = self._KS
449
if ks._roman==1 and ks._n%2==0 and ks._starred:
450
return 2
451
return 4
452
453
def bad_reduction_type(self):
454
r"""
455
Return the type of bad reduction of this reduction data.
456
457
OUTPUT:
458
459
(int or ``None``):
460
461
- +1 for split multiplicative reduction
462
- -1 for non-split multiplicative reduction
463
- 0 for additive reduction
464
- ``None`` for good reduction
465
466
EXAMPLES::
467
468
sage: E=EllipticCurve('14a1')
469
sage: [(p,E.local_data(p).bad_reduction_type()) for p in prime_range(15)]
470
[(2, -1), (3, None), (5, None), (7, 1), (11, None), (13, None)]
471
472
sage: K.<a>=NumberField(x^3-2)
473
sage: P17a, P17b = [P for P,e in K.factor(17)]
474
sage: E = EllipticCurve([0,0,0,0,2*a+1])
475
sage: [(p,E.local_data(p).bad_reduction_type()) for p in [P17a,P17b]]
476
[(Fractional ideal (4*a^2 - 2*a + 1), None), (Fractional ideal (2*a + 1), 0)]
477
"""
478
return self._reduction_type
479
480
def has_good_reduction(self):
481
r"""
482
Return True if there is good reduction.
483
484
EXAMPLES::
485
486
sage: E = EllipticCurve('14a1')
487
sage: [(p,E.local_data(p).has_good_reduction()) for p in prime_range(15)]
488
[(2, False), (3, True), (5, True), (7, False), (11, True), (13, True)]
489
490
sage: K.<a> = NumberField(x^3-2)
491
sage: P17a, P17b = [P for P,e in K.factor(17)]
492
sage: E = EllipticCurve([0,0,0,0,2*a+1])
493
sage: [(p,E.local_data(p).has_good_reduction()) for p in [P17a,P17b]]
494
[(Fractional ideal (4*a^2 - 2*a + 1), True),
495
(Fractional ideal (2*a + 1), False)]
496
"""
497
return self._reduction_type is None
498
499
def has_bad_reduction(self):
500
r"""
501
Return True if there is bad reduction.
502
503
EXAMPLES::
504
505
sage: E = EllipticCurve('14a1')
506
sage: [(p,E.local_data(p).has_bad_reduction()) for p in prime_range(15)]
507
[(2, True), (3, False), (5, False), (7, True), (11, False), (13, False)]
508
509
::
510
511
sage: K.<a> = NumberField(x^3-2)
512
sage: P17a, P17b = [P for P,e in K.factor(17)]
513
sage: E = EllipticCurve([0,0,0,0,2*a+1])
514
sage: [(p,E.local_data(p).has_bad_reduction()) for p in [P17a,P17b]]
515
[(Fractional ideal (4*a^2 - 2*a + 1), False),
516
(Fractional ideal (2*a + 1), True)]
517
"""
518
return not self._reduction_type is None
519
520
def has_multiplicative_reduction(self):
521
r"""
522
Return True if there is multiplicative reduction.
523
524
.. note::
525
526
See also ``has_split_multiplicative_reduction()`` and
527
``has_nonsplit_multiplicative_reduction()``.
528
529
EXAMPLES::
530
531
sage: E = EllipticCurve('14a1')
532
sage: [(p,E.local_data(p).has_multiplicative_reduction()) for p in prime_range(15)]
533
[(2, True), (3, False), (5, False), (7, True), (11, False), (13, False)]
534
535
::
536
537
sage: K.<a> = NumberField(x^3-2)
538
sage: P17a, P17b = [P for P,e in K.factor(17)]
539
sage: E = EllipticCurve([0,0,0,0,2*a+1])
540
sage: [(p,E.local_data(p).has_multiplicative_reduction()) for p in [P17a,P17b]]
541
[(Fractional ideal (4*a^2 - 2*a + 1), False), (Fractional ideal (2*a + 1), False)]
542
"""
543
return self._reduction_type in (-1,+1)
544
545
def has_split_multiplicative_reduction(self):
546
r"""
547
Return True if there is split multiplicative reduction.
548
549
EXAMPLES::
550
551
sage: E = EllipticCurve('14a1')
552
sage: [(p,E.local_data(p).has_split_multiplicative_reduction()) for p in prime_range(15)]
553
[(2, False), (3, False), (5, False), (7, True), (11, False), (13, False)]
554
555
::
556
557
sage: K.<a> = NumberField(x^3-2)
558
sage: P17a, P17b = [P for P,e in K.factor(17)]
559
sage: E = EllipticCurve([0,0,0,0,2*a+1])
560
sage: [(p,E.local_data(p).has_split_multiplicative_reduction()) for p in [P17a,P17b]]
561
[(Fractional ideal (4*a^2 - 2*a + 1), False),
562
(Fractional ideal (2*a + 1), False)]
563
"""
564
return self._reduction_type == +1
565
566
def has_nonsplit_multiplicative_reduction(self):
567
r"""
568
Return True if there is non-split multiplicative reduction.
569
570
EXAMPLES::
571
572
sage: E = EllipticCurve('14a1')
573
sage: [(p,E.local_data(p).has_nonsplit_multiplicative_reduction()) for p in prime_range(15)]
574
[(2, True), (3, False), (5, False), (7, False), (11, False), (13, False)]
575
576
::
577
578
sage: K.<a> = NumberField(x^3-2)
579
sage: P17a, P17b = [P for P,e in K.factor(17)]
580
sage: E = EllipticCurve([0,0,0,0,2*a+1])
581
sage: [(p,E.local_data(p).has_nonsplit_multiplicative_reduction()) for p in [P17a,P17b]]
582
[(Fractional ideal (4*a^2 - 2*a + 1), False), (Fractional ideal (2*a + 1), False)]
583
"""
584
return self._reduction_type == -1
585
586
def has_additive_reduction(self):
587
r"""
588
Return True if there is additive reduction.
589
590
EXAMPLES::
591
592
sage: E = EllipticCurve('27a1')
593
sage: [(p,E.local_data(p).has_additive_reduction()) for p in prime_range(15)]
594
[(2, False), (3, True), (5, False), (7, False), (11, False), (13, False)]
595
596
::
597
598
sage: K.<a> = NumberField(x^3-2)
599
sage: P17a, P17b = [P for P,e in K.factor(17)]
600
sage: E = EllipticCurve([0,0,0,0,2*a+1])
601
sage: [(p,E.local_data(p).has_additive_reduction()) for p in [P17a,P17b]]
602
[(Fractional ideal (4*a^2 - 2*a + 1), False),
603
(Fractional ideal (2*a + 1), True)]
604
"""
605
return self._reduction_type == 0
606
607
def _tate(self, proof = None):
608
r"""
609
Tate's algorithm for an elliptic curve over a number field.
610
611
Computes both local reduction data at a prime ideal and a
612
local minimal model.
613
614
The model is not required to be integral on input. If `P` is
615
principal, uses a generator as uniformizer, so it will not
616
affect integrality or minimality at other primes. If `P` is not
617
principal, the minimal model returned will preserve
618
integrality at other primes, but not minimality.
619
620
.. note::
621
622
Called only by ``EllipticCurveLocalData.__init__()``.
623
624
OUTPUT:
625
626
(tuple) ``(Emin, p, val_disc, fp, KS, cp)`` where:
627
628
- ``Emin`` (EllipticCurve) is a model (integral and) minimal at P
629
- ``p`` (int) is the residue characteristic
630
- ``val_disc`` (int) is the valuation of the local minimal discriminant
631
- ``fp`` (int) is the valuation of the conductor
632
- ``KS`` (string) is the Kodaira symbol
633
- ``cp`` (int) is the Tamagawa number
634
635
636
EXAMPLES (this raised a type error in sage prior to 4.4.4, see ticket #7930) ::
637
638
sage: E = EllipticCurve('99d1')
639
640
sage: R.<X> = QQ[]
641
sage: K.<t> = NumberField(X^3 + X^2 - 2*X - 1)
642
sage: L.<s> = NumberField(X^3 + X^2 - 36*X - 4)
643
644
sage: EK = E.base_extend(K)
645
sage: toK = EK.torsion_order()
646
sage: da = EK.local_data() # indirect doctest
647
648
sage: EL = E.base_extend(L)
649
sage: da = EL.local_data() # indirect doctest
650
651
EXAMPLES:
652
653
The following example shows that the bug at #9324 is fixed::
654
655
sage: K.<a> = NumberField(x^2-x+6)
656
sage: E = EllipticCurve([0,0,0,-53160*a-43995,-5067640*a+19402006])
657
sage: E.conductor() # indirect doctest
658
Fractional ideal (18, 6*a)
659
660
The following example shows that the bug at #9417 is fixed::
661
662
sage: K.<a> = NumberField(x^2+18*x+1)
663
sage: E = EllipticCurve(K, [0, -36, 0, 320, 0])
664
sage: E.tamagawa_number(K.ideal(2))
665
4
666
667
"""
668
E = self._curve
669
P = self._prime
670
K = E.base_ring()
671
OK = K.maximal_order()
672
t = verbose("Running Tate's algorithm with P = %s"%P, level=1)
673
F = OK.residue_field(P)
674
p = F.characteristic()
675
676
# In case P is not principal we mostly use a uniformiser which
677
# is globally integral (with positive valuation at some other
678
# primes); for this to work, it is essential that we can
679
# reduce (mod P) elements of K which are not integral (but are
680
# P-integral). However, if the model is non-minimal and we
681
# end up dividing a_i by pi^i then at that point we use a
682
# uniformiser pi which has non-positive valuation at all other
683
# primes, so that we can divide by it without losing
684
# integrality at other primes.
685
686
principal_flag = P.is_principal()
687
if principal_flag:
688
pi = P.gens_reduced()[0]
689
verbose("P is principal, generator pi = %s"%pi, t, 1)
690
else:
691
pi = K.uniformizer(P, 'positive')
692
verbose("P is not principal, uniformizer pi = %s"%pi, t, 1)
693
pi2 = pi*pi; pi3 = pi*pi2; pi4 = pi*pi3
694
pi_neg = None
695
prime = pi if K is QQ else P
696
697
pval = lambda x: x.valuation(prime)
698
pdiv = lambda x: x.is_zero() or pval(x) > 0
699
# Since ResidueField is cached in a way that
700
# does not care much about embeddings of number
701
# fields, it can happen that F.p.ring() is different
702
# from K. This is a problem: If F.p.ring() has no
703
# embedding but K has, then there is no coercion
704
# from F.p.ring().maximal_order() to K. But it is
705
# no problem to do an explicit conversion in that
706
# case (Simon King, trac ticket #8800).
707
708
from sage.categories.pushout import pushout, CoercionException
709
try:
710
if hasattr(F.p.ring(), 'maximal_order'): # it is not ZZ
711
_tmp_ = pushout(F.p.ring().maximal_order(),K)
712
pinv = lambda x: F.lift(~F(x))
713
proot = lambda x,e: F.lift(F(x).nth_root(e, extend = False, all = True)[0])
714
preduce = lambda x: F.lift(F(x))
715
except CoercionException: # the pushout does not exist, we need conversion
716
pinv = lambda x: K(F.lift(~F(x)))
717
proot = lambda x,e: K(F.lift(F(x).nth_root(e, extend = False, all = True)[0]))
718
preduce = lambda x: K(F.lift(F(x)))
719
720
def _pquadroots(a, b, c):
721
r"""
722
Local function returning True iff `ax^2 + bx + c` has roots modulo `P`
723
"""
724
(a, b, c) = (F(a), F(b), F(c))
725
if a == 0:
726
return (b != 0) or (c == 0)
727
elif p == 2:
728
return len(PolynomialRing(F, "x")([c,b,a]).roots()) > 0
729
else:
730
return (b**2 - 4*a*c).is_square()
731
def _pcubicroots(b, c, d):
732
r"""
733
Local function returning the number of roots of `x^3 +
734
b*x^2 + c*x + d` modulo `P`, counting multiplicities
735
"""
736
737
return sum([rr[1] for rr in PolynomialRing(F, 'x')([F(d), F(c), F(b), F(1)]).roots()],0)
738
739
if p == 2:
740
halfmodp = OK(Integer(0))
741
else:
742
halfmodp = pinv(Integer(2))
743
744
A = E.a_invariants()
745
A = [0, A[0], A[1], A[2], A[3], 0, A[4]]
746
indices = [1,2,3,4,6]
747
if min([pval(a) for a in A if a != 0]) < 0:
748
verbose("Non-integral model at P: valuations are %s; making integral"%([pval(a) for a in A if a != 0]), t, 1)
749
e = 0
750
for i in range(7):
751
if A[i] != 0:
752
e = max(e, (-pval(A[i])/i).ceil())
753
pie = pi**e
754
for i in range(7):
755
if A[i] != 0:
756
A[i] *= pie**i
757
verbose("P-integral model is %s, with valuations %s"%([A[i] for i in indices], [pval(A[i]) for i in indices]), t, 1)
758
759
split = None # only relevant for multiplicative reduction
760
761
(a1, a2, a3, a4, a6) = (A[1], A[2], A[3], A[4], A[6])
762
while True:
763
C = EllipticCurve([a1, a2, a3, a4, a6]);
764
(b2, b4, b6, b8) = C.b_invariants()
765
(c4, c6) = C.c_invariants()
766
delta = C.discriminant()
767
val_disc = pval(delta)
768
769
if val_disc == 0:
770
## Good reduction already
771
cp = 1
772
fp = 0
773
KS = KodairaSymbol("I0")
774
break #return
775
776
# Otherwise, we change coordinates so that p | a3, a4, a6
777
if p == 2:
778
if pdiv(b2):
779
r = proot(a4, 2)
780
t = proot(((r + a2)*r + a4)*r + a6, 2)
781
else:
782
temp = pinv(a1)
783
r = temp * a3
784
t = temp * (a4 + r*r)
785
elif p == 3:
786
if pdiv(b2):
787
r = proot(-b6, 3)
788
else:
789
r = -pinv(b2) * b4
790
t = a1 * r + a3
791
else:
792
if pdiv(c4):
793
r = -pinv(12) * b2
794
else:
795
r = -pinv(12*c4) * (c6 + b2 * c4)
796
t = -halfmodp * (a1 * r + a3)
797
r = preduce(r)
798
t = preduce(t)
799
verbose("Before first transform C = %s"%C)
800
verbose("[a1,a2,a3,a4,a6] = %s"%([a1, a2, a3, a4, a6]))
801
C = C.rst_transform(r, 0, t)
802
(a1, a2, a3, a4, a6) = C.a_invariants()
803
(b2, b4, b6, b8) = C.b_invariants()
804
if min([pval(a) for a in (a1, a2, a3, a4, a6) if a != 0]) < 0:
805
raise RuntimeError, "Non-integral model after first transform!"
806
verbose("After first transform %s\n, [a1,a2,a3,a4,a6] = %s\n, valuations = %s"%([r, 0, t], [a1, a2, a3, a4, a6], [pval(a1), pval(a2), pval(a3), pval(a4), pval(a6)]), t, 2)
807
if pval(a3) == 0:
808
raise RuntimeError, "p does not divide a3 after first transform!"
809
if pval(a4) == 0:
810
raise RuntimeError, "p does not divide a4 after first transform!"
811
if pval(a6) == 0:
812
raise RuntimeError, "p does not divide a6 after first transform!"
813
814
# Now we test for Types In, II, III, IV
815
# NB the c invariants never change.
816
817
if not pdiv(c4):
818
# Multiplicative reduction: Type In (n = val_disc)
819
split = False
820
if _pquadroots(1, a1, -a2):
821
cp = val_disc
822
split = True
823
elif Integer(2).divides(val_disc):
824
cp = 2
825
else:
826
cp = 1
827
KS = KodairaSymbol("I%s"%val_disc)
828
fp = 1
829
break #return
830
831
# Additive reduction
832
833
if pval(a6) < 2:
834
## Type II
835
KS = KodairaSymbol("II")
836
fp = val_disc
837
cp = 1
838
break #return
839
if pval(b8) < 3:
840
## Type III
841
KS = KodairaSymbol("III")
842
fp = val_disc - 1
843
cp = 2
844
break #return
845
if pval(b6) < 3:
846
## Type IV
847
cp = 1
848
a3t = preduce(a3/pi)
849
a6t = preduce(a6/pi2)
850
if _pquadroots(1, a3t, -a6t): cp = 3
851
KS = KodairaSymbol("IV")
852
fp = val_disc - 2
853
break #return
854
855
# If our curve is none of these types, we change coords so that
856
# p | a1, a2; p^2 | a3, a4; p^3 | a6
857
if p == 2:
858
s = proot(a2, 2) # so s^2=a2 (mod pi)
859
t = pi*proot(a6/pi2, 2) # so t^2=a6 (mod pi^3)
860
elif p == 3:
861
s = a1 # so a1'=2s+a1=3a1=0 (mod pi)
862
t = a3 # so a3'=2t+a3=3a3=0 (mod pi^2)
863
else:
864
s = -a1*halfmodp # so a1'=2s+a1=0 (mod pi)
865
t = -a3*halfmodp # so a3'=2t+a3=0 (mod pi^2)
866
C = C.rst_transform(0, s, t)
867
(a1, a2, a3, a4, a6) = C.a_invariants()
868
(b2, b4, b6, b8) = C.b_invariants()
869
verbose("After second transform %s\n[a1, a2, a3, a4, a6] = %s\nValuations: %s"%([0, s, t], [a1,a2,a3,a4,a6],[pval(a1),pval(a2),pval(a3),pval(a4),pval(a6)]), t, 2)
870
if pval(a1) == 0:
871
raise RuntimeError, "p does not divide a1 after second transform!"
872
if pval(a2) == 0:
873
raise RuntimeError, "p does not divide a2 after second transform!"
874
if pval(a3) < 2:
875
raise RuntimeError, "p^2 does not divide a3 after second transform!"
876
if pval(a4) < 2:
877
raise RuntimeError, "p^2 does not divide a4 after second transform!"
878
if pval(a6) < 3:
879
raise RuntimeError, "p^3 does not divide a6 after second transform!"
880
if min(pval(a1), pval(a2), pval(a3), pval(a4), pval(a6)) < 0:
881
raise RuntimeError, "Non-integral model after second transform!"
882
883
# Analyze roots of the cubic T^3 + bT^2 + cT + d = 0 mod P, where
884
# b = a2/p, c = a4/p^2, d = a6/p^3
885
b = preduce(a2/pi)
886
c = preduce(a4/pi2)
887
d = preduce(a6/pi3)
888
bb = b*b
889
cc = c*c
890
bc = b*c
891
w = 27*d*d - bb*cc + 4*b*bb*d - 18*bc*d + 4*c*cc
892
x = 3*c - bb
893
if pdiv(w):
894
if pdiv(x):
895
sw = 3
896
else:
897
sw = 2
898
else:
899
sw = 1
900
verbose("Analyzing roots of cubic T^3 + %s*T^2 + %s*T + %s, case %s"%(b, c, d, sw), t, 1)
901
if sw == 1:
902
## Three distinct roots - Type I*0
903
verbose("Distinct roots", t, 1)
904
KS = KodairaSymbol("I0*")
905
cp = 1 + _pcubicroots(b, c, d)
906
fp = val_disc - 4
907
break #return
908
elif sw == 2:
909
## One double root - Type I*m for some m
910
verbose("One double root", t, 1)
911
## Change coords so that the double root is T = 0 mod p
912
if p == 2:
913
r = proot(c, 2)
914
elif p == 3:
915
r = c * pinv(b)
916
else:
917
r = (bc - 9*d)*pinv(2*x)
918
r = pi * preduce(r)
919
C = C.rst_transform(r, 0, 0)
920
(a1, a2, a3, a4, a6) = C.a_invariants()
921
(b2, b4, b6, b8) = C.b_invariants()
922
# The rest of this branch is just to compute cp, fp, KS.
923
# We use pi to keep transforms integral.
924
ix = 3; iy = 3; mx = pi2; my = mx
925
while True:
926
a2t = preduce(a2 / pi)
927
a3t = preduce(a3 / my)
928
a4t = preduce(a4 / (pi*mx))
929
a6t = preduce(a6 / (mx*my))
930
if pdiv(a3t*a3t + 4*a6t):
931
if p == 2:
932
t = my*proot(a6t, 2)
933
else:
934
t = my*preduce(-a3t*halfmodp)
935
C = C.rst_transform(0, 0, t)
936
(a1, a2, a3, a4, a6) = C.a_invariants()
937
(b2, b4, b6, b8) = C.b_invariants()
938
my *= pi
939
iy += 1
940
a2t = preduce(a2 / pi)
941
a3t = preduce(a3/my)
942
a4t = preduce(a4/(pi*mx))
943
a6t = preduce(a6/(mx*my))
944
if pdiv(a4t*a4t - 4*a6t*a2t):
945
if p == 2:
946
r = mx*proot(a6t*pinv(a2t), 2)
947
else:
948
r = mx*preduce(-a4t*pinv(2*a2t))
949
C = C.rst_transform(r, 0, 0)
950
(a1, a2, a3, a4, a6) = C.a_invariants()
951
(b2, b4, b6, b8) = C.b_invariants()
952
mx *= pi
953
ix += 1 # and stay in loop
954
else:
955
if _pquadroots(a2t, a4t, a6t):
956
cp = 4
957
else:
958
cp = 2
959
break # exit loop
960
else:
961
if _pquadroots(1, a3t, -a6t):
962
cp = 4
963
else:
964
cp = 2
965
break
966
KS = KodairaSymbol("I%s*"%(ix+iy-5))
967
fp = val_disc - ix - iy + 1
968
break #return
969
else: # sw == 3
970
## The cubic has a triple root
971
verbose("Triple root", t, 1)
972
## First we change coordinates so that T = 0 mod p
973
if p == 2:
974
r = b
975
elif p == 3:
976
r = proot(-d, 3)
977
else:
978
r = -b * pinv(3)
979
r = pi*preduce(r)
980
C = C.rst_transform(r, 0, 0)
981
(a1, a2, a3, a4, a6) = C.a_invariants()
982
(b2, b4, b6, b8) = C.b_invariants()
983
verbose("After third transform %s\n[a1,a2,a3,a4,a6] = %s\nValuations: %s"%([r,0,0],[a1,a2,a3,a4,a6],[pval(ai) for ai in [a1,a2,a3,a4,a6]]), t, 2)
984
if min(pval(ai) for ai in [a1,a2,a3,a4,a6]) < 0:
985
raise RuntimeError, "Non-integral model after third transform!"
986
if pval(a2) < 2 or pval(a4) < 3 or pval(a6) < 4:
987
raise RuntimeError, "Cubic after transform does not have a triple root at 0"
988
a3t = preduce(a3/pi2)
989
a6t = preduce(a6/pi4)
990
# We test for Type IV*
991
if not pdiv(a3t*a3t + 4*a6t):
992
cp = 3 if _pquadroots(1, a3t, -a6t) else 1
993
KS = KodairaSymbol("IV*")
994
fp = val_disc - 6
995
break #return
996
# Now change coordinates so that p^3|a3, p^5|a6
997
if p==2:
998
t = -pi2*proot(a6t, 2)
999
else:
1000
t = pi2*preduce(-a3t*halfmodp)
1001
C = C.rst_transform(0, 0, t)
1002
(a1, a2, a3, a4, a6) = C.a_invariants()
1003
(b2, b4, b6, b8) = C.b_invariants()
1004
# We test for types III* and II*
1005
if pval(a4) < 4:
1006
## Type III*
1007
KS = KodairaSymbol("III*")
1008
fp = val_disc - 7
1009
cp = 2
1010
break #return
1011
if pval(a6) < 6:
1012
## Type II*
1013
KS = KodairaSymbol("II*")
1014
fp = val_disc - 8
1015
cp = 1
1016
break #return
1017
if pi_neg is None:
1018
if principal_flag:
1019
pi_neg = pi
1020
else:
1021
pi_neg = K.uniformizer(P, 'negative')
1022
pi_neg2 = pi_neg*pi_neg
1023
pi_neg3 = pi_neg*pi_neg2
1024
pi_neg4 = pi_neg*pi_neg3
1025
pi_neg6 = pi_neg4*pi_neg2
1026
a1 /= pi_neg
1027
a2 /= pi_neg2
1028
a3 /= pi_neg3
1029
a4 /= pi_neg4
1030
a6 /= pi_neg6
1031
verbose("Non-minimal equation, dividing out...\nNew model is %s"%([a1, a2, a3, a4, a6]), t, 1)
1032
return (C, p, val_disc, fp, KS, cp, split)
1033
1034
1035
def check_prime(K,P):
1036
r"""
1037
Function to check that `P` determines a prime of `K`, and return that ideal.
1038
1039
INPUT:
1040
1041
- ``K`` -- a number field (including `\QQ`).
1042
1043
- ``P`` -- an element of ``K`` or a (fractional) ideal of ``K``.
1044
1045
OUTPUT:
1046
1047
- If ``K`` is `\QQ`: the prime integer equal to or which generates `P`.
1048
1049
- If ``K`` is not `\QQ`: the prime ideal equal to or generated by `P`.
1050
1051
.. note::
1052
1053
If `P` is not a prime and does not generate a prime, a TypeError is raised.
1054
1055
EXAMPLES::
1056
1057
sage: from sage.schemes.elliptic_curves.ell_local_data import check_prime
1058
sage: check_prime(QQ,3)
1059
3
1060
sage: check_prime(QQ,ZZ.ideal(31))
1061
31
1062
sage: K.<a>=NumberField(x^2-5)
1063
sage: check_prime(K,a)
1064
Fractional ideal (a)
1065
sage: check_prime(K,a+1)
1066
Fractional ideal (a + 1)
1067
sage: [check_prime(K,P) for P in K.primes_above(31)]
1068
[Fractional ideal (5/2*a + 1/2), Fractional ideal (5/2*a - 1/2)]
1069
"""
1070
if K is QQ:
1071
if isinstance(P, (int,long,Integer)):
1072
P = Integer(P)
1073
if P.is_prime():
1074
return P
1075
else:
1076
raise TypeError, "%s is not prime"%P
1077
else:
1078
if is_Ideal(P) and P.base_ring() is ZZ and P.is_prime():
1079
return P.gen()
1080
raise TypeError, "%s is not a prime ideal of %s"%(P,ZZ)
1081
1082
if not is_NumberField(K):
1083
raise TypeError, "%s is not a number field"%K
1084
1085
if is_NumberFieldFractionalIdeal(P):
1086
if P.is_prime():
1087
return P
1088
else:
1089
raise TypeError, "%s is not a prime ideal of %s"%(P,K)
1090
1091
if is_NumberFieldElement(P):
1092
if P in K:
1093
P = K.ideal(P)
1094
else:
1095
raise TypeError, "%s is not an element of %s"%(P,K)
1096
if P.is_prime():
1097
return P
1098
else:
1099
raise TypeError, "%s is not a prime ideal of %s"%(P,K)
1100
1101
raise TypeError, "%s is not a valid prime of %s"%(P,K)
1102
1103