Book a Demo!
CoCalc Logo Icon
StoreFeaturesDocsShareSupportNewsAboutPoliciesSign UpSign In
sagemath
GitHub Repository: sagemath/sagelib
Path: blob/master/sage/schemes/elliptic_curves/ell_finite_field.py
4156 views
1
"""
2
Elliptic curves over finite fields
3
4
AUTHORS:
5
6
- William Stein (2005): Initial version
7
8
- Robert Bradshaw et al....
9
10
- John Cremona (2008-02): Point counting and group structure for
11
non-prime fields, Frobenius endomorphism and order, elliptic logs
12
13
- Mariah Lenox (2011-03): Added set_order method
14
"""
15
16
#*****************************************************************************
17
# Copyright (C) 2005 William Stein <[email protected]>
18
#
19
# Distributed under the terms of the GNU General Public License (GPL)
20
#
21
# This code is distributed in the hope that it will be useful,
22
# but WITHOUT ANY WARRANTY; without even the implied warranty of
23
# MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the GNU
24
# General Public License for more details.
25
#
26
# The full text of the GPL is available at:
27
#
28
# http://www.gnu.org/licenses/
29
#*****************************************************************************
30
31
from sage.misc.randstate import current_randstate
32
33
from sage.schemes.plane_curves.projective_curve import Hasse_bounds
34
from ell_field import EllipticCurve_field
35
from constructor import EllipticCurve, EllipticCurve_from_j
36
from sage.schemes.hyperelliptic_curves.hyperelliptic_finite_field import HyperellipticCurve_finite_field
37
import sage.rings.ring as ring
38
from sage.rings.all import Integer, ZZ, PolynomialRing, GF, polygen
39
from sage.rings.finite_rings.all import is_FiniteFieldElement
40
import sage.groups.generic as generic
41
import ell_point
42
from sage.rings.arith import gcd, lcm
43
from sage.structure.sequence import Sequence
44
45
import sage.plot.all as plot
46
47
import sage.libs.pari
48
pari = sage.libs.pari.all.pari
49
50
class EllipticCurve_finite_field(EllipticCurve_field, HyperellipticCurve_finite_field):
51
"""
52
Elliptic curve over a finite field.
53
"""
54
def __init__(self, x, y=None):
55
"""
56
Special constructor for elliptic curves over a finite field
57
58
EXAMPLES::
59
60
sage: EllipticCurve(GF(101),[2,3])
61
Elliptic Curve defined by y^2 = x^3 + 2*x + 3 over Finite Field of size 101
62
63
::
64
65
sage: F=GF(101^2, 'a')
66
sage: EllipticCurve([F(2),F(3)])
67
Elliptic Curve defined by y^2 = x^3 + 2*x + 3 over Finite Field in a of size 101^2
68
69
Elliptic curves over `\ZZ/N\ZZ` with `N` prime are of type
70
"elliptic curve over a finite field":
71
72
sage: F = Zmod(101)
73
sage: EllipticCurve(F, [2, 3])
74
Elliptic Curve defined by y^2 = x^3 + 2*x + 3 over Ring of integers modulo 101
75
sage: E = EllipticCurve([F(2), F(3)])
76
sage: type(E)
77
<class 'sage.schemes.elliptic_curves.ell_finite_field.EllipticCurve_finite_field_with_category'>
78
sage: E.category()
79
Category of schemes over Ring of integers modulo 101
80
81
Elliptic curves over `\ZZ/N\ZZ` with `N` composite are of type
82
"generic elliptic curve"::
83
84
sage: F = Zmod(95)
85
sage: EllipticCurve(F, [2, 3])
86
Elliptic Curve defined by y^2 = x^3 + 2*x + 3 over Ring of integers modulo 95
87
sage: E = EllipticCurve([F(2), F(3)])
88
sage: type(E)
89
<class 'sage.schemes.elliptic_curves.ell_generic.EllipticCurve_generic_with_category'>
90
sage: E.category()
91
Category of schemes over Ring of integers modulo 95
92
sage: TestSuite(E).run(skip=["_test_elements"])
93
"""
94
if isinstance(x, list):
95
seq = Sequence(x)
96
else:
97
seq = Sequence(y, universe=x)
98
ainvs = list(seq)
99
field = seq.universe()
100
if not isinstance(field, ring.Ring):
101
raise TypeError
102
103
EllipticCurve_field.__init__(self, ainvs)
104
105
self._point = ell_point.EllipticCurvePoint_finite_field
106
107
def plot(self, *args, **kwds):
108
"""
109
Draw a graph of this elliptic curve over a prime finite field.
110
111
INPUT:
112
113
114
- ``*args, **kwds`` - all other options are passed
115
to the circle graphing primitive.
116
117
118
EXAMPLES::
119
120
sage: E = EllipticCurve(FiniteField(17), [0,1])
121
sage: P = plot(E, rgbcolor=(0,0,1))
122
"""
123
R = self.base_ring()
124
if not R.is_prime_field():
125
raise NotImplementedError
126
127
G = plot.Graphics()
128
G += plot.points([P[0:2] for P in self.points() if not P.is_zero()], *args, **kwds)
129
130
return G
131
132
def _points_via_group_structure(self):
133
"""
134
Return a list of all the points on the curve, for prime fields only
135
(see points() for the general case)
136
137
EXAMPLES::
138
139
sage: S=EllipticCurve(GF(97),[2,3])._points_via_group_structure()
140
sage: len(S)
141
100
142
143
See trac #4687, where the following example did not work::
144
145
sage: E=EllipticCurve(GF(2),[0, 0, 1, 1, 1])
146
sage: E.points()
147
[(0 : 1 : 0)]
148
149
::
150
151
sage: E=EllipticCurve(GF(2),[0, 0, 1, 0, 1])
152
sage: E.points()
153
[(0 : 1 : 0), (1 : 0 : 1), (1 : 1 : 1)]
154
155
::
156
157
sage: E=EllipticCurve(GF(4,'a'),[0, 0, 1, 0, 1])
158
sage: E.points()
159
[(0 : 1 : 0), (0 : a : 1), (0 : a + 1 : 1), (1 : 0 : 1), (1 : 1 : 1), (a : 0 : 1), (a : 1 : 1), (a + 1 : 0 : 1), (a + 1 : 1 : 1)]
160
"""
161
# TODO, eliminate when polynomial calling is fast
162
G = self.abelian_group()
163
pts = [x.element() for x in G.gens()]
164
165
ni = G.generator_orders()
166
ngens = G.ngens()
167
168
H0=[self(0)]
169
if ngens == 0: # trivial group
170
return H0
171
for m in range(1,ni[0]): H0.append(H0[-1]+pts[0])
172
if ngens == 1: # cyclic group
173
return H0
174
175
# else noncyclic group
176
H1=[self(0)]
177
for m in range(1,ni[1]): H1.append(H1[-1]+pts[1])
178
return [P+Q for P in H0 for Q in H1]
179
180
def points(self):
181
r"""
182
All the points on this elliptic curve. The list of points is cached
183
so subsequent calls are free.
184
185
EXAMPLES::
186
187
sage: p = 5
188
sage: F = GF(p)
189
sage: E = EllipticCurve(F, [1, 3])
190
sage: a_sub_p = E.change_ring(QQ).ap(p); a_sub_p
191
2
192
193
::
194
195
sage: len(E.points())
196
4
197
sage: p + 1 - a_sub_p
198
4
199
sage: E.points()
200
[(0 : 1 : 0), (1 : 0 : 1), (4 : 1 : 1), (4 : 4 : 1)]
201
202
::
203
204
sage: K = GF(p**2,'a')
205
sage: E = E.change_ring(K)
206
sage: len(E.points())
207
32
208
sage: (p + 1)**2 - a_sub_p**2
209
32
210
sage: w = E.points(); w
211
[(0 : 1 : 0), (0 : 2*a + 4 : 1), (0 : 3*a + 1 : 1), (1 : 0 : 1), (2 : 2*a + 4 : 1), (2 : 3*a + 1 : 1), (3 : 2*a + 4 : 1), (3 : 3*a + 1 : 1), (4 : 1 : 1), (4 : 4 : 1), (a : 1 : 1), (a : 4 : 1), (a + 2 : a + 1 : 1), (a + 2 : 4*a + 4 : 1), (a + 3 : a : 1), (a + 3 : 4*a : 1), (a + 4 : 0 : 1), (2*a : 2*a : 1), (2*a : 3*a : 1), (2*a + 4 : a + 1 : 1), (2*a + 4 : 4*a + 4 : 1), (3*a + 1 : a + 3 : 1), (3*a + 1 : 4*a + 2 : 1), (3*a + 2 : 2*a + 3 : 1), (3*a + 2 : 3*a + 2 : 1), (4*a : 0 : 1), (4*a + 1 : 1 : 1), (4*a + 1 : 4 : 1), (4*a + 3 : a + 3 : 1), (4*a + 3 : 4*a + 2 : 1), (4*a + 4 : a + 4 : 1), (4*a + 4 : 4*a + 1 : 1)]
212
213
Note that the returned list is an immutable sorted Sequence::
214
215
sage: w[0] = 9
216
Traceback (most recent call last):
217
...
218
ValueError: object is immutable; please change a copy instead.
219
"""
220
try:
221
return self.__points
222
except AttributeError: pass
223
224
from sage.structure.sequence import Sequence
225
k = self.base_ring()
226
if k.is_prime_field() and k.order()>50:
227
v = self._points_via_group_structure()
228
else:
229
v =self._points_fast_sqrt()
230
v.sort()
231
self.__points = Sequence(v, immutable=True)
232
return self.__points
233
234
rational_points = points
235
236
def count_points(self, n=1):
237
"""
238
Returns the cardinality of this elliptic curve over the base field or extensions.
239
240
INPUT:
241
242
- ``n`` (int) -- a positive integer
243
244
OUTPUT:
245
246
If `n=1`, returns the cardinality of the curve over its base field.
247
248
If `n>1`, returns a list `[c_1, c_2, ..., c_n]` where `c_d` is
249
the cardinality of the curve over the extension of degree `d`
250
of its base field.
251
252
EXAMPLES::
253
254
sage: p = 101
255
sage: F = GF(p)
256
sage: E = EllipticCurve(F, [2,3])
257
sage: E.count_points(1)
258
96
259
sage: E.count_points(5)
260
[96, 10368, 1031904, 104053248, 10509895776]
261
262
::
263
264
sage: F.<a> = GF(p^2)
265
sage: E = EllipticCurve(F, [a,a])
266
sage: E.cardinality()
267
10295
268
sage: E.count_points()
269
10295
270
sage: E.count_points(1)
271
10295
272
sage: E.count_points(5)
273
[10295, 104072155, 1061518108880, 10828567126268595, 110462212555439192375]
274
275
"""
276
try:
277
n = Integer(n)
278
except TypeError:
279
raise TypeError, "n must be a positive integer"
280
281
if n<1:
282
raise ValueError, "n must be a positive integer"
283
284
if n==1:
285
return self.cardinality()
286
else:
287
return [ self.cardinality(extension_degree=i) for i in range(1,n+1)]
288
289
def random_element(self):
290
"""
291
Returns a random point on this elliptic curve.
292
293
If `q` is small, finds all points and returns one at random.
294
Otherwise, returns the point at infinity with probability
295
`1/(q+1)` where the base field has cardinality `q`, and then
296
picks random `x`-coordinates from the base field until one
297
gives a rational point.
298
299
EXAMPLES::
300
301
sage: k = GF(next_prime(7^5))
302
sage: E = EllipticCurve(k,[2,4])
303
sage: P = E.random_element(); P
304
(16740 : 12486 : 1)
305
sage: type(P)
306
<class 'sage.schemes.elliptic_curves.ell_point.EllipticCurvePoint_finite_field'>
307
sage: P in E
308
True
309
310
::
311
312
sage: k.<a> = GF(7^5)
313
sage: E = EllipticCurve(k,[2,4])
314
sage: P = E.random_element(); P
315
(2*a^4 + 3*a^3 + 5*a^2 + 6*a + 4 : 6*a^4 + 4*a^3 + a + 6 : 1)
316
sage: type(P)
317
<class 'sage.schemes.elliptic_curves.ell_point.EllipticCurvePoint_finite_field'>
318
sage: P in E
319
True
320
321
::
322
323
sage: k.<a> = GF(2^5)
324
sage: E = EllipticCurve(k,[a^2,a,1,a+1,1])
325
sage: P = E.random_element(); P
326
(a^4 + a^2 + 1 : a^3 + a : 1)
327
sage: type(P)
328
<class 'sage.schemes.elliptic_curves.ell_point.EllipticCurvePoint_finite_field'>
329
sage: P in E
330
True
331
332
Ensure that the entire point set is reachable::
333
334
sage: E = EllipticCurve(GF(11), [2,1])
335
sage: len(set(E.random_element() for _ in range(100)))
336
16
337
sage: E.cardinality()
338
16
339
340
TESTS:
341
342
See trac #8311::
343
344
sage: E = EllipticCurve(GF(3), [0,0,0,2,2])
345
sage: E.random_element()
346
(0 : 1 : 0)
347
sage: E.cardinality()
348
1
349
350
sage: E = EllipticCurve(GF(2), [0,0,1,1,1])
351
sage: E.random_point()
352
(0 : 1 : 0)
353
sage: E.cardinality()
354
1
355
356
sage: F.<a> = GF(4)
357
sage: E = EllipticCurve(F, [0, 0, 1, 0, a])
358
sage: E.random_point()
359
(0 : 1 : 0)
360
sage: E.cardinality()
361
1
362
363
"""
364
random = current_randstate().c_rand_double
365
k = self.base_field()
366
q = k.order()
367
368
# For small fields we find all the rational points and pick
369
# one at random. Note that the group can be trivial for
370
# q=2,3,4 only (see #8311) so these cases need special
371
# treatment.
372
373
if q < 5:
374
pts = self.points() # will be cached
375
return pts[ZZ.random_element(len(pts))]
376
377
378
# The following allows the origin self(0) to be picked
379
if random() <= 1/float(q+1):
380
return self(0)
381
382
while True:
383
v = self.lift_x(k.random_element(), all=True)
384
if v:
385
return v[int(random() * len(v))]
386
387
random_point = random_element
388
389
390
def trace_of_frobenius(self):
391
r"""
392
Return the trace of Frobenius acting on this elliptic curve.
393
394
.. note::
395
396
This computes the curve cardinality, which may be
397
time-consuming.
398
399
EXAMPLES::
400
401
sage: E=EllipticCurve(GF(101),[2,3])
402
sage: E.trace_of_frobenius()
403
6
404
sage: E=EllipticCurve(GF(11^5,'a'),[2,5])
405
sage: E.trace_of_frobenius()
406
802
407
408
The following shows that the issue from trac #2849 is fixed::
409
410
sage: E=EllipticCurve(GF(3^5,'a'),[-1,-1])
411
sage: E.trace_of_frobenius()
412
-27
413
"""
414
return 1 + self.base_field().order() - self.cardinality()
415
416
def _cardinality_with_j_invariant_1728(self):
417
r"""
418
Special function to compute cardinality when j=1728.
419
420
EXAMPLES: An example with q=p=1 (mod 4)
421
422
::
423
424
sage: F=GF(10009)
425
sage: [EllipticCurve(F,[0,0,0,11^i,0])._cardinality_with_j_invariant_1728() for i in range(4)]
426
[10016, 10210, 10004, 9810]
427
428
An example with q=p=3 (mod 4)
429
430
::
431
432
sage: F=GF(10007)
433
sage: [EllipticCurve(F,[0,0,0,5^i,0])._cardinality_with_j_invariant_1728() for i in range(4)]
434
[10008, 10008, 10008, 10008]
435
436
An example with `q=p^2`, p=3 (mod 4)
437
438
::
439
440
sage: F.<a>=GF(10007^2,'a')
441
sage: [EllipticCurve(F,[0,0,0,a^i,0])._cardinality_with_j_invariant_1728() for i in range(4)]
442
[100160064, 100140050, 100120036, 100140050]
443
444
Examples with `q=2^d`, d odd (3 isomorphism classes)::
445
446
sage: F.<a> = GF(2**15,'a')
447
sage: ais = [[0,0,1,0,0],[0,0,1,1,0],[0,0,1,1,1]]
448
sage: curves=[EllipticCurve(F,ai) for ai in ais]
449
sage: all([all([e1==e2 or not e1.is_isomorphic(e2) for e1 in curves]) for e2 in curves])
450
True
451
sage: [e._cardinality_with_j_invariant_1728() for e in curves]
452
[32769, 33025, 32513]
453
454
Examples with `q=2^d`, d even (7 isomorphism classes)::
455
456
sage: F.<a> = GF(2**16,'a')
457
sage: b = a^11 # trace 1
458
sage: ais = [[0,0,1,0,0],[0,0,1,0,b],[0,0,1,b,0],[0,0,a,0,0],[0,0,a,0,a^2*b],[0,0,a^2,0,0],[0,0,a^2,0,a^4*b]]
459
sage: curves=[EllipticCurve(F,ai) for ai in ais]
460
sage: all([all([e1==e2 or not e1.is_isomorphic(e2) for e1 in curves]) for e2 in curves])
461
True
462
sage: [e._cardinality_with_j_invariant_1728() for e in curves]
463
[65025, 66049, 65537, 65793, 65281, 65793, 65281]
464
465
Examples with `q=3^d`, d odd (4 isomorphism classes)::
466
467
sage: F.<a> = GF(3**15,'a')
468
sage: b=a^7 # has trace 1
469
sage: ais=[[0,0,0,1,0],[0,0,0,-1,0],[0,0,0,-1,b],[0,0,0,-1,-b]]
470
sage: curves=[EllipticCurve(F,ai) for ai in ais]
471
sage: all([all([e1==e2 or not e1.is_isomorphic(e2) for e1 in curves]) for e2 in curves])
472
True
473
sage: [e._cardinality_with_j_invariant_1728() for e in curves]
474
[14348908, 14348908, 14342347, 14355469]
475
476
Examples with `q=3^d`, d even (6 isomorphism classes)::
477
478
sage: F.<g>=GF(3^18,'g')
479
sage: i=F(-1).sqrt()
480
sage: a=g^8 # has trace 1
481
sage: ais= [[0,0,0,1,0],[0,0,0,1,i*a],[0,0,0,g,0],[0,0,0,g^3,0],[0,0,0,g^2,0], [0,0,0,g^2,i*a*g^3]]
482
sage: curves=[EllipticCurve(F,ai) for ai in ais]
483
sage: all([all([e1==e2 or not e1.is_isomorphic(e2) for e1 in curves]) for e2 in curves])
484
True
485
sage: [E._cardinality_with_j_invariant_1728() for E in curves]
486
[387459856, 387400807, 387420490, 387420490, 387381124, 387440173]
487
"""
488
try:
489
return self._order
490
except AttributeError:
491
pass
492
493
k = self.base_ring()
494
assert self.j_invariant()==k(1728)
495
q = k.cardinality()
496
p = k.characteristic()
497
d = k.degree()
498
x=polygen(ZZ)
499
500
# p=2, j=0=1728
501
#
502
# Number of isomorphism classes is 3 (in odd degree) or 7 (in even degree)
503
#
504
if p==2:
505
if d%2==1:
506
# The 3 classes are represented, independently of d,
507
# by [0,0,1,0,0], [0,0,1,1,0], [0,0,1,1,1]
508
E=EllipticCurve(k,[0,0,1,0,0])
509
if self.is_isomorphic(E):
510
N = q+1
511
else:
512
n = (d+1)//2
513
t = 2**n
514
n = n%4
515
if n==0 or n==1: t=-t
516
E=EllipticCurve(k,[0,0,1,1,1])
517
if self.is_isomorphic(E): t=-t
518
N = q+1-t
519
else:
520
# The 7 classes are represented by E1=[0,0,1,0,0],
521
# E2=[0,0,1,0,b], E3=[0,0,1,b,0], E4=[0,0,a,0,0],
522
# E4=[0,0,a,0,a^2*b], E6=[0,0,a^2,0,0],
523
# E7=[0,0,a^2,0,a^4*b], where a is a non-cube and b
524
# has trace 1. E1's Frobenius is pi=(-2)**(d//2); the
525
# Frobeniuses are then pi, -pi, 0; w*pi, -w*pi;
526
# w^2*pi, -w^2*pi where w is either cube root of
527
# unity, so the traces are 2*pi, -2*pi, 0, -pi, +pi;
528
# -pi, +pi.
529
delta = self.discriminant()
530
discube = (delta**((q-1)//3) == k(1))
531
pi = (-2)**(d//2)
532
if discube:
533
a = k.gen()
534
b = a
535
while b.trace()==0: b*=a
536
if self.is_isomorphic(EllipticCurve(k,[0,0,1,b,0])):
537
t = 0
538
else:
539
t = 2*pi
540
if not self.is_isomorphic(EllipticCurve(k,[0,0,1,0,0])):
541
t = -t
542
543
else:
544
t = pi
545
if self.is_isomorphic(EllipticCurve(k,[0,0,delta,0,0])):
546
t = -t
547
N = q+1-t
548
549
550
# p=3, j=0=1728
551
#
552
# Number of isomorphism classes is 4 (odd degree) or 6 (even degree)
553
#
554
elif p==3:
555
if d%2==1:
556
# The 4 classes are represented by [0,0,0,1,0],
557
# [0,0,0,-1,0], [0,0,0,-1,a], [0,0,0,-1,-a] where a
558
# has trace 1
559
delta = self.discriminant()
560
if (-delta).is_square():
561
t = 0
562
else:
563
u = delta.sqrt()
564
if not u.is_square(): u=-u
565
tr = ((self.a3()**2+self.a6())/u).trace()
566
if tr==0:
567
t = 0
568
else:
569
d2 = (d+1)//2
570
t = 3**d2
571
if d2%2==1: t = -t
572
if tr==-1: t = -t
573
N = q+1-t
574
else:
575
# The 6 classes are represented by: [0,0,0,1,0],
576
# [0,0,0,1,i*a]; [0,0,0,g,0], [0,0,0,g^3,0];
577
# [0,0,0,g^2,0], [0,0,0,g^2,i*a*g^3]; where g
578
# generates the multiplicative group modulo 4th
579
# powers, and a has nonzero trace.
580
A4 = self.a4()-self.a1()*self.a3()
581
i = k(-1).sqrt()
582
t = 0
583
if A4.is_square():
584
u = A4.sqrt()
585
t = (-3)**(d//2)
586
if u.is_square():
587
A6 = (self.a3()**2+self.a6())/u
588
if (i*A6).trace()==0:
589
t = 2*t
590
else:
591
t = -t
592
else:
593
A6 = (self.a3()**2+self.a6())/(u*A4)
594
if (i*A6).trace()==0:
595
t = -2*t
596
N = q+1-t
597
598
# p>3, j=1728
599
#
600
# Number of isomorphism classes is 4 if q=1 (mod 4), else 2
601
#
602
elif p%4==3:
603
if d%2==1:
604
t = 0
605
else:
606
t = (-p)**(d//2)
607
w = (self.c4()/k(48))**((q-1)//4)
608
if w==1: t = 2*t
609
elif w==-1: t = -2*t
610
else: t = 0
611
612
N = q+1-t
613
614
# p=1 (mod 4). First find Frobenius pi=a+b*i for [0,0,0,-1,0] over GF(p):
615
# N(pi)=p and N(pi-1)=0 (mod 8).
616
#
617
else: # p%4==1
618
R = ZZ.extension(x**2+1,'i')
619
i = R.gen(1)
620
pi = R.fraction_field().factor(p)[0][0].gens_reduced()[0]
621
a,b = pi.list()
622
if a%2==0:
623
a,b = -b,a
624
if (a+b+1)%4==0:
625
a,b = -a,-b
626
pi = a+b*i # Now pi=a+b*i with (a,b)=(1,0),(3,2) mod 4
627
628
# Lift to Frobenius for [0,0,0,-1,0] over GF(p^d):
629
630
if d>1:
631
pi = pi**d
632
a,b = pi.list()
633
634
# Compute appropriate quartic twist:
635
636
w = (self.c4()/k(48))**((q-1)//4)
637
if w==1:
638
t = 2*a
639
elif w==-1:
640
t = -2*a
641
elif k(b)==w*k(a):
642
t = 2*b
643
else:
644
t = -2*b
645
N = q+1-t
646
647
self._order = Integer(N)
648
return self._order
649
650
def _cardinality_with_j_invariant_0(self):
651
r"""
652
Special function to compute cardinality when j=0.
653
654
EXAMPLES: An example with q=p=1 (mod 6)
655
656
::
657
658
sage: F=GF(1009)
659
sage: [EllipticCurve(F,[0,0,0,0,11^i])._cardinality_with_j_invariant_0() for i in range(6)]
660
[948, 967, 1029, 1072, 1053, 991]
661
662
An example with q=p=5 (mod 6)
663
664
::
665
666
sage: F=GF(1013)
667
sage: [EllipticCurve(F,[0,0,0,0,3^i])._cardinality_with_j_invariant_0() for i in range(6)]
668
[1014, 1014, 1014, 1014, 1014, 1014]
669
670
An example with `q=p^2`, p=5 (mod 6)
671
672
::
673
674
sage: F.<a>=GF(1013^2,'a')
675
sage: [EllipticCurve(F,[0,0,0,0,a^i])._cardinality_with_j_invariant_0() for i in range(6)]
676
[1028196, 1027183, 1025157, 1024144, 1025157, 1027183]
677
678
For examples in characteristic 2 and 3, see the function
679
_cardinality_with_j_invariant_1728()
680
"""
681
682
try:
683
return self._order
684
except AttributeError:
685
pass
686
687
k = self.base_ring()
688
assert self.j_invariant()==k(0)
689
p = k.characteristic()
690
if p==2 or p==3: # then 0==1728
691
return self._cardinality_with_j_invariant_1728()
692
693
q = k.cardinality()
694
d = k.degree()
695
x=polygen(ZZ)
696
697
# p>3, j=0
698
#
699
# Number of isomorphism classes is 4 if q=1 (mod 4), else 2
700
#
701
if p%6==5:
702
if d%2==1:
703
t = 0
704
else:
705
t = (-p)**(d//2)
706
w = (self.c6()/k(-864))**((q-1)//6)
707
if w==1: t = 2*t
708
elif w==-1: t = -2*t
709
elif w**3==1: t = -t
710
711
N = q+1-t
712
713
# p=1 (mod 6). First find Frobenius pi=a+b*w for [0,0,0,0,1] over GF(p):
714
# N(pi)=p and N(pi-1)=0 (mod 12).
715
#
716
else: # p%6==1
717
R = ZZ.extension(x**2-x+1,'zeta6')
718
zeta6 = R.gen(1)
719
pi = R.fraction_field().factor(p)[0][0].gens_reduced()[0]
720
while (pi-1).norm()%12 !=0: pi*=zeta6
721
a,b = pi.list()
722
z = k(-b)/k(a) # a *specific* 6th root of unity in k
723
724
# Now pi=a+b*zeta6 with N(pi-1)=0 (mod 12)
725
726
# Lift to Frobenius for [0,0,0,0,1] over GF(p^d):
727
728
if d>1:
729
pi = pi**d
730
a,b = pi.list()
731
732
# Compute appropriate sextic twist:
733
734
w = (self.c6()/k(-864))**((q-1)//6)
735
736
if w==1: t = 2*a+b # = Trace(pi)
737
elif w==-1: t = -2*a-b # = Trace(-pi)
738
elif w==z: t = a-b # = Trace(pi*zeta6)
739
elif w==z**2: t = -a-2*b # = Trace(pi*zeta6**2)
740
elif w==z**4: t = b-a # = Trace(pi*zeta6**4)
741
elif w==z**5: t = a+2*b # = Trace(pi*zeta6**5)
742
743
N = q+1-t
744
745
self._order = Integer(N)
746
return self._order
747
748
def cardinality(self, algorithm='heuristic', extension_degree=1):
749
r"""
750
Return the number of points on this elliptic curve over an
751
extension field (default: the base field).
752
753
INPUT:
754
755
756
- ``algorithm`` - string (default: 'heuristic'), used
757
only for point counting over prime fields
758
759
- ``'heuristic'`` - use a heuristic to choose between
760
``'pari'`` and ``'bsgs'``.
761
762
- ``'pari'`` - use the baby step giant step or SEA methods
763
as implemented in PARI via the C-library function ellap.
764
765
- ``'bsgs'`` - use the baby step giant step method as
766
implemented in Sage, with the Cremona-Sutherland version
767
of Mestre's trick.
768
769
- ``'all'`` - (over prime fields only) compute cardinality
770
with all of PARI and bsgs; return result if they agree
771
or raise a RuntimeError if they do not.
772
773
- ``extension_degree`` - int (default: 1); if the
774
base field is `k=GF(p^n)` and extension_degree=d, returns
775
the cardinality of `E(GF(p^{n d}))`.
776
777
778
OUTPUT: an integer
779
780
The cardinality is cached.
781
782
Over prime fields, one of the above algorithms is used. Over
783
non-prime fields, the serious point counting is done on a standard
784
curve with the same j-invariant over the field GF(p)(j), then
785
lifted to the base_field, and finally account is taken of twists.
786
787
For j=0 and j=1728 special formulas are used instead.
788
789
EXAMPLES::
790
791
sage: EllipticCurve(GF(4,'a'),[1,2,3,4,5]).cardinality()
792
8
793
sage: k.<a> = GF(3^3)
794
sage: l = [a^2 + 1, 2*a^2 + 2*a + 1, a^2 + a + 1, 2, 2*a]
795
sage: EllipticCurve(k,l).cardinality()
796
29
797
798
::
799
800
sage: l = [1, 1, 0, 2, 0]
801
sage: EllipticCurve(k,l).cardinality()
802
38
803
804
An even bigger extension (which we check against Magma)::
805
806
sage: EllipticCurve(GF(3^100,'a'),[1,2,3,4,5]).cardinality()
807
515377520732011331036459693969645888996929981504
808
sage: magma.eval("Order(EllipticCurve([GF(3^100)|1,2,3,4,5]))") # optional - magma
809
'515377520732011331036459693969645888996929981504'
810
811
::
812
813
sage: EllipticCurve(GF(10007),[1,2,3,4,5]).cardinality()
814
10076
815
sage: EllipticCurve(GF(10007),[1,2,3,4,5]).cardinality(algorithm='pari')
816
10076
817
sage: EllipticCurve(GF(next_prime(10**20)),[1,2,3,4,5]).cardinality()
818
100000000011093199520
819
820
The cardinality is cached::
821
822
sage: E = EllipticCurve(GF(3^100,'a'),[1,2,3,4,5])
823
sage: E.cardinality() is E.cardinality()
824
True
825
sage: E=EllipticCurve(GF(11^2,'a'),[3,3])
826
sage: E.cardinality()
827
128
828
sage: EllipticCurve(GF(11^100,'a'),[3,3]).cardinality()
829
137806123398222701841183371720896367762643312000384671846835266941791510341065565176497846502742959856128
830
831
TESTS::
832
833
sage: EllipticCurve(GF(10007),[1,2,3,4,5]).cardinality(algorithm='foobar')
834
Traceback (most recent call last):
835
...
836
ValueError: Algorithm is not known
837
"""
838
if extension_degree>1:
839
# A recursive call to cardinality() with
840
# extension_degree=1, which will cache the cardinality, is
841
# made by the call to frobenius_order() here:
842
R=self.frobenius_order()
843
if R.degree()==1:
844
return (self.frobenius()**extension_degree-1)**2
845
else:
846
return (self.frobenius()**extension_degree-1).norm()
847
848
# Now extension_degree==1
849
if algorithm != 'all':
850
try:
851
return self._order
852
except AttributeError:
853
pass
854
855
k = self.base_ring()
856
q = k.cardinality()
857
858
if q < 50:
859
return self.cardinality_exhaustive()
860
861
# use special code for j=0, 1728 (for any field)
862
j = self.j_invariant()
863
if j==k(0):
864
return self._cardinality_with_j_invariant_0()
865
if j==k(1728):
866
return self._cardinality_with_j_invariant_1728()
867
868
N = 0
869
p = k.characteristic()
870
d = k.degree()
871
872
# Over prime fields, we have a variety of algorithms to choose from:
873
874
if d == 1:
875
if algorithm == 'heuristic':
876
algorithm = 'pari'
877
if algorithm == 'pari':
878
N = self.cardinality_pari()
879
elif algorithm == 'sea':
880
N = self.cardinality_pari() # purely for backwards compatibility
881
elif algorithm == 'bsgs':
882
N = self.cardinality_bsgs()
883
elif algorithm == 'all':
884
N1 = self.cardinality_pari()
885
N2 = self.cardinality_bsgs()
886
if N1 == N2:
887
N = N1
888
else:
889
raise RuntimeError, "BUG! Cardinality with pari=%s but with bsgs=%s"%(N1, N2)
890
else:
891
raise ValueError, "Algorithm is not known"
892
self._order = Integer(N)
893
return self._order
894
895
# now k is not a prime field and j is not 0, 1728
896
897
# we count points on a standard curve with the same
898
# j-invariant, defined over the field it generates, then lift
899
# to the curve's own field, and finally allow for twists
900
901
# Since j is not 0, 1728 the only twists are quadratic
902
903
j_pol=j.minimal_polynomial()
904
j_deg=j_pol.degree()
905
906
# if not possible to work over a smaller field:
907
if d==j_deg:
908
self._order = self.cardinality_bsgs()
909
return self._order
910
911
kj=GF(p**j_deg,name='a',modulus=j_pol)
912
jkj=kj.gen() if j_deg>1 else j_pol.roots(multiplicities=False)[0]
913
914
# recursive call which will do all the real work:
915
Ej = EllipticCurve_from_j(jkj)
916
N=Ej.cardinality(extension_degree=d//j_deg)
917
918
# if curve ia a (quadratic) twist of the "standard" one:
919
if not self.is_isomorphic(EllipticCurve_from_j(j)): N=2*(q+1)-N
920
921
self._order = N
922
return self._order
923
924
order = cardinality # alias
925
926
def frobenius_polynomial(self):
927
r"""
928
Return the characteristic polynomial of Frobenius.
929
930
The Frobenius endomorphism of the elliptic curve has quadratic
931
characteristic polynomial. In most cases this is irreducible and
932
defines an imaginary quadratic order; for some supersingular
933
curves, Frobenius is an integer a and the polynomial is
934
`(x-a)^2`.
935
936
.. note::
937
938
This computes the curve cardinality, which may be
939
time-consuming.
940
941
EXAMPLES::
942
943
sage: E=EllipticCurve(GF(11),[3,3])
944
sage: E.frobenius_polynomial()
945
x^2 - 4*x + 11
946
947
For some supersingular curves, Frobenius is in Z and the polynomial
948
is a square::
949
950
sage: E=EllipticCurve(GF(25,'a'),[0,0,0,0,1])
951
sage: E.frobenius_polynomial().factor()
952
(x + 5)^2
953
"""
954
x=polygen(ZZ)
955
return x**2-self.trace_of_frobenius()*x+self.base_field().cardinality()
956
957
def frobenius_order(self):
958
r"""
959
Return the quadratic order Z[phi] where phi is the Frobenius
960
endomorphism of the elliptic curve
961
962
.. note::
963
964
This computes the curve cardinality, which may be
965
time-consuming.
966
967
EXAMPLES::
968
969
sage: E=EllipticCurve(GF(11),[3,3])
970
sage: E.frobenius_order()
971
Order in Number Field in phi with defining polynomial x^2 - 4*x + 11
972
973
For some supersingular curves, Frobenius is in Z and the Frobenius
974
order is Z::
975
976
sage: E=EllipticCurve(GF(25,'a'),[0,0,0,0,1])
977
sage: R=E.frobenius_order()
978
sage: R
979
Order in Number Field in phi with defining polynomial x + 5
980
sage: R.degree()
981
1
982
"""
983
f = self.frobenius_polynomial().factor()[0][0]
984
return ZZ.extension(f,names='phi')
985
986
def frobenius(self):
987
r"""
988
Return the frobenius of self as an element of a quadratic order
989
990
.. note::
991
992
This computes the curve cardinality, which may be
993
time-consuming.
994
995
Frobenius is only determined up to conjugacy.
996
997
EXAMPLES::
998
999
sage: E=EllipticCurve(GF(11),[3,3])
1000
sage: E.frobenius()
1001
phi
1002
sage: E.frobenius().minpoly()
1003
x^2 - 4*x + 11
1004
1005
For some supersingular curves, Frobenius is in Z::
1006
1007
sage: E=EllipticCurve(GF(25,'a'),[0,0,0,0,1])
1008
sage: E.frobenius()
1009
-5
1010
"""
1011
R = self.frobenius_order()
1012
if R.degree()==1:
1013
return self.frobenius_polynomial().roots(multiplicities=False)[0]
1014
else:
1015
return R.gen(1)
1016
1017
def cardinality_exhaustive(self):
1018
r"""
1019
Return the cardinality of self over the base field. Simply adds up
1020
the number of points with each x-coordinate: only used for small
1021
field sizes!
1022
1023
EXAMPLES::
1024
1025
sage: p=next_prime(10^3)
1026
sage: E=EllipticCurve(GF(p),[3,4])
1027
sage: E.cardinality_exhaustive()
1028
1020
1029
sage: E=EllipticCurve(GF(3^4,'a'),[1,1])
1030
sage: E.cardinality_exhaustive()
1031
64
1032
"""
1033
self._order = Integer(1+sum([len(self.lift_x(x,all=True)) for x in self.base_field()]))
1034
return self._order
1035
1036
def cardinality_pari(self):
1037
r"""
1038
Return the cardinality of self over the (prime) base field using PARI.
1039
1040
The result is not cached.
1041
1042
EXAMPLES::
1043
1044
sage: p=next_prime(10^3)
1045
sage: E=EllipticCurve(GF(p),[3,4])
1046
sage: E.cardinality_pari()
1047
1020
1048
sage: K=GF(next_prime(10^6))
1049
sage: E=EllipticCurve(K,[1,0,0,1,1])
1050
sage: E.cardinality_pari()
1051
999945
1052
1053
TESTS::
1054
1055
sage: K.<a>=GF(3^20)
1056
sage: E=EllipticCurve(K,[1,0,0,1,a])
1057
sage: E.cardinality_pari()
1058
Traceback (most recent call last):
1059
...
1060
ValueError: cardinality_pari() only works over prime fields.
1061
sage: E.cardinality()
1062
3486794310
1063
1064
"""
1065
k = self.base_ring()
1066
p = k.characteristic()
1067
if k.degree()==1:
1068
return ZZ(p + 1 - int(self._pari_().ellap(p)))
1069
else:
1070
raise ValueError, "cardinality_pari() only works over prime fields."
1071
1072
def cardinality_bsgs(self, verbose=False):
1073
r"""
1074
Return the cardinality of self over the base field. Will be called
1075
by user function cardinality only when necessary, i.e. when the
1076
j_invariant is not in the prime field.
1077
1078
ALGORITHM: A variant of "Mestre's trick" extended to all finite
1079
fields by Cremona and Sutherland, 2008.
1080
1081
.. note::
1082
1083
1. The Mestre-Schoof-Cremona-Sutherland algorithm may fail for
1084
a small finite number of curves over `F_q` for `q` at most 49, so
1085
for `q<50` we use an exhaustive count.
1086
1087
2. Quadratic twists are not implemented in characteristic 2
1088
when `j=0 (=1728)`; but this case is treated separately.
1089
1090
EXAMPLES::
1091
1092
sage: p=next_prime(10^3)
1093
sage: E=EllipticCurve(GF(p),[3,4])
1094
sage: E.cardinality_bsgs()
1095
1020
1096
sage: E=EllipticCurve(GF(3^4,'a'),[1,1])
1097
sage: E.cardinality_bsgs()
1098
64
1099
sage: F.<a>=GF(101^3,'a')
1100
sage: E=EllipticCurve([2*a^2 + 48*a + 27, 89*a^2 + 76*a + 24])
1101
sage: E.cardinality_bsgs()
1102
1031352
1103
"""
1104
E1 = self
1105
k = self.base_field()
1106
q = k.order()
1107
if q<50:
1108
if verbose:
1109
print "q=",q,"< 50 so using exhaustive count"
1110
return self.cardinality_exhaustive()
1111
1112
# Construct the quadratic twist:
1113
E2 = E1.quadratic_twist()
1114
if verbose:
1115
print "Quadratic twist is ",E2.ainvs()
1116
1117
bounds = Hasse_bounds(q)
1118
lower, upper = bounds
1119
B = upper-q-1 # = floor(2*sqrt(q))
1120
a = ZZ(0)
1121
N1 = N2 = M = ZZ(1)
1122
kmin = -B
1123
kmax = B
1124
q1 = q+1
1125
# Throughout, we have #E=q+1-t where |t|<=B and t=a+k*M = a
1126
# (mod M) where kmin <= k <= kmax.
1127
1128
# M is the lcm of the orders of all the points found on E1 and
1129
# E2, which will eventually exceed 2*B, at which point
1130
# kmin=kmax.
1131
1132
if q > 2**10:
1133
N1 = ZZ(2)**sum([e for P,e in E1._p_primary_torsion_basis(2)])
1134
N2 = ZZ(2)**sum([e for P,e in E2._p_primary_torsion_basis(2)])
1135
if q > 2**20:
1136
N1 *= ZZ(3)**sum([e for P,e in E1._p_primary_torsion_basis(3)])
1137
N2 *= ZZ(3)**sum([e for P,e in E2._p_primary_torsion_basis(3)])
1138
if q > 2**40:
1139
N1 *= ZZ(5)**sum([e for P,e in E1._p_primary_torsion_basis(5)])
1140
N2 *= ZZ(5)**sum([e for P,e in E2._p_primary_torsion_basis(5)])
1141
# We now know that t=q+1 (mod N1) and t=-(q+1) (mod N2)
1142
a = q1
1143
M = N1
1144
g,u,v = M.xgcd(N2) # g==u*M+v*N2
1145
if N2>g:
1146
a = (a*v*N2-q1*u*M)//g
1147
M *= (N2//g) # = lcm(M,N2)
1148
a = a%M
1149
if verbose:
1150
print "(a,M)=",(a,M)
1151
kmin = ((-B-a)/M).ceil()
1152
kmax = ((B-a)/M).floor()
1153
if kmin==kmax:
1154
self._order = q1-a-kmin*M
1155
if verbose: print "no random points were needed"
1156
return self._order
1157
if verbose: print "(2,3,5)-torsion subgroup gives M=",M
1158
1159
# N1, N2 are divisors of the orders of E1, E2 separately,
1160
# which are used to speed up the computation of the orders of
1161
# points on each curve. For large q it is worth initializing
1162
# these with the full order of the (2,3,5)-torsion which are
1163
# often non-trivial.
1164
1165
while kmax!=kmin:
1166
# Get a random point on E1 and find its order, using the
1167
# Hasse bounds and the fact that we know that the group
1168
# order is a multiple of N1:
1169
n = generic.order_from_bounds(E1.random_point(),bounds,N1,operation='+')
1170
if verbose: print "New point on E has order ",n
1171
# update N1 and M
1172
N1 = N1.lcm(n)
1173
g,u,v = M.xgcd(n) # g==u*M+v*n
1174
if n>g:
1175
# update congruence a (mod M) with q+1 (mod n)
1176
a = (a*v*n+q1*u*M)//g
1177
M *= (n//g) # = lcm(M,n)
1178
a = a%M
1179
if verbose: print "(a,M)=",(a,M)
1180
kmin = ((-B-a)/M).ceil()
1181
kmax = ((B-a)/M).floor()
1182
if kmin==kmax:
1183
self._order = q1-a-kmin*M
1184
return self._order
1185
if verbose: print "number of possibilities is now ",kmax-kmin+1
1186
1187
# Get a random point on E2 and find its order, using the
1188
# Hasse bounds and the fact that we know that the group
1189
# order is a multiple of N2:
1190
n = generic.order_from_bounds(E2.random_point(),bounds,N2,operation='+')
1191
if verbose: print "New point on E' has order ",n
1192
# update N2 and M
1193
N2 = N2.lcm(n)
1194
g,u,v = M.xgcd(n) # g==u*M+v*n
1195
if n>g:
1196
# update congruence a (mod M) with -(q+1) (mod n)
1197
a = (a*v*n-q1*u*M)//g
1198
M *= (n//g) # = lcm(M,n)
1199
a = a%M
1200
if verbose: print "(a,M)=",(a,M)
1201
kmin = ((-B-a)/M).ceil()
1202
kmax = ((B-a)/M).floor()
1203
if kmin==kmax:
1204
self._order = q1-a-kmin*M
1205
return self._order
1206
if verbose: print "number of possibilities is now ",kmax-kmin+1
1207
1208
def gens(self):
1209
"""
1210
Returns a tuple of length up to 2 of points which generate the
1211
abelian group of points on this elliptic curve. See
1212
abelian_group() for limitations.
1213
1214
The algorithm uses random points on the curve, and hence the
1215
generators are likely to differ from one run to another; but they
1216
are cached so will be consistent in any one run of Sage.
1217
1218
AUTHORS:
1219
1220
- John Cremona
1221
1222
EXAMPLES::
1223
1224
sage: E=EllipticCurve(GF(11),[2,5])
1225
sage: E.gens() # random output
1226
((0 : 7 : 1),)
1227
sage: EllipticCurve(GF(41),[2,5]).gens() # random output
1228
((21 : 1 : 1), (8 : 0 : 1))
1229
sage: F.<a>=GF(3^6,'a')
1230
sage: E=EllipticCurve([a,a+1])
1231
sage: pts=E.gens()
1232
sage: len(pts)
1233
1
1234
sage: pts[0].order()==E.cardinality()
1235
True
1236
"""
1237
try:
1238
G = self.abelian_group()
1239
return [x.element() for x in G.gens()]
1240
except AttributeError:
1241
pass
1242
1243
def __iter__(self):
1244
"""
1245
Return an iterator through the points of this elliptic curve.
1246
1247
EXAMPLES::
1248
1249
sage: E = EllipticCurve(GF(11), [1,2])
1250
sage: for P in E: print P, P.order()
1251
(0 : 1 : 0) 1
1252
(1 : 2 : 1) 4
1253
(1 : 9 : 1) 4
1254
(2 : 1 : 1) 8
1255
...
1256
(10 : 0 : 1) 2
1257
"""
1258
for P in self.points():
1259
yield P
1260
1261
def __getitem__(self, n):
1262
"""
1263
Return the n'th point in self's __points list. This enables users
1264
to iterate over the curve's point set.
1265
1266
EXAMPLE::
1267
1268
sage: E=EllipticCurve(GF(97),[2,3])
1269
sage: S=E.points()
1270
sage: E[10]
1271
(10 : 76 : 1)
1272
sage: E[15]
1273
(17 : 10 : 1)
1274
sage: for P in E: print P.order()
1275
1
1276
50
1277
50
1278
50
1279
50
1280
5
1281
5
1282
50
1283
...
1284
"""
1285
return self.points()[n]
1286
1287
def abelian_group(self, debug=False):
1288
r"""
1289
Returns the abelian group structure of the group of points on this
1290
elliptic curve.
1291
1292
.. warning::
1293
1294
The algorithm is definitely *not* intended for use with
1295
*large* finite fields! The factorization of the orders of
1296
elements must be feasible. Also, baby-step-giant-step
1297
methods are used which have space and time requirements
1298
which are `O(\sqrt{q})`.
1299
1300
Also, the algorithm uses random points on the curve and hence the
1301
generators are likely to differ from one run to another; but the
1302
group is cached so the generators will not change in any one run of
1303
Sage.
1304
1305
INPUT:
1306
1307
1308
- ``debug`` - (default: False): if True, print
1309
debugging messages
1310
1311
1312
OUTPUT:
1313
1314
- an abelian group
1315
1316
- tuple of images of each of the generators of the abelian
1317
group as points on this curve
1318
1319
AUTHORS:
1320
1321
- John Cremona
1322
1323
EXAMPLES::
1324
1325
sage: E=EllipticCurve(GF(11),[2,5])
1326
sage: E.abelian_group()
1327
Additive abelian group isomorphic to Z/10 embedded in Abelian group of points on Elliptic Curve defined by y^2 = x^3 + 2*x + 5 over Finite Field of size 11
1328
1329
::
1330
1331
sage: E=EllipticCurve(GF(41),[2,5])
1332
sage: E.abelian_group()
1333
Additive abelian group isomorphic to Z/2 + Z/22 ...
1334
1335
::
1336
1337
sage: F.<a>=GF(3^6,'a')
1338
sage: E=EllipticCurve([a^4 + a^3 + 2*a^2 + 2*a, 2*a^5 + 2*a^3 + 2*a^2 + 1])
1339
sage: E.abelian_group()
1340
Additive abelian group isomorphic to Z/26 + Z/26 ...
1341
1342
::
1343
1344
sage: F.<a>=GF(101^3,'a')
1345
sage: E=EllipticCurve([2*a^2 + 48*a + 27, 89*a^2 + 76*a + 24])
1346
sage: E.abelian_group()
1347
Additive abelian group isomorphic to Z/1031352 ...
1348
1349
The group can be trivial::
1350
1351
sage: E=EllipticCurve(GF(2),[0,0,1,1,1])
1352
sage: E.abelian_group()
1353
Trivial group embedded in Abelian group of points on Elliptic Curve defined by y^2 + y = x^3 + x + 1 over Finite Field of size 2
1354
1355
Of course, there are plenty of points if we extend the field::
1356
1357
sage: E.cardinality(extension_degree=100)
1358
1267650600228231653296516890625
1359
1360
This tests the patch for trac #3111, using 10 primes randomly
1361
selected::
1362
1363
sage: E = EllipticCurve('389a')
1364
sage: for p in [5927, 2297, 1571, 1709, 3851, 127, 3253, 5783, 3499, 4817]:
1365
... G = E.change_ring(GF(p)).abelian_group()
1366
sage: for p in prime_range(10000): # long time (19s on sage.math, 2011)
1367
... if p != 389:
1368
... G = E.change_ring(GF(p)).abelian_group()
1369
1370
This tests that the bug reported in trac #3926 has been fixed::
1371
1372
sage: K.<i> = QuadraticField(-1)
1373
sage: OK = K.ring_of_integers()
1374
sage: P=K.factor(10007)[0][0]
1375
sage: OKmodP = OK.residue_field(P)
1376
sage: E = EllipticCurve([0,0,0,i,i+3])
1377
sage: Emod = E.change_ring(OKmodP); Emod
1378
Elliptic Curve defined by y^2 = x^3 + ibar*x + (ibar+3) over Residue field in ibar of Fractional ideal (10007)
1379
sage: Emod.abelian_group() #random generators
1380
(Multiplicative Abelian Group isomorphic to C50067594 x C2,
1381
((3152*ibar + 7679 : 7330*ibar + 7913 : 1), (8466*ibar + 1770 : 0 : 1)))
1382
"""
1383
if not debug:
1384
# if we're in debug mode, always recalculate
1385
try:
1386
return self.__abelian_group
1387
except AttributeError:
1388
pass
1389
1390
k = self.base_field()
1391
q = k.order()
1392
p = k.characteristic()
1393
d = k.degree()
1394
j = self.j_invariant()
1395
if d>1:
1396
d = j.minimal_polynomial().degree()
1397
1398
1399
# Before computing the group structure we compute the
1400
# cardinality. While this is not strictly necessary, it makes
1401
# the code simpler and also makes computation of orders of
1402
# points faster.
1403
1404
# j=0,1728
1405
1406
if j==k(0):
1407
N = self._cardinality_with_j_invariant_0()
1408
if j==k(1728):
1409
N = self._cardinality_with_j_invariant_1728()
1410
1411
bounds = Hasse_bounds(q)
1412
lower, upper = bounds
1413
if debug:
1414
print "Lower and upper bounds on group order: [",lower,",",upper,"]"
1415
1416
try:
1417
N=self._order
1418
if debug:
1419
print "Group order already known to be ",N
1420
except:
1421
if (q<50):
1422
if debug:
1423
print "Computing group order naively"
1424
N=self.cardinality_exhaustive()
1425
elif d==1:
1426
if debug:
1427
print "Computing group order using PARI"
1428
N=self.cardinality()
1429
else:
1430
if debug:
1431
print "Computing group order using bsgs"
1432
N=self.cardinality_bsgs()
1433
if debug:
1434
print "... group order = ",N
1435
1436
self._order=N
1437
plist = N.prime_factors()
1438
P1=self(0)
1439
P2=self(0)
1440
n1= Integer(1)
1441
n2= Integer(1)
1442
P1._order=n1
1443
P2._order=n2
1444
npts = 0
1445
1446
# At all stages the current subgroup is generated by P1, P2 with
1447
# orders n1,n2 which are disjoint. We stop when n1*n2 == N
1448
1449
if debug:
1450
"About to start generating random points"
1451
1452
while n1*n2 != N:
1453
if debug:
1454
"Getting a new random point"
1455
Q = self.random_point()
1456
while Q.is_zero(): Q = self.random_point()
1457
npts += 1
1458
if debug:
1459
print "Q = ",Q,":",
1460
print " Order(Q) = ", Q.order()
1461
1462
Q1=n1*Q;
1463
1464
if Q1.is_zero() and npts>=10: # then P1,n1 will not change but we may increase n2
1465
if debug: print "Case 2: n2 may increase"
1466
n1a = 1; n1b = n1
1467
P1a = P1
1468
n1a = n1.prime_to_m_part(N//n1)
1469
n1b = n1//n1a
1470
Q = n1a*Q # has order | n1b
1471
P1a = n1a*P1 # has order = n1b
1472
if debug: print "n1a=",n1a
1473
a = None
1474
for m in n1b.divisors():
1475
try:
1476
a = generic.bsgs(m*P1a,m*Q,(0,(n1b//m)-1),operation='+')
1477
break
1478
except ValueError:
1479
pass
1480
assert a != None
1481
a *= (m*n1a)
1482
if debug: print "linear relation gives m=",m,", a=",a
1483
if debug: assert m*Q==a*P1
1484
if m>1: # else Q is in <P1>
1485
Q=Q-(a//m)*P1; # has order m and is disjoint from P1
1486
if debug: assert Q.order()==m
1487
Q._order=m
1488
if n2==1: # this is our first nontrivial P2
1489
P2=Q
1490
n2=m
1491
if debug:
1492
print "Adding second generator ",P2," of order ",n2
1493
print "Subgroup order now ",n1*n2,"=",n1,"*",n2
1494
else: # we must merge P2 and Q:
1495
oldn2=n2 # holds old value
1496
P2,n2=generic.merge_points((P2,n2),(Q,m),operation='+', check=debug)
1497
if debug: assert P2.order()==n2
1498
P2._order=n2
1499
if debug:
1500
if n2>oldn2:
1501
print "Replacing second generator by ",P2,
1502
print " of order ",n2, " gaining index ",n2//oldn2
1503
print "Subgroup order now ",n1*n2,"=",n1,"*",n2
1504
elif not Q1.is_zero(): # Q1 nonzero: n1 will increase
1505
if debug: print "Case 1: n1 may increase"
1506
oldn1=n1
1507
if n2>1:
1508
P3=(n1//n2)*P1 # so P2,P3 are a basis for n2-torsion
1509
if debug: assert P3.order()==n2
1510
P3._order=n2
1511
if debug: print "storing generator ",P3," of ",n2,"-torsion"
1512
m = generic.order_from_multiple(Q,N,plist,operation='+', check=debug)
1513
P1,n1=generic.merge_points((P1,n1),(Q,m), check=debug)
1514
if debug: assert P1.order()==n1
1515
P1._order=n1
1516
if debug:
1517
print "Replacing first generator by ",P1," of order ",
1518
print n1,", gaining index ",n1//oldn1
1519
print "Subgroup order now ",n1*n2,"=",n1,"*",n2
1520
1521
# Now replace P2 by a point of order n2 s.t. it and
1522
# (n1//n2)*P1 are still a basis for n2-torsion:
1523
if n2>1:
1524
a,m = generic.linear_relation(P1,P3,operation='+')
1525
if debug: print "linear relation gives m=",m,", a=",a
1526
P3=P3-(a//m)*P1
1527
if debug: assert P3.order()==m
1528
P3._order=m
1529
if debug: print "First P2 component =",P3
1530
if m==n2:
1531
P2=P3
1532
else:
1533
a,m = generic.linear_relation(P1,P2,operation='+')
1534
if debug: print "linear relation gives m=",m,", a=",a
1535
P2=P2-(a//m)*P1;
1536
if debug: assert P2.order()==m
1537
P2._order=m
1538
if debug: print "Second P2 component =",P2
1539
P2,n2=generic.merge_points((P2,n2),(P3,m), check=debug)
1540
if debug: assert P2.order()==n2
1541
P2._order=n2
1542
if debug: print "Combined P2 component =",P2
1543
1544
if debug:
1545
if P1.order()!=n1:
1546
print "Generator P1 = ",P1," has order ",P1.order(),
1547
print " and not ",n1
1548
raise ValueError
1549
if P2.order()!=n2:
1550
print "Generator P2 = ",P2," has order ",P2.order()
1551
print " and not ",n2
1552
raise ValueError
1553
if n2>1:
1554
if generic.linear_relation(P1,P2,operation='+')[1]!=n2:
1555
print "Generators not independent!"
1556
raise ValueError
1557
print "Generators: P1 = ",P1," of order ",n1,
1558
print ", P2 = ",P2," of order ",n2
1559
print "Subgroup order is now ",n1*n2,"=",n1,"*",n2
1560
1561
# Finished: record group order, structure and generators
1562
1563
from sage.groups.additive_abelian.additive_abelian_wrapper import AdditiveAbelianGroupWrapper
1564
self._order = n1*n2
1565
if n1==1:
1566
self.__abelian_group = AdditiveAbelianGroupWrapper(self.point_homset(), [], [])
1567
else:
1568
if n2==1:
1569
self.__abelian_group = AdditiveAbelianGroupWrapper(self.point_homset(), [P1], [n1])
1570
else:
1571
self.__abelian_group = AdditiveAbelianGroupWrapper(self.point_homset(), [P1, P2], [n1, n2])
1572
return self.__abelian_group
1573
1574
def is_isogenous(self, other, field=None, proof=True):
1575
"""
1576
Returns whether or not self is isogenous to other
1577
1578
INPUT:
1579
1580
- ``other`` -- another elliptic curve.
1581
1582
- ``field`` (default None) -- a field containing the base
1583
fields of the two elliptic curves into which the two curves
1584
may be extended to test if they are isogenous over this
1585
field. By default is_isogenous will not try to find this
1586
field unless one of the curves can be extended into the base
1587
field of the other, in which case it will test over the
1588
larger base field.
1589
1590
- ``proof`` (default True) -- this parameter is here only to
1591
be consistent with versions for other types of elliptic
1592
curves.
1593
1594
OUTPUT:
1595
1596
(bool) True if there is an isogeny from curve ``self`` to
1597
curve ``other`` defined over ``field``.
1598
1599
EXAMPLES::
1600
1601
sage: E1 = EllipticCurve(GF(11^2,'a'),[2,7]); E1
1602
Elliptic Curve defined by y^2 = x^3 + 2*x + 7 over Finite Field in a of size 11^2
1603
sage: E1.is_isogenous(5)
1604
Traceback (most recent call last):
1605
...
1606
ValueError: Second argument is not an Elliptic Curve.
1607
sage: E1.is_isogenous(E1)
1608
True
1609
1610
sage: E2 = EllipticCurve(GF(7^3,'b'),[3,1]); E2
1611
Elliptic Curve defined by y^2 = x^3 + 3*x + 1 over Finite Field in b of size 7^3
1612
sage: E1.is_isogenous(E2)
1613
Traceback (most recent call last):
1614
...
1615
ValueError: The base fields must have the same characteristic.
1616
1617
sage: E3 = EllipticCurve(GF(11^2,'c'),[4,3]); E3
1618
Elliptic Curve defined by y^2 = x^3 + 4*x + 3 over Finite Field in c of size 11^2
1619
sage: E1.is_isogenous(E3)
1620
False
1621
1622
sage: E4 = EllipticCurve(GF(11^6,'d'),[6,5]); E4
1623
Elliptic Curve defined by y^2 = x^3 + 6*x + 5 over Finite Field in d of size 11^6
1624
sage: E1.is_isogenous(E4)
1625
True
1626
1627
sage: E5 = EllipticCurve(GF(11^7,'e'),[4,2]); E5
1628
Elliptic Curve defined by y^2 = x^3 + 4*x + 2 over Finite Field in e of size 11^7
1629
sage: E1.is_isogenous(E5)
1630
Traceback (most recent call last):
1631
...
1632
ValueError: Curves have different base fields: use the field parameter.
1633
1634
When the field is given:
1635
1636
sage: E1 = EllipticCurve(GF(13^2,'a'),[2,7]); E1
1637
Elliptic Curve defined by y^2 = x^3 + 2*x + 7 over Finite Field in a of size 13^2
1638
sage: E1.is_isogenous(5,GF(13^6,'f'))
1639
Traceback (most recent call last):
1640
...
1641
ValueError: Second argument is not an Elliptic Curve.
1642
sage: E6 = EllipticCurve(GF(11^3,'g'),[9,3]); E6
1643
Elliptic Curve defined by y^2 = x^3 + 9*x + 3 over Finite Field in g of size 11^3
1644
sage: E1.is_isogenous(E6,QQ)
1645
Traceback (most recent call last):
1646
...
1647
ValueError: The base fields must have the same characteristic.
1648
sage: E7 = EllipticCurve(GF(13^5,'h'),[2,9]); E7
1649
Elliptic Curve defined by y^2 = x^3 + 2*x + 9 over Finite Field in h of size 13^5
1650
sage: E1.is_isogenous(E7,GF(13^4,'i'))
1651
Traceback (most recent call last):
1652
...
1653
ValueError: Field must be an extension of the base fields of both curves
1654
sage: E1.is_isogenous(E7,GF(13^10,'j'))
1655
False
1656
sage: E1.is_isogenous(E7,GF(13^30,'j'))
1657
False
1658
"""
1659
from ell_generic import is_EllipticCurve
1660
if not is_EllipticCurve(other):
1661
raise ValueError, "Second argument is not an Elliptic Curve."
1662
if self.is_isomorphic(other):
1663
return True
1664
elif self.base_field().characteristic() != other.base_field().characteristic():
1665
raise ValueError, "The base fields must have the same characteristic."
1666
elif field==None:
1667
if self.base_field().degree() == other.base_field().degree():
1668
if self.cardinality() == other.cardinality():
1669
return True
1670
else:
1671
return False
1672
elif self.base_field().degree() == gcd(self.base_field().degree(),other.base_field().degree()):
1673
if self.cardinality(extension_degree=other.base_field().degree()//self.base_field().degree()) == other.cardinality():
1674
return True
1675
else:
1676
return False
1677
elif other.base_field().degree() == gcd(self.base_field().degree(),other.base_field().degree()):
1678
if other.cardinality(extension_degree=self.base_field().degree()//other.base_field().degree()) == self.cardinality():
1679
return True
1680
else:
1681
return False
1682
else:
1683
raise ValueError, "Curves have different base fields: use the field parameter."
1684
else:
1685
if not lcm(self.base_field().degree(), other.base_field().degree()).divides(field.degree()):
1686
raise ValueError, "Field must be an extension of the base fields of both curves"
1687
else:
1688
if \
1689
self.cardinality(extension_degree=field.degree()//self.base_field().degree())\
1690
== other.cardinality(extension_degree=field.degree()//other.base_field().degree()):
1691
return True
1692
else:
1693
return False
1694
1695
def is_supersingular(self, proof=True):
1696
r"""
1697
Return True if this elliptic curve is supersingular, else False.
1698
1699
INPUT:
1700
1701
- ``proof`` (boolean, default True) -- If True, returns a
1702
proved result. If False, then a return value of False is
1703
certain but a return value of True may be based on a
1704
probabilistic test. See the documentaion of the function
1705
:meth:`is_j_supersingular` for more details.
1706
1707
EXAMPLES::
1708
1709
sage: F = GF(101)
1710
sage: EllipticCurve(j=F(0)).is_supersingular()
1711
True
1712
sage: EllipticCurve(j=F(1728)).is_supersingular()
1713
False
1714
sage: EllipticCurve(j=F(66)).is_supersingular()
1715
True
1716
sage: EllipticCurve(j=F(99)).is_supersingular()
1717
False
1718
1719
TESTS::
1720
1721
sage: from sage.schemes.elliptic_curves.ell_finite_field import supersingular_j_polynomial, is_j_supersingular
1722
sage: F = GF(103)
1723
sage: ssjlist = [F(1728)] + supersingular_j_polynomial(103).roots(multiplicities=False)
1724
sage: Set([j for j in F if is_j_supersingular(j)]) == Set(ssjlist)
1725
True
1726
1727
"""
1728
return is_j_supersingular(self.j_invariant(), proof=proof)
1729
1730
def is_ordinary(self, proof=True):
1731
r"""
1732
Return True if this elliptic curve is ordinary, else False.
1733
1734
INPUT:
1735
1736
- ``proof`` (boolean, default True) -- If True, returns a
1737
proved result. If False, then a return value of True is
1738
certain but a return value of False may be based on a
1739
probabilistic test. See the documentaion of the function
1740
:meth:`is_j_supersingular` for more details.
1741
1742
EXAMPLES::
1743
1744
sage: F = GF(101)
1745
sage: EllipticCurve(j=F(0)).is_ordinary()
1746
False
1747
sage: EllipticCurve(j=F(1728)).is_ordinary()
1748
True
1749
sage: EllipticCurve(j=F(66)).is_ordinary()
1750
False
1751
sage: EllipticCurve(j=F(99)).is_ordinary()
1752
True
1753
1754
"""
1755
return not is_j_supersingular(self.j_invariant(), proof=proof)
1756
1757
def set_order(self, value, num_checks=8):
1758
r"""
1759
Set the value of self._order to value.
1760
1761
Use this when you know a priori the order of the curve to
1762
avoid a potentially expensive order calculation.
1763
1764
INPUT:
1765
1766
- ``value`` - Integer in the Hasse-Weil range for this
1767
curve.
1768
1769
- ``num_checks`` - Integer (default: 8) number of times to
1770
check whether value*(a random point on this curve) is
1771
equal to the identity.
1772
1773
1774
OUTPUT:
1775
1776
None
1777
1778
EXAMPLES:
1779
1780
This example illustrates basic usage.
1781
1782
::
1783
1784
sage: E = EllipticCurve(GF(7), [0, 1]) # This curve has order 6
1785
sage: E.set_order(6)
1786
sage: E.order()
1787
6
1788
sage: E.order() * E.random_point()
1789
(0 : 1 : 0)
1790
1791
We now give a more interesting case, the NIST-P521 curve. Its
1792
order is too big to calculate with Sage, and takes a long time
1793
using other packages, so it is very useful here.
1794
1795
::
1796
1797
sage: p = 2^521 - 1
1798
sage: prev_proof_state = proof.arithmetic()
1799
sage: proof.arithmetic(False) # turn off primality checking
1800
sage: F = GF(p)
1801
sage: A = p - 3
1802
sage: B = 1093849038073734274511112390766805569936207598951683748994586394495953116150735016013708737573759623248592132296706313309438452531591012912142327488478985984
1803
sage: q = 6864797660130609714981900799081393217269435300143305409394463459185543183397655394245057746333217197532963996371363321113864768612440380340372808892707005449
1804
sage: E = EllipticCurve([F(A), F(B)])
1805
sage: E.set_order(q)
1806
sage: G = E.random_point()
1807
sage: E.order() * G # This takes practically no time.
1808
(0 : 1 : 0)
1809
sage: proof.arithmetic(prev_proof_state) # restore state
1810
1811
It is an error to pass a value which is not an integer in the
1812
Hasse-Weil range::
1813
1814
sage: E = EllipticCurve(GF(7), [0, 1]) # This curve has order 6
1815
sage: E.set_order("hi")
1816
Traceback (most recent call last):
1817
...
1818
ValueError: Value hi illegal (not an integer in the Hasse range)
1819
sage: E.set_order(3.14159)
1820
Traceback (most recent call last):
1821
...
1822
ValueError: Value 3.14159000000000 illegal (not an integer in the Hasse range)
1823
sage: E.set_order(0)
1824
Traceback (most recent call last):
1825
...
1826
ValueError: Value 0 illegal (not an integer in the Hasse range)
1827
sage: E.set_order(1000)
1828
Traceback (most recent call last):
1829
...
1830
ValueError: Value 1000 illegal (not an integer in the Hasse range)
1831
1832
It is also very likely an error to pass a value which is not
1833
the actual order of this curve. How unlikely is determined by
1834
num_checks, the factorization of the actual order, and the
1835
actual group structure::
1836
1837
sage: E = EllipticCurve(GF(7), [0, 1]) # This curve has order 6
1838
sage: E.set_order(11)
1839
Traceback (most recent call last):
1840
...
1841
ValueError: Value 11 illegal (multiple of random point not the identity)
1842
1843
However, set_order can be fooled, though it's not likely in
1844
"real cases of interest". For instance, the order can be set
1845
to a multiple of the actual order::
1846
1847
sage: E = EllipticCurve(GF(7), [0, 1]) # This curve has order 6
1848
sage: E.set_order(12) # 12 just fits in the Hasse range
1849
sage: E.order()
1850
12
1851
1852
Or, the order can be set incorrectly along with num_checks set
1853
too small::
1854
1855
sage: E = EllipticCurve(GF(7), [0, 1]) # This curve has order 6
1856
sage: E.set_order(4, num_checks=0)
1857
WARNING: No checking done in set_order
1858
sage: E.order()
1859
4
1860
1861
The value of num_checks must be an integer. Negative values
1862
are interpreted as zero, which means don't do any checking::
1863
1864
sage: E = EllipticCurve(GF(7), [0, 1]) # This curve has order 6
1865
sage: E.set_order(4, num_checks=-12)
1866
WARNING: No checking done in set_order
1867
sage: E.order()
1868
4
1869
1870
NOTES:
1871
1872
The implementation is based on the fact that orders of elliptic curves
1873
are cached in the (pseudo-private) _order slot.
1874
1875
AUTHORS:
1876
1877
- Mariah Lenox (2011-02-16)
1878
"""
1879
# Is value in the Hasse range?
1880
q = self.base_field().order()
1881
a,b = Hasse_bounds(q,1)
1882
#a = q + 1 - 2*q.isqrt()
1883
#b = q + 1 + 2*q.isqrt()
1884
if not value in ZZ:
1885
raise ValueError('Value %s illegal (not an integer in the Hasse range)'%value)
1886
if not a <= value <= b:
1887
raise ValueError('Value %s illegal (not an integer in the Hasse range)'%value)
1888
# Is value*random == identity?
1889
for i in range(num_checks):
1890
G = self.random_point()
1891
if value * G != self(0):
1892
raise ValueError('Value %s illegal (multiple of random point not the identity)'%value)
1893
if(num_checks <= 0):
1894
print 'WARNING: No checking done in set_order'
1895
self._order = value
1896
1897
def supersingular_j_polynomial(p):
1898
"""
1899
Return a polynomial whose roots are the supersingular `j`-invariants
1900
in characteristic `p`, other than 0, 1728.
1901
1902
INPUT:
1903
1904
- `p` (integer) -- a prime number.
1905
1906
ALGORITHM:
1907
1908
First compute H(X) whose roots are the Legendre
1909
`\lambda`-invariants of supersingular curves (Silverman V.4.1(b))
1910
in charactersitic `p`. Then, using a resultant computation with
1911
the polynomial relating `\lambda` and `j` (Silverman III.1.7(b)),
1912
we recover the polynomial (in variable ``j``) whose roots are the
1913
`j`-invariants. Factors of `j` and `j-1728` are removed if
1914
present.
1915
1916
EXAMPLES::
1917
1918
sage: from sage.schemes.elliptic_curves.ell_finite_field import supersingular_j_polynomial
1919
sage: f = supersingular_j_polynomial(67); f
1920
j^5 + 53*j^4 + 4*j^3 + 47*j^2 + 36*j + 8
1921
sage: f.factor()
1922
(j + 1) * (j^2 + 8*j + 45) * (j^2 + 44*j + 24)
1923
1924
::
1925
1926
sage: [supersingular_j_polynomial(p) for p in prime_range(30)]
1927
[1, 1, 1, 1, 1, j + 8, j + 9, j + 12, j + 4, j^2 + 2*j + 21]
1928
1929
TESTS::
1930
1931
sage: supersingular_j_polynomial(6)
1932
Traceback (most recent call last):
1933
...
1934
ValueError: p (=6) should be a prime number
1935
1936
"""
1937
try:
1938
p = ZZ(p)
1939
except TypeError:
1940
raise ValueError, "p (=%s) should be a prime number"%p
1941
if not p.is_prime():
1942
raise ValueError, "p (=%s) should be a prime number"%p
1943
1944
J = polygen(GF(p),'j')
1945
if p<13:
1946
return J.parent().one()
1947
from sage.rings.all import binomial
1948
from sage.misc.all import prod
1949
m=(p-1)//2
1950
X,T = PolynomialRing(GF(p),2,names=['X','T']).gens()
1951
H = sum([binomial(m,i)**2 * T**i for i in xrange(m+1)])
1952
F = T**2 * (T-1)**2 * X - 256*(T**2-T+1)**3
1953
R = F.resultant(H,T)
1954
R = prod([fi for fi,e in R([J,0]).factor()])
1955
if R(0)==0:
1956
R = R//J
1957
if R(1728)==0:
1958
R = R//(J-1728)
1959
return R
1960
1961
# For p in [13..300] we have precomputed these polynomials and store
1962
# them (as lists of their coefficients in ZZ) in a dict:
1963
1964
supersingular_j_polynomials = dict()
1965
1966
supersingular_j_polynomials[13] = [8, 1]
1967
supersingular_j_polynomials[17] = [9, 1]
1968
supersingular_j_polynomials[19] = [12, 1]
1969
supersingular_j_polynomials[23] = [4, 1]
1970
supersingular_j_polynomials[29] = [21, 2, 1]
1971
supersingular_j_polynomials[31] = [8, 25, 1]
1972
supersingular_j_polynomials[37] = [11, 5, 23, 1]
1973
supersingular_j_polynomials[41] = [18, 10, 19, 1]
1974
supersingular_j_polynomials[43] = [32, 11, 21, 1]
1975
supersingular_j_polynomials[47] = [35, 33, 31, 1]
1976
supersingular_j_polynomials[53] = [24, 9, 30, 7, 1]
1977
supersingular_j_polynomials[59] = [39, 31, 35, 39, 1]
1978
supersingular_j_polynomials[61] = [60, 21, 27, 8, 60, 1]
1979
supersingular_j_polynomials[67] = [8, 36, 47, 4, 53, 1]
1980
supersingular_j_polynomials[71] = [18, 54, 28, 33, 1, 1]
1981
supersingular_j_polynomials[73] = [7, 39, 38, 9, 68, 60, 1]
1982
supersingular_j_polynomials[79] = [10, 25, 1, 63, 57, 55, 1]
1983
supersingular_j_polynomials[83] = [43, 72, 81, 81, 62, 11, 1]
1984
supersingular_j_polynomials[89] = [42, 79, 23, 22, 37, 86, 60, 1]
1985
supersingular_j_polynomials[97] = [19, 28, 3, 72, 2, 96, 10, 60, 1]
1986
supersingular_j_polynomials[101] = [9, 76, 45, 79, 1, 68, 87, 60, 1]
1987
supersingular_j_polynomials[103] = [64, 15, 24, 58, 70, 83, 84, 100, 1]
1988
supersingular_j_polynomials[107] = [6, 18, 72, 59, 43, 19, 17, 68, 1]
1989
supersingular_j_polynomials[109] = [107, 22, 39, 83, 30, 34, 108, 104, 60, 1]
1990
supersingular_j_polynomials[113] = [86, 71, 75, 6, 47, 97, 100, 4, 60, 1]
1991
supersingular_j_polynomials[127] = [32, 31, 5, 50, 115, 122, 114, 67, 38, 35, 1]
1992
supersingular_j_polynomials[131] = [65, 64, 10, 34, 129, 35, 94, 127, 7, 7, 1]
1993
supersingular_j_polynomials[137] = [104, 83, 3, 82, 112, 23, 77, 135, 18, 50, 60, 1]
1994
supersingular_j_polynomials[139] = [87, 79, 109, 21, 138, 9, 104, 130, 61, 118, 90, 1]
1995
supersingular_j_polynomials[149] = [135, 55, 80, 86, 87, 74, 32, 60, 130, 80, 146, 60, 1]
1996
supersingular_j_polynomials[151] = [94, 125, 8, 6, 93, 21, 114, 80, 107, 58, 42, 18, 1]
1997
supersingular_j_polynomials[157] = [14, 95, 22, 58, 110, 23, 71, 51, 47, 5, 147, 59, 60, 1]
1998
supersingular_j_polynomials[163] = [102, 26, 74, 95, 112, 151, 98, 107, 27, 37, 25, 111, 109, 1]
1999
supersingular_j_polynomials[167] = [14, 9, 27, 109, 97, 55, 51, 74, 145, 125, 36, 113, 89, 1]
2000
supersingular_j_polynomials[173] = [152, 73, 56, 12, 18, 96, 98, 49, 30, 43, 52, 79, 163, 60, 1]
2001
supersingular_j_polynomials[179] = [110, 51, 3, 94, 123, 90, 156, 90, 88, 119, 158, 27, 71, 29, 1]
2002
supersingular_j_polynomials[181] = [7, 65, 77, 29, 139, 34, 65, 84, 164, 73, 51, 136, 7, 141, 60, 1]
2003
supersingular_j_polynomials[191] = [173, 140, 144, 3, 135, 80, 182, 84, 93, 75, 83, 17, 22, 42, 160, 1]
2004
supersingular_j_polynomials[193] = [23, 48, 26, 15, 108, 141, 124, 44, 132, 49, 72, 173, 126, 101, 22, 60, 1]
2005
supersingular_j_polynomials[197] = [14, 111, 64, 170, 193, 32, 124, 91, 112, 163, 14, 112, 167, 191, 183, 60, 1]
2006
supersingular_j_polynomials[199] = [125, 72, 65, 30, 63, 45, 10, 177, 91, 102, 28, 27, 5, 150, 51, 128, 1]
2007
supersingular_j_polynomials[211] = [27, 137, 128, 90, 102, 141, 5, 77, 131, 144, 83, 108, 23, 105, 98, 13, 80, 1]
2008
supersingular_j_polynomials[223] = [56, 183, 46, 133, 191, 94, 20, 8, 92, 100, 57, 200, 166, 67, 59, 218, 28, 32, 1]
2009
supersingular_j_polynomials[227] = [79, 192, 142, 66, 11, 114, 100, 208, 57, 147, 32, 5, 144, 93, 185, 147, 92, 16, 1]
2010
supersingular_j_polynomials[229] = [22, 55, 182, 130, 228, 172, 63, 25, 108, 99, 100, 101, 220, 111, 205, 199, 91, 163, 60, 1]
2011
supersingular_j_polynomials[233] = [101, 148, 85, 113, 226, 68, 71, 103, 61, 44, 173, 175, 5, 225, 227, 99, 146, 170, 60, 1]
2012
supersingular_j_polynomials[239] = [225, 81, 47, 26, 133, 182, 238, 2, 144, 154, 234, 178, 165, 130, 35, 61, 144, 112, 207, 1]
2013
supersingular_j_polynomials[241] = [224, 51, 227, 139, 134, 186, 187, 152, 161, 175, 213, 59, 105, 88, 87, 124, 202, 40, 15, 60, 1]
2014
supersingular_j_polynomials[251] = [30, 183, 80, 127, 40, 56, 230, 168, 192, 48, 226, 61, 214, 54, 165, 147, 105, 88, 38, 171, 1]
2015
supersingular_j_polynomials[257] = [148, 201, 140, 146, 169, 147, 220, 4, 205, 224, 35, 42, 198, 97, 127, 7, 110, 229, 118, 202, 60, 1]
2016
supersingular_j_polynomials[263] = [245, 126, 72, 213, 14, 64, 152, 83, 169, 114, 9, 128, 138, 231, 103, 85, 114, 211, 173, 249, 135, 1]
2017
supersingular_j_polynomials[269] = [159, 32, 69, 95, 201, 266, 190, 176, 76, 151, 212, 21, 106, 49, 263, 105, 136, 194, 215, 181, 237, 60, 1]
2018
supersingular_j_polynomials[271] = [169, 87, 179, 109, 133, 101, 31, 167, 208, 99, 127, 120, 83, 62, 36, 23, 61, 50, 69, 263, 265, 111, 1]
2019
supersingular_j_polynomials[277] = [251, 254, 171, 72, 190, 237, 12, 231, 123, 217, 263, 151, 270, 183, 29, 228, 85, 4, 67, 101, 29, 169, 60, 1]
2020
supersingular_j_polynomials[281] = [230, 15, 146, 69, 41, 23, 142, 232, 18, 80, 58, 134, 270, 62, 272, 70, 247, 189, 118, 255, 274, 159, 60, 1]
2021
supersingular_j_polynomials[283] = [212, 4, 42, 155, 38, 1, 270, 175, 172, 256, 264, 232, 50, 82, 244, 127, 148, 46, 249, 72, 59, 124, 75, 1]
2022
supersingular_j_polynomials[293] = [264, 66, 165, 144, 243, 25, 163, 210, 18, 107, 160, 153, 70, 255, 91, 211, 22, 7, 256, 50, 150, 94, 225, 60, 1]
2023
2024
2025
def is_j_supersingular(j, proof=True):
2026
r"""
2027
Return True if `j` is a supersingular `j`-invariant.
2028
2029
INPUT:
2030
2031
- ``j`` (finite field element) -- an element of a finite field
2032
2033
- ``proof`` (boolean, default True) -- If True, returns a proved
2034
result. If False, then a return value of False is certain but a
2035
return value of True may be based on a probabilistic test. See
2036
the ALGORITHM section below for more details.
2037
2038
OUTPUT:
2039
2040
(boolean) True if `j` is supersingular, else False.
2041
2042
ALGORITHM:
2043
2044
For small characteristics `p` we check whether the `j`-invariant
2045
is in a precomputed list of supersingular values. Otherwise we
2046
next check the `j`-invariant. If `j=0`, the curve is
2047
supersingular if and only if `p=2` or `p\equiv3\pmod{4}`; if
2048
`j=1728`, the curve is supersingular if and only if `p=3` or
2049
`p\equiv2\pmod{3}`. Next, if the base field is the prime field
2050
`{\rm GF}(p)`, we check that `(p+1)P=0` for several random points
2051
`P`, returning False if any fail: supersingular curves over `{\rm
2052
GF}(p)` have cardinality `p+1`. If Proof is false we now return
2053
True. Otherwise we compute the cardinality and return True if and
2054
only if it is divisible by `p`.
2055
2056
EXAMPLES::
2057
2058
sage: from sage.schemes.elliptic_curves.ell_finite_field import is_j_supersingular, supersingular_j_polynomials
2059
sage: [(p,[j for j in GF(p) if is_j_supersingular(j)]) for p in prime_range(30)]
2060
[(2, [0]), (3, [0]), (5, [0]), (7, [6]), (11, [0, 1]), (13, [5]), (17, [0, 8]), (19, [7, 18]), (23, [0, 3, 19]), (29, [0, 2, 25])]
2061
2062
sage: [j for j in GF(109) if is_j_supersingular(j)]
2063
[17, 41, 43]
2064
sage: PolynomialRing(GF(109),'j')(supersingular_j_polynomials[109]).roots()
2065
[(43, 1), (41, 1), (17, 1)]
2066
2067
sage: [p for p in prime_range(100) if is_j_supersingular(GF(p)(0))]
2068
[2, 3, 5, 11, 17, 23, 29, 41, 47, 53, 59, 71, 83, 89]
2069
sage: [p for p in prime_range(100) if is_j_supersingular(GF(p)(1728))]
2070
[2, 3, 7, 11, 19, 23, 31, 43, 47, 59, 67, 71, 79, 83]
2071
sage: [p for p in prime_range(100) if is_j_supersingular(GF(p)(123456))]
2072
[2, 3, 59, 89]
2073
2074
"""
2075
if not is_FiniteFieldElement(j):
2076
raise ValueError, "%s must be an element of a finite field"%j
2077
2078
F = j.parent()
2079
p = F.characteristic()
2080
d = F.degree()
2081
2082
if j.is_zero():
2083
return p==3 or p%3==2
2084
2085
if (j-1728).is_zero():
2086
return p==2 or p%4==3
2087
2088
# From now on we know that j != 0, 1728
2089
2090
if p in (2,3,5,7,11):
2091
return False # since j=0, 1728 are the only s.s. invariants
2092
2093
# supersingular j-invariants have degree at most 2:
2094
2095
jpol = j.minimal_polynomial()
2096
degj = jpol.degree()
2097
if degj > 2:
2098
return False
2099
2100
# if p occurs in the precomputed list, use that:
2101
2102
try:
2103
coeffs = supersingular_j_polynomials[p]
2104
return PolynomialRing(F,'x')(coeffs)(j).is_zero()
2105
except KeyError:
2106
pass
2107
2108
# Over GF(p), supersingular elliptic curves have cardinality
2109
# exactly p+1, so we check some random points in order to detect
2110
# non-supersingularity. Over GF(p^2) (for p at least 5) the
2111
# cardinality is either (p-1)^2 or (p+1)^2, and the group has
2112
# exponent p+1 or p-1, so we can do a similar random check: unless
2113
# (p+1)*P=0 for all the random points, or (p-1)*P=0 for all of
2114
# them, we can certainly return False.
2115
2116
# First we replace j by an element of GF(p) or GF(p^2) (since F
2117
# might be a proper extension of these):
2118
2119
if degj==1:
2120
j = -jpol(0) # = j, but in GF(p)
2121
elif d>2:
2122
F = GF(p^2,'a')
2123
j = jpol.roots(F,multiplicities=False)[0] # j, but in GF(p^2)
2124
2125
E = EllipticCurve(j=j)
2126
if degj==1:
2127
for i in range(10):
2128
P = E.random_element()
2129
if not ((p+1)*P).is_zero():
2130
return False
2131
else:
2132
n = None # will hold either p+1 or p-1 later
2133
for i in range(10):
2134
P = E.random_element()
2135
# avoid 2-torsion; we know that a1=a3=0 and #E>4!
2136
while P[2].is_zero() or P[1].is_zero():
2137
P = E.random_element()
2138
2139
if n is None: # not yet decided between p+1 and p-1
2140
pP = p*P
2141
if not pP[0]==P[0]: # i.e. pP is neither P nor -P
2142
return False
2143
if pP[1]==P[1]: # then p*P == P != -P
2144
n=p-1
2145
else: # then p*P == -P != P
2146
n=p+1
2147
else:
2148
if not (n*P).is_zero():
2149
return False
2150
2151
2152
# when proof is False we return True for any curve which passes
2153
# the probabilistic test:
2154
2155
if not proof:
2156
return True
2157
2158
# otherwise we check the trace of Frobenius (which could be
2159
# expensive since it involves counting the number of points on E):
2160
2161
return E.trace_of_frobenius() % p == 0
2162
2163
2164
2165