Book a Demo!
CoCalc Logo Icon
StoreFeaturesDocsShareSupportNewsAboutPoliciesSign UpSign In
sagemath
GitHub Repository: sagemath/sagelib
Path: blob/master/sage/schemes/elliptic_curves/ell_point.py
4128 views
1
# -*- coding: utf-8 -*-
2
r"""
3
Points on elliptic curves
4
5
The base class ``EllipticCurvePoint_field``, derived from
6
``AdditiveGroupElement``, provides support for points on elliptic
7
curves defined over general fields. The derived classes
8
``EllipticCurvePoint_number_field`` and
9
``EllipticCurvePoint_finite_field`` provide further support for point
10
on curves defined over number fields (including the rational field
11
`\QQ`) and over finite fields. Although there is no special
12
class for points over `\QQ`, there is currently greater
13
functionality implemented over `\QQ` than over other number
14
fields.
15
16
The class ``EllipticCurvePoint``, which is based on
17
``SchemeMorphism_point_projective_ring``, currently has little extra
18
functionality.
19
20
EXAMPLES:
21
22
An example over `\QQ`::
23
24
sage: E = EllipticCurve('389a1')
25
sage: P = E(-1,1); P
26
(-1 : 1 : 1)
27
sage: Q = E(0,-1); Q
28
(0 : -1 : 1)
29
sage: P+Q
30
(4 : 8 : 1)
31
sage: P-Q
32
(1 : 0 : 1)
33
sage: 3*P-5*Q
34
(328/361 : -2800/6859 : 1)
35
36
An example over a number field::
37
38
sage: K.<i> = QuadraticField(-1)
39
sage: E = EllipticCurve(K,[1,0,0,0,-1])
40
sage: P = E(0,i); P
41
(0 : i : 1)
42
sage: P.order()
43
+Infinity
44
sage: 101*P-100*P==P
45
True
46
47
An example over a finite field::
48
49
sage: K.<a> = GF(101^3)
50
sage: E = EllipticCurve(K,[1,0,0,0,-1])
51
sage: P = E(40*a^2 + 69*a + 84 , 58*a^2 + 73*a + 45)
52
sage: P.order()
53
1032210
54
sage: E.cardinality()
55
1032210
56
57
Arithmetic with a point over an extension of a finite field::
58
59
sage: k.<a> = GF(5^2)
60
sage: E = EllipticCurve(k,[1,0]); E
61
Elliptic Curve defined by y^2 = x^3 + x over Finite Field in a of size 5^2
62
sage: P = E([a,2*a+4])
63
sage: 5*P
64
(2*a + 3 : 2*a : 1)
65
sage: P*5
66
(2*a + 3 : 2*a : 1)
67
sage: P + P + P + P + P
68
(2*a + 3 : 2*a : 1)
69
70
::
71
72
sage: F = Zmod(3)
73
sage: E = EllipticCurve(F,[1,0]);
74
sage: P = E([2,1])
75
sage: import sys
76
sage: n = sys.maxint
77
sage: P*(n+1)-P*n == P
78
True
79
80
Arithmetic over `\ZZ/N\ZZ` with composite `N` is supported. When an
81
operation tries to invert a non-invertible element, a
82
ZeroDivisionError is raised and a factorization of the modulus appears
83
in the error message::
84
85
sage: N = 1715761513
86
sage: E = EllipticCurve(Integers(N),[3,-13])
87
sage: P = E(2,1)
88
sage: LCM([2..60])*P
89
Traceback (most recent call last):
90
...
91
ZeroDivisionError: Inverse of 1520944668 does not exist (characteristic = 1715761513 = 26927*63719)
92
93
94
AUTHORS:
95
96
- William Stein (2005) -- Initial version
97
98
- Robert Bradshaw et al....
99
100
- John Cremona (Feb 2008) -- Point counting and group structure for
101
non-prime fields, Frobenius endomorphism and order, elliptic logs
102
103
- John Cremona (Aug 2008) -- Introduced ``EllipticCurvePoint_number_field`` class
104
105
- Tobias Nagel, Michael Mardaus, John Cremona (Dec 2008) -- `p`-adic elliptic logarithm over `\QQ`
106
107
- David Hansen (Jan 2009) -- Added ``weil_pairing`` function to ``EllipticCurvePoint_finite_field`` class
108
109
- Mariah Lenox (March 2011) -- Added ``tate_pairing`` and ``ate_pairing``
110
functions to ``EllipticCurvePoint_finite_field`` class
111
112
"""
113
114
#*****************************************************************************
115
# Copyright (C) 2005 William Stein <[email protected]>
116
#
117
# Distributed under the terms of the GNU General Public License (GPL)
118
# as published by the Free Software Foundation; either version 2 of
119
# the License, or (at your option) any later version.
120
# http://www.gnu.org/licenses/
121
#*****************************************************************************
122
123
124
import math
125
126
from sage.structure.element import AdditiveGroupElement
127
import sage.plot.all as plot
128
129
from sage.rings.padics.factory import Qp
130
from sage.rings.padics.precision_error import PrecisionError
131
132
import sage.rings.all as rings
133
from sage.groups.all import AbelianGroup
134
import sage.groups.generic as generic
135
from sage.libs.pari.all import pari, PariError
136
from sage.structure.sequence import Sequence
137
138
from sage.schemes.plane_curves.projective_curve import Hasse_bounds
139
from sage.schemes.generic.morphism import (SchemeMorphism_point_projective_ring,
140
SchemeMorphism_point_abelian_variety_field,
141
is_SchemeMorphism)
142
143
from constructor import EllipticCurve
144
145
oo = rings.infinity # infinity
146
147
class EllipticCurvePoint(SchemeMorphism_point_projective_ring):
148
"""
149
A point on an elliptic curve.
150
"""
151
def __cmp__(self, other):
152
"""
153
Standard comparison function for points on elliptic curves, to
154
allow sorting and equality testing.
155
156
EXAMPLES:
157
sage: E=EllipticCurve(QQ,[1,1])
158
sage: P=E(0,1)
159
sage: P.order()
160
+Infinity
161
sage: Q=P+P
162
sage: P==Q
163
False
164
sage: Q+Q == 4*P
165
True
166
"""
167
assert isinstance(other, (int, long, rings.Integer)) and other == 0
168
if self.is_zero():
169
return 0
170
else:
171
return -1
172
173
174
class EllipticCurvePoint_field(SchemeMorphism_point_abelian_variety_field):
175
"""
176
A point on an elliptic curve over a field. The point has coordinates
177
in the base field.
178
179
EXAMPLES::
180
181
sage: E = EllipticCurve('37a')
182
sage: E([0,0])
183
(0 : 0 : 1)
184
sage: E(0,0) # brackets are optional
185
(0 : 0 : 1)
186
sage: E([GF(5)(0), 0]) # entries are coerced
187
(0 : 0 : 1)
188
189
sage: E(0.000, 0)
190
(0 : 0 : 1)
191
192
sage: E(1,0,0)
193
Traceback (most recent call last):
194
...
195
TypeError: Coordinates [1, 0, 0] do not define a point on
196
Elliptic Curve defined by y^2 + y = x^3 - x over Rational Field
197
198
::
199
200
sage: E = EllipticCurve([0,0,1,-1,0])
201
sage: S = E(QQ); S
202
Abelian group of points on Elliptic Curve defined by y^2 + y = x^3 - x over Rational Field
203
204
sage: K.<i>=NumberField(x^2+1)
205
sage: E=EllipticCurve(K,[0,1,0,-160,308])
206
sage: P=E(26,-120)
207
sage: Q=E(2+12*i,-36+48*i)
208
sage: P.order() == Q.order() == 4 # long time (3s)
209
True
210
sage: 2*P==2*Q
211
False
212
213
::
214
215
sage: K.<t>=FractionField(PolynomialRing(QQ,'t'))
216
sage: E=EllipticCurve([0,0,0,0,t^2])
217
sage: P=E(0,t)
218
sage: P,2*P,3*P
219
((0 : t : 1), (0 : -t : 1), (0 : 1 : 0))
220
221
222
TESTS::
223
224
sage: loads(S.dumps()) == S
225
True
226
sage: E = EllipticCurve('37a')
227
sage: P = E(0,0); P
228
(0 : 0 : 1)
229
sage: loads(P.dumps()) == P
230
True
231
sage: T = 100*P
232
sage: loads(T.dumps()) == T
233
True
234
235
Test pickling an elliptic curve that has known points on it::
236
237
sage: e = EllipticCurve([0, 0, 1, -1, 0]); g = e.gens(); loads(dumps(e)) == e
238
True
239
"""
240
def __init__(self, curve, v, check=True):
241
"""
242
Constructor for a point on an elliptic curve.
243
244
INPUT:
245
246
- curve -- an elliptic curve
247
- v -- data determining a point (another point, the integer
248
0, or a tuple of coordinates)
249
250
EXAMPLE::
251
252
sage: E = EllipticCurve('43a')
253
sage: P = E([2, -4, 2]); P
254
(1 : -2 : 1)
255
sage: P == E([1,-2])
256
True
257
sage: P = E(0); P
258
(0 : 1 : 0)
259
sage: P=E(2, -4, 2); P
260
(1 : -2 : 1)
261
"""
262
point_homset = curve.point_homset()
263
if is_SchemeMorphism(v) or isinstance(v,EllipticCurvePoint_field):
264
v = list(v)
265
elif v == 0:
266
# some of the code assumes that E(0) has integral entries
267
# irregardless of the base ring...
268
#R = self.base_ring()
269
#v = (R.zero(),R.one(),R.zero())
270
v = (0,1,0)
271
if check:
272
# mostly from SchemeMorphism_point_projective_field
273
d = point_homset.codomain().ambient_space().ngens()
274
if not isinstance(v,(list,tuple)):
275
raise TypeError, \
276
"Argument v (= %s) must be a scheme point, list, or tuple."%str(v)
277
if len(v) != d and len(v) != d-1:
278
raise TypeError, "v (=%s) must have %s components"%(v, d)
279
v = Sequence(v, point_homset.value_ring())
280
if len(v) == d-1: # very common special case
281
v.append(v.universe()(1))
282
283
n = len(v)
284
all_zero = True
285
for i in range(n):
286
c = v[n-1-i]
287
if c:
288
all_zero = False
289
if c == 1:
290
break
291
for j in range(n-i):
292
v[j] /= c
293
break
294
if all_zero:
295
raise ValueError, "%s does not define a valid point since all entries are 0"%repr(v)
296
297
x, y, z = v
298
if z == 0:
299
test = x
300
else:
301
a1, a2, a3, a4, a6 = curve.ainvs()
302
test = y**2 + (a1*x+a3)*y - (((x+a2)*x+a4)*x+a6)
303
if not test == 0:
304
raise TypeError, "Coordinates %s do not define a point on %s"%(list(v),curve)
305
306
SchemeMorphism_point_abelian_variety_field.__init__(self, point_homset, v, check=False)
307
#AdditiveGroupElement.__init__(self, point_homset)
308
309
def _repr_(self):
310
"""
311
Return a string representation of this point.
312
313
EXAMPLE::
314
315
sage: E = EllipticCurve('39a')
316
sage: P = E([-2, 1, 1])
317
sage: P._repr_()
318
'(-2 : 1 : 1)'
319
"""
320
return self.codomain().ambient_space()._repr_generic_point(self._coords)
321
322
def _latex_(self):
323
"""
324
Return a LaTeX representation of this point.
325
326
EXAMPLE::
327
328
sage: E = EllipticCurve('40a')
329
sage: P = E([3, 0])
330
sage: P._latex_()
331
'\\left(3 : 0 : 1\\right)'
332
"""
333
return self.codomain().ambient_space()._latex_generic_point(self._coords)
334
335
def __getitem__(self, n):
336
"""
337
Return the n'th coordinate of this point.
338
339
EXAMPLE::
340
341
sage: E = EllipticCurve('42a')
342
sage: P = E([-17, -51, 17])
343
sage: [P[i] for i in [2,1,0]]
344
[1, -3, -1]
345
"""
346
return self._coords[n]
347
348
def __iter__(self):
349
"""
350
Return the coordinates of this point as a list.
351
352
EXAMPLE::
353
354
sage: E = EllipticCurve('37a')
355
sage: list(E([0,0]))
356
[0, 0, 1]
357
"""
358
return iter(self._coords)
359
360
def __tuple__(self):
361
"""
362
Return the coordinates of this point as a tuple.
363
364
EXAMPLE::
365
366
sage: E = EllipticCurve('44a')
367
sage: P = E([1, -2, 1])
368
sage: P.__tuple__()
369
(1, -2, 1)
370
"""
371
return tuple(self._coords) # Warning: _coords is a list!
372
373
def __cmp__(self, other):
374
"""
375
Comparison function for points to allow sorting and equality testing.
376
377
EXAMPLES::
378
379
sage: E = EllipticCurve('45a')
380
sage: P = E([2, -1, 1])
381
sage: P == E(0)
382
False
383
sage: P+P == E(0)
384
True
385
"""
386
if not isinstance(other, EllipticCurvePoint_field):
387
try:
388
other = self.codomain().ambient_space()(other)
389
except TypeError:
390
return -1
391
return cmp(self._coords, other._coords)
392
393
def _pari_(self):
394
r"""
395
Converts this point to PARI format.
396
397
EXAMPLES::
398
399
sage: E = EllipticCurve([0,0,0,3,0])
400
sage: O = E(0)
401
sage: P = E.point([1,2])
402
sage: O._pari_()
403
[0]
404
sage: P._pari_()
405
[1, 2]
406
407
The following implicitly calls O._pari_() and P._pari_()::
408
409
sage: pari(E).elladd(O,P)
410
[1, 2]
411
412
TESTS::
413
414
Try the same over a finite field::
415
416
sage: E = EllipticCurve(GF(11), [0,0,0,3,0])
417
sage: O = E(0)
418
sage: P = E.point([1,2])
419
sage: O._pari_()
420
[0]
421
sage: P._pari_()
422
[Mod(1, 11), Mod(2, 11)]
423
sage: pari(E).elladd(O,P)
424
[Mod(1, 11), Mod(2, 11)]
425
"""
426
if self[2]:
427
return pari([self[0]/self[2], self[1]/self[2]])
428
else:
429
return pari([0])
430
431
def scheme(self):
432
"""
433
Return the scheme of this point, i.e., the curve it is on.
434
This is synonymous with curve() which is perhaps more
435
intuitive.
436
437
EXAMPLES::
438
439
sage: E=EllipticCurve(QQ,[1,1])
440
sage: P=E(0,1)
441
sage: P.scheme()
442
Elliptic Curve defined by y^2 = x^3 + x + 1 over Rational Field
443
sage: P.scheme() == P.curve()
444
True
445
sage: K.<a>=NumberField(x^2-3,'a')
446
sage: P=E.base_extend(K)(1,a)
447
sage: P.scheme()
448
Elliptic Curve defined by y^2 = x^3 + x + 1 over Number Field in a with defining polynomial x^2 - 3
449
"""
450
#The following text is just not true: it applies to the class
451
#EllipticCurvePoint, which appears to be never used, but does
452
#not apply to EllipticCurvePoint_field which is simply derived
453
#from AdditiveGroupElement.
454
#
455
#"Technically, points on curves in Sage are scheme maps from
456
# the domain Spec(F) where F is the base field of the curve to
457
# the codomain which is the curve. See also domain() and
458
# codomain()."
459
460
return self.codomain()
461
462
def domain(self):
463
"""
464
Return the domain of this point, which is `Spec(F)` where `F` is
465
the field of definition.
466
467
EXAMPLES::
468
469
sage: E=EllipticCurve(QQ,[1,1])
470
sage: P=E(0,1)
471
sage: P.domain()
472
Spectrum of Rational Field
473
sage: K.<a>=NumberField(x^2-3,'a')
474
sage: P=E.base_extend(K)(1,a)
475
sage: P.domain()
476
Spectrum of Number Field in a with defining polynomial x^2 - 3
477
"""
478
return self.parent().domain()
479
480
def codomain(self):
481
"""
482
Return the codomain of this point, which is the curve it is
483
on. Synonymous with curve() which is perhaps more intuitive.
484
485
EXAMPLES::
486
487
sage: E=EllipticCurve(QQ,[1,1])
488
sage: P=E(0,1)
489
sage: P.domain()
490
Spectrum of Rational Field
491
sage: K.<a>=NumberField(x^2-3,'a')
492
sage: P=E.base_extend(K)(1,a)
493
sage: P.codomain()
494
Elliptic Curve defined by y^2 = x^3 + x + 1 over Number Field in a with defining polynomial x^2 - 3
495
sage: P.codomain() == P.curve()
496
True
497
"""
498
return self.parent().codomain()
499
500
def order(self):
501
r"""
502
Return the order of this point on the elliptic curve.
503
504
If the point is zero, returns 1, otherwise raise a
505
NotImplementedError.
506
507
For curves over number fields and finite fields, see below.
508
509
.. note::
510
511
:meth:`additive_order` is a synonym for :meth:`order`
512
513
EXAMPLE::
514
515
sage: K.<t>=FractionField(PolynomialRing(QQ,'t'))
516
sage: E=EllipticCurve([0,0,0,-t^2,0])
517
sage: P=E(t,0)
518
sage: P.order()
519
Traceback (most recent call last):
520
...
521
NotImplementedError: Computation of order of a point not implemented over general fields.
522
sage: E(0).additive_order()
523
1
524
sage: E(0).order() == 1
525
True
526
527
"""
528
if self.is_zero():
529
return rings.Integer(1)
530
raise NotImplementedError, "Computation of order of a point not implemented over general fields."
531
532
additive_order = order
533
534
def curve(self):
535
"""
536
Return the curve that this point is on.
537
538
EXAMPLES::
539
540
sage: E = EllipticCurve('389a')
541
sage: P = E([-1,1])
542
sage: P.curve()
543
Elliptic Curve defined by y^2 + y = x^3 + x^2 - 2*x over Rational Field
544
"""
545
return self.scheme()
546
547
def __nonzero__(self):
548
"""
549
Return True if this is not the zero point on the curve.
550
551
EXAMPLES::
552
553
sage: E = EllipticCurve('37a')
554
sage: P = E(0); P
555
(0 : 1 : 0)
556
sage: P.is_zero()
557
True
558
sage: P = E.gens()[0]
559
sage: P.is_zero()
560
False
561
"""
562
return bool(self[2])
563
564
def has_finite_order(self):
565
"""
566
Return True if this point has finite additive order as an element
567
of the group of points on this curve.
568
569
For fields other than number fields and finite fields, this is
570
NotImplemented unless self.is_zero().
571
572
EXAMPLES::
573
574
sage: K.<t>=FractionField(PolynomialRing(QQ,'t'))
575
sage: E=EllipticCurve([0,0,0,-t^2,0])
576
sage: P = E(0)
577
sage: P.has_finite_order()
578
True
579
sage: P=E(t,0)
580
sage: P.has_finite_order()
581
Traceback (most recent call last):
582
...
583
NotImplementedError: Computation of order of a point not implemented over general fields.
584
sage: (2*P).is_zero()
585
True
586
"""
587
if self.is_zero(): return True
588
return self.order() != oo
589
590
is_finite_order = has_finite_order # for backward compatibility
591
592
def has_infinite_order(self):
593
"""
594
Return True if this point has infinite additive order as an element
595
of the group of points on this curve.
596
597
For fields other than number fields and finite fields, this is
598
NotImplemented unless self.is_zero().
599
600
EXAMPLES::
601
602
sage: K.<t>=FractionField(PolynomialRing(QQ,'t'))
603
sage: E=EllipticCurve([0,0,0,-t^2,0])
604
sage: P = E(0)
605
sage: P.has_infinite_order()
606
False
607
sage: P=E(t,0)
608
sage: P.has_infinite_order()
609
Traceback (most recent call last):
610
...
611
NotImplementedError: Computation of order of a point not implemented over general fields.
612
sage: (2*P).is_zero()
613
True
614
"""
615
if self.is_zero(): return False
616
return self.order() == oo
617
618
def plot(self, **args):
619
"""
620
Plot this point on an elliptic curve.
621
622
INPUT:
623
624
- ``**args`` -- all arguments get passed directly onto the point
625
plotting function.
626
627
EXAMPLES::
628
629
sage: E = EllipticCurve('389a')
630
sage: P = E([-1,1])
631
sage: P.plot(pointsize=30, rgbcolor=(1,0,0))
632
"""
633
if self.is_zero():
634
return plot.text("$\\infty$", (-3,3), **args)
635
636
else:
637
return plot.point((self[0], self[1]), **args)
638
639
def _add_(self, right):
640
"""
641
Add self to right.
642
643
EXAMPLES::
644
645
sage: E = EllipticCurve('389a')
646
sage: P = E([-1,1]); Q = E([0,0])
647
sage: P + Q
648
(1 : 0 : 1)
649
sage: P._add_(Q) == P + Q
650
True
651
652
Example to show that bug \#4820 is fixed::
653
654
sage: [type(c) for c in 2*EllipticCurve('37a1').gen(0)]
655
[<type 'sage.rings.rational.Rational'>,
656
<type 'sage.rings.rational.Rational'>,
657
<type 'sage.rings.rational.Rational'>]
658
"""
659
# Use Prop 7.1.7 of Cohen "A Course in Computational Algebraic Number Theory"
660
if self.is_zero():
661
return right
662
if right.is_zero():
663
return self
664
E = self.curve()
665
a1, a2, a3, a4, a6 = E.ainvs()
666
x1, y1 = self[0], self[1]
667
x2, y2 = right[0], right[1]
668
if x1 == x2 and y1 == -y2 - a1*x2 - a3:
669
return E(0) # point at infinity
670
671
if x1 == x2 and y1 == y2:
672
try:
673
m = (3*x1*x1 + 2*a2*x1 + a4 - a1*y1) / (2*y1 + a1*x1 + a3)
674
except ZeroDivisionError:
675
R = E.base_ring()
676
if R.is_finite():
677
N = R.characteristic()
678
from sage.rings.all import ZZ
679
N1 = N.gcd(ZZ(2*y1 + a1*x1 + a3))
680
N2 = N//N1
681
raise ZeroDivisionError, "Inverse of %s does not exist (characteristic = %s = %s*%s)"%(2*y1 + a1*x1 + a3, N,N1,N2)
682
else:
683
raise ZeroDivisionError, "Inverse of %s does not exist"%(2*y1 + a1*x1 + a3)
684
else:
685
try:
686
m = (y1-y2)/(x1-x2)
687
except ZeroDivisionError:
688
R = E.base_ring()
689
if R.is_finite():
690
N = R.characteristic()
691
from sage.rings.all import ZZ
692
N1 = N.gcd(ZZ(x1-x2))
693
N2 = N//N1
694
raise ZeroDivisionError, "Inverse of %s does not exist (characteristic = %s = %s*%s)"%(x1-x2, N,N1,N2)
695
else:
696
raise ZeroDivisionError, "Inverse of %s does not exist"%(x1-x2)
697
698
x3 = -x1 - x2 - a2 + m*(m+a1)
699
y3 = -y1 - a3 - a1*x3 + m*(x1-x3)
700
# See \#4820 for why we need to coerce 1 into the base ring here:
701
return E.point([x3, y3, E.base_ring()(1)], check=False)
702
703
def _sub_(self, right):
704
"""
705
Subtract right from self.
706
707
EXAMPLES::
708
709
sage: E = EllipticCurve('389a')
710
sage: P = E([-1,1]); Q = E([0,0])
711
sage: P - Q
712
(4 : 8 : 1)
713
sage: P - Q == P._sub_(Q)
714
True
715
sage: (P - Q) + Q
716
(-1 : 1 : 1)
717
sage: P
718
(-1 : 1 : 1)
719
"""
720
return self + (-right)
721
722
def __neg__(self):
723
"""
724
Return the additive inverse of this point.
725
726
EXAMPLES::
727
728
sage: E = EllipticCurve('389a')
729
sage: P = E([-1,1])
730
sage: Q = -P; Q
731
(-1 : -2 : 1)
732
sage: Q + P
733
(0 : 1 : 0)
734
735
Example to show that bug \#4820 is fixed::
736
737
sage: [type(c) for c in -EllipticCurve('37a1').gen(0)]
738
[<type 'sage.rings.rational.Rational'>,
739
<type 'sage.rings.rational.Rational'>,
740
<type 'sage.rings.rational.Rational'>]
741
"""
742
if self.is_zero():
743
return self
744
E, x, y = self.curve(), self[0], self[1]
745
# See \#4820 for why we need to coerce 1 into the base ring here:
746
return E.point([x, -y - E.a1()*x - E.a3(), E.base_ring()(1)], check=False)
747
748
def xy(self):
749
"""
750
Return the `x` and `y` coordinates of this point, as a 2-tuple.
751
If this is the point at infinity a ZeroDivisionError is raised.
752
753
EXAMPLES::
754
755
sage: E = EllipticCurve('389a')
756
sage: P = E([-1,1])
757
sage: P.xy()
758
(-1, 1)
759
sage: Q = E(0); Q
760
(0 : 1 : 0)
761
sage: Q.xy()
762
Traceback (most recent call last):
763
...
764
ZeroDivisionError: Rational division by zero
765
"""
766
if self[2] == 1:
767
return self[0], self[1]
768
else:
769
return self[0]/self[2], self[1]/self[2]
770
771
def is_divisible_by(self, m):
772
"""
773
Return True if there exists a point `Q` defined over the same
774
field as self such that `mQ` == self.
775
776
INPUT:
777
778
- ``m`` -- a positive integer.
779
780
OUTPUT:
781
782
(bool) -- True if there is a solution, else False.
783
784
.. warning::
785
786
This function usually triggers the computation of the
787
`m`-th division polynomial of the associated elliptic
788
curve, which will be expensive if `m` is large, though it
789
will be cached for subsequent calls with the same `m`.
790
791
EXAMPLES::
792
793
sage: E = EllipticCurve('389a')
794
sage: Q = 5*E(0,0); Q
795
(-2739/1444 : -77033/54872 : 1)
796
sage: Q.is_divisible_by(4)
797
False
798
sage: Q.is_divisible_by(5)
799
True
800
801
A finite field example::
802
803
sage: E = EllipticCurve(GF(101),[23,34])
804
sage: E.cardinality().factor()
805
2 * 53
806
sage: Set([T.order() for T in E.points()])
807
{1, 106, 2, 53}
808
sage: len([T for T in E.points() if T.is_divisible_by(2)])
809
53
810
sage: len([T for T in E.points() if T.is_divisible_by(3)])
811
106
812
813
TESTS:
814
815
This shows that the bug reported at #10076 is fixed::
816
817
sage: K = QuadraticField(8,'a')
818
sage: E = EllipticCurve([K(0),0,0,-1,0])
819
sage: P = E([-1,0])
820
sage: P.is_divisible_by(2)
821
False
822
sage: P.division_points(2)
823
[]
824
825
Note that it is not sufficient to test that
826
``self.division_points(m,poly_only=True)`` has roots::
827
828
sage: P.division_points(2, poly_only=True).roots()
829
[(1/2*a - 1, 1), (-1/2*a - 1, 1)]
830
831
sage: tor = E.torsion_points(); len(tor)
832
8
833
sage: [T.order() for T in tor]
834
[4, 2, 4, 2, 4, 2, 4, 1]
835
sage: all([T.is_divisible_by(3) for T in tor])
836
True
837
sage: Set([T for T in tor if T.is_divisible_by(2)])
838
{(0 : 1 : 0), (1 : 0 : 1)}
839
sage: Set([2*T for T in tor])
840
{(0 : 1 : 0), (1 : 0 : 1)}
841
842
843
"""
844
# Coerce the input m to an integer
845
m = rings.Integer(m)
846
847
# Check for trivial cases of m = 1, -1 and 0.
848
if m == 1 or m == -1:
849
return True
850
if m == 0:
851
return self == 0 # then m*self=self for all m!
852
m = m.abs()
853
854
# Now the following line would of course be correct, but we
855
# work harder to be more efficient:
856
# return len(self.division_points(m)) > 0
857
858
P = self
859
860
# If P has finite order n and gcd(m,n)=1 then the result is
861
# True. However, over general fields computing P.order() is
862
# not implemented.
863
864
try:
865
n = P.order()
866
if not n == oo:
867
if m.gcd(n)==1:
868
return True
869
except NotImplementedError:
870
pass
871
872
P_is_2_torsion = (P==-P)
873
g = P.division_points(m, poly_only=True)
874
875
if not P_is_2_torsion:
876
# In this case deg(g)=m^2, and each root in K lifts to two
877
# points Q,-Q both in E(K), of which exactly one is a
878
# solution. So we just check the existence of roots:
879
return len(g.roots())>0
880
881
# Now 2*P==0
882
883
if m%2==1:
884
return True # P itself is a solution when m is odd
885
886
# Now m is even and 2*P=0. Roots of g in K may or may not
887
# lift to solutions in E(K), so we fall back to the default.
888
# Note that division polynomials are cached so this is not
889
# inefficient:
890
891
return len(self.division_points(m)) > 0
892
893
def division_points(self, m, poly_only=False):
894
r"""
895
Return a list of all points `Q` such that `mQ=P` where `P` = self.
896
897
Only points on the elliptic curve containing self and defined
898
over the base field are included.
899
900
INPUT:
901
902
- ``m`` -- a positive integer
903
- ``poly_only`` -- bool (default: False); if True return polynomial whose roots give all possible `x`-coordinates of `m`-th roots of self.
904
905
OUTPUT:
906
907
(list) -- a (possibly empty) list of solutions `Q` to `mQ=P`, where `P` = self.
908
909
EXAMPLES:
910
911
We find the five 5-torsion points on an elliptic curve::
912
913
sage: E = EllipticCurve('11a'); E
914
Elliptic Curve defined by y^2 + y = x^3 - x^2 - 10*x - 20 over Rational Field
915
sage: P = E(0); P
916
(0 : 1 : 0)
917
sage: P.division_points(5)
918
[(0 : 1 : 0), (5 : -6 : 1), (5 : 5 : 1), (16 : -61 : 1), (16 : 60 : 1)]
919
920
Note above that 0 is included since [5]*0 = 0.
921
922
We create a curve of rank 1 with no torsion and do a consistency check::
923
924
sage: E = EllipticCurve('11a').quadratic_twist(-7)
925
sage: Q = E([44,-270])
926
sage: (4*Q).division_points(4)
927
[(44 : -270 : 1)]
928
929
We create a curve over a non-prime finite field with group of order `18`::
930
931
sage: k.<a> = GF(25)
932
sage: E = EllipticCurve(k, [1,2+a,3,4*a,2])
933
sage: P = E([3,3*a+4])
934
sage: factor(E.order())
935
2 * 3^2
936
sage: P.order()
937
9
938
939
We find the `1`-division points as a consistency check -- there
940
is just one, of course::
941
942
sage: P.division_points(1)
943
[(3 : 3*a + 4 : 1)]
944
945
The point `P` has order coprime to 2 but divisible by 3, so::
946
947
sage: P.division_points(2)
948
[(2*a + 1 : 3*a + 4 : 1), (3*a + 1 : a : 1)]
949
950
We check that each of the 2-division points works as claimed::
951
952
sage: [2*Q for Q in P.division_points(2)]
953
[(3 : 3*a + 4 : 1), (3 : 3*a + 4 : 1)]
954
955
Some other checks::
956
957
sage: P.division_points(3)
958
[]
959
sage: P.division_points(4)
960
[(0 : 3*a + 2 : 1), (1 : 0 : 1)]
961
sage: P.division_points(5)
962
[(1 : 1 : 1)]
963
964
An example over a number field (see trac #3383)::
965
966
sage: E = EllipticCurve('19a1')
967
sage: K.<t> = NumberField(x^9-3*x^8-4*x^7+16*x^6-3*x^5-21*x^4+5*x^3+7*x^2-7*x+1)
968
sage: EK = E.base_extend(K)
969
sage: E(0).division_points(3)
970
[(0 : 1 : 0), (5 : -10 : 1), (5 : 9 : 1)]
971
sage: EK(0).division_points(3)
972
[(0 : 1 : 0), (5 : 9 : 1), (5 : -10 : 1)]
973
sage: E(0).division_points(9)
974
[(0 : 1 : 0), (5 : -10 : 1), (5 : 9 : 1)]
975
sage: EK(0).division_points(9)
976
[(0 : 1 : 0), (5 : 9 : 1), (5 : -10 : 1), (-150/121*t^8 + 414/121*t^7 + 1481/242*t^6 - 2382/121*t^5 - 103/242*t^4 + 629/22*t^3 - 367/242*t^2 - 1307/121*t + 625/121 : 35/484*t^8 - 133/242*t^7 + 445/242*t^6 - 799/242*t^5 + 373/484*t^4 + 113/22*t^3 - 2355/484*t^2 - 753/242*t + 1165/484 : 1), (-150/121*t^8 + 414/121*t^7 + 1481/242*t^6 - 2382/121*t^5 - 103/242*t^4 + 629/22*t^3 - 367/242*t^2 - 1307/121*t + 625/121 : -35/484*t^8 + 133/242*t^7 - 445/242*t^6 + 799/242*t^5 - 373/484*t^4 - 113/22*t^3 + 2355/484*t^2 + 753/242*t - 1649/484 : 1), (-1383/484*t^8 + 970/121*t^7 + 3159/242*t^6 - 5211/121*t^5 + 37/484*t^4 + 654/11*t^3 - 909/484*t^2 - 4831/242*t + 6791/484 : 927/121*t^8 - 5209/242*t^7 - 8187/242*t^6 + 27975/242*t^5 - 1147/242*t^4 - 1729/11*t^3 + 1566/121*t^2 + 12873/242*t - 10871/242 : 1), (-1383/484*t^8 + 970/121*t^7 + 3159/242*t^6 - 5211/121*t^5 + 37/484*t^4 + 654/11*t^3 - 909/484*t^2 - 4831/242*t + 6791/484 : -927/121*t^8 + 5209/242*t^7 + 8187/242*t^6 - 27975/242*t^5 + 1147/242*t^4 + 1729/11*t^3 - 1566/121*t^2 - 12873/242*t + 10629/242 : 1), (-4793/484*t^8 + 6791/242*t^7 + 10727/242*t^6 - 18301/121*t^5 + 2347/484*t^4 + 2293/11*t^3 - 7311/484*t^2 - 17239/242*t + 26767/484 : 30847/484*t^8 - 21789/121*t^7 - 34605/121*t^6 + 117164/121*t^5 - 10633/484*t^4 - 29437/22*t^3 + 39725/484*t^2 + 55428/121*t - 176909/484 : 1), (-4793/484*t^8 + 6791/242*t^7 + 10727/242*t^6 - 18301/121*t^5 + 2347/484*t^4 + 2293/11*t^3 - 7311/484*t^2 - 17239/242*t + 26767/484 : -30847/484*t^8 + 21789/121*t^7 + 34605/121*t^6 - 117164/121*t^5 + 10633/484*t^4 + 29437/22*t^3 - 39725/484*t^2 - 55428/121*t + 176425/484 : 1)]
977
978
"""
979
# Coerce the input m to an integer
980
m = rings.Integer(m)
981
# Check for trivial cases of m = 1, -1 and 0.
982
if m == 1 or m == -1:
983
return [self]
984
if m == 0:
985
if self == 0: # then every point Q is a solution, but...
986
return [self]
987
else:
988
return []
989
990
# ans will contain the list of division points.
991
ans = []
992
993
# We compute a polynomial g whose roots give all possible x
994
# coordinates of the m-division points. The number of
995
# solutions (over the algebraic closure) is m^2, assuming that
996
# the characteristic does not divide m.
997
998
E = self.curve()
999
P = self
1000
nP = -P
1001
P_is_2_torsion = (P==nP)
1002
1003
# If self is the 0, then self is a solution, and the correct
1004
# poly is the m'th division polynomial
1005
if P == 0:
1006
ans.append(P)
1007
g = E.division_polynomial(m)
1008
else:
1009
# The poly g here is 0 at x(Q) iff x(m*Q) = x(P).
1010
g = E._multiple_x_numerator(m) - P[0]*E._multiple_x_denominator(m)
1011
1012
# When 2*P=0, then -Q is a solution iff Q is. For even m,
1013
# no 2-torsion point is a solution, so that g is the
1014
# square of a polynomial g1 of degree m^2/2, and each root
1015
# of g1 leads to a pair of solutions Q, -Q to m*Q=P. For
1016
# odd m, P itself is the only 2-torsion solution, so g has
1017
# the form (x-x(P))*g1(x)^2 where g1 has degree (m^2-1)/2
1018
# and each root of g1 leads to a pair Q, -Q.
1019
1020
if P_is_2_torsion:
1021
if m%2==0:
1022
# This computes g.sqrt() which is not implemented
1023
g = g.gcd(g.derivative())*g.leading_coefficient().sqrt()
1024
1025
# When 2*P!=0, then for each solution Q to m*Q=P, -Q is
1026
# not a solution (and points of order 2 are not
1027
# solutions). Hence the roots of g are distinct and each
1028
# gives rise to precisely one solution Q.
1029
1030
else:
1031
g0 = g.variables()[0] - P[0]
1032
g = g // g0
1033
g = g.gcd(g.derivative())*g.leading_coefficient().sqrt()
1034
g = g0*g
1035
1036
if poly_only:
1037
return g
1038
1039
for x in g.roots(multiplicities=False):
1040
if E.is_x_coord(x):
1041
# Make a point on the curve with this x coordinate.
1042
Q = E.lift_x(x)
1043
nQ = -Q
1044
mQ = m*Q
1045
# if P==-P then Q works iff -Q works, so we include
1046
# both unless they are equal:
1047
if P_is_2_torsion:
1048
if mQ == P:
1049
ans.append(Q)
1050
if nQ != Q:
1051
ans.append(nQ)
1052
else:
1053
# P is not 2-torsion so at most one of Q, -Q works
1054
# and we must try both:
1055
if mQ == P:
1056
ans.append(Q)
1057
elif mQ == nP:
1058
ans.append(nQ)
1059
1060
# Finally, sort and return
1061
ans.sort()
1062
return ans
1063
1064
def _divide_out(self,p):
1065
r"""
1066
Return `(Q,k)` where `p^kQ` == self and `Q` cannot be divided by `p`.
1067
1068
..WARNING:
1069
1070
It is up to the caller to make sure that this does not loop
1071
endlessly. It is used in
1072
``EllipticCurve_generic._p_primary_torsion_basis()``, when
1073
self will always have (finite) order which is a power of `p`,
1074
so that the order of `Q` increases by a factor of `p` at each
1075
stage.
1076
1077
Since it will clearly be in danger of looping when
1078
self.is_zero(), this case is caught, but otherwise caveat
1079
user.
1080
1081
EXAMPLES::
1082
1083
sage: E = EllipticCurve('37a1')
1084
sage: P = E([0, 0])
1085
sage: R = 12*P
1086
sage: R._divide_out(2)
1087
((-1 : -1 : 1), 2)
1088
sage: R._divide_out(3)
1089
((2 : -3 : 1), 1)
1090
sage: R._divide_out(5)
1091
((1357/841 : 28888/24389 : 1), 0)
1092
sage: R._divide_out(12)
1093
Traceback (most recent call last):
1094
...
1095
ValueError: p (=12) should be prime.
1096
"""
1097
p = rings.Integer(p)
1098
if not p.is_prime():
1099
raise ValueError, "p (=%s) should be prime."%p
1100
1101
if self.is_zero():
1102
raise ValueError, "self must not be 0."
1103
1104
k=0; Q=self
1105
pts = Q.division_points(p)
1106
while len(pts) > 0:
1107
Q = pts[0]
1108
k += 1
1109
pts = Q.division_points(p)
1110
return (Q,k)
1111
1112
def set_order(self, value):
1113
r"""
1114
Set the value of self._order to value.
1115
1116
Use this when you know a priori the order of this point to avoid a
1117
potentially expensive order calculation.
1118
1119
INPUT:
1120
1121
- ``value`` - positive Integer
1122
1123
OUTPUT:
1124
1125
None
1126
1127
EXAMPLES:
1128
1129
This example illustrates basic usage.
1130
1131
::
1132
1133
sage: E = EllipticCurve(GF(7), [0, 1]) # This curve has order 6
1134
sage: G = E(5, 0)
1135
sage: G.set_order(2)
1136
sage: 2*G
1137
(0 : 1 : 0)
1138
1139
We now give a more interesting case, the NIST-P521 curve. Its
1140
order is too big to calculate with SAGE, and takes a long time
1141
using other packages, so it is very useful here.
1142
1143
::
1144
1145
sage: p = 2^521 - 1
1146
sage: prev_proof_state = proof.arithmetic()
1147
sage: proof.arithmetic(False) # turn off primality checking
1148
sage: F = GF(p)
1149
sage: A = p - 3
1150
sage: B = 1093849038073734274511112390766805569936207598951683748994586394495953116150735016013708737573759623248592132296706313309438452531591012912142327488478985984
1151
sage: q = 6864797660130609714981900799081393217269435300143305409394463459185543183397655394245057746333217197532963996371363321113864768612440380340372808892707005449
1152
sage: E = EllipticCurve([F(A), F(B)])
1153
sage: G = E.random_point()
1154
sage: G.set_order(q)
1155
sage: G.order() * G # This takes practically no time.
1156
(0 : 1 : 0)
1157
sage: proof.arithmetic(prev_proof_state) # restore state
1158
1159
It is an error to pass a `value` equal to `0`::
1160
1161
sage: E = EllipticCurve(GF(7), [0, 1]) # This curve has order 6
1162
sage: G = E.random_point()
1163
sage: G.set_order(0)
1164
Traceback (most recent call last):
1165
...
1166
ValueError: Value 0 illegal for point order
1167
sage: G.set_order(1000)
1168
Traceback (most recent call last):
1169
...
1170
ValueError: Value 1000 illegal: outside max Hasse bound
1171
1172
It is also very likely an error to pass a value which is not the actual
1173
order of this point. How unlikely is determined by the factorization of
1174
the actual order, and the actual group structure::
1175
1176
sage: E = EllipticCurve(GF(7), [0, 1]) # This curve has order 6
1177
sage: G = E(5, 0) # G has order 2
1178
sage: G.set_order(11)
1179
Traceback (most recent call last):
1180
...
1181
ValueError: Value 11 illegal: 11 * (5 : 0 : 1) is not the identity
1182
1183
However, set_order can be fooled, though it's not likely in "real cases
1184
of interest". For instance, the order can be set to a multiple the
1185
actual order::
1186
1187
sage: E = EllipticCurve(GF(7), [0, 1]) # This curve has order 6
1188
sage: G = E(5, 0) # G has order 2
1189
sage: G.set_order(8)
1190
sage: G.order()
1191
8
1192
1193
NOTES:
1194
1195
The implementation is based of the fact that orders of elliptic curve
1196
points are cached in the (pseudo-private) _order slot.
1197
1198
AUTHORS:
1199
1200
- Mariah Lenox (2011-02-16)
1201
"""
1202
E = self.curve()
1203
q = E.base_ring().order()
1204
O = E(0)
1205
if value == 0:
1206
raise ValueError('Value 0 illegal for point order')
1207
if value == 1:
1208
if self == O:
1209
self._order = 1
1210
return
1211
else:
1212
raise ValueError('Value 1 illegal order for non-identity')
1213
low, hi = Hasse_bounds(q)
1214
if value > hi:
1215
raise ValueError('Value %s illegal: outside max Hasse bound'%value)
1216
if value * self != O:
1217
raise ValueError('Value %s illegal: %s * %s is not the identity'%(value, value, self))
1218
if (value - 1) * self == O:
1219
raise ValueError('Value %s illegal: %s * %s is the identity'%(value, value-1, self))
1220
self._order = value
1221
1222
############################## end ################################
1223
1224
def _line_(self,R,Q):
1225
r"""
1226
Computes the value at `Q` of a straight line through points self and `R`.
1227
1228
INPUT:
1229
1230
- ``R, Q`` -- points on self.curve() with ``Q`` nonzero.
1231
1232
OUTPUT:
1233
1234
An element of the base field self.curve().base_field().
1235
A ValueError is raised if ``Q`` is zero.
1236
1237
EXAMPLES::
1238
1239
sage: F.<a>=GF(2^5)
1240
sage: E=EllipticCurve(F,[0,0,1,1,1])
1241
sage: P = E(a^4 + 1, a^3)
1242
sage: Q = E(a^4, a^4 + a^3)
1243
sage: O = E(0)
1244
sage: P._line_(P,-2*P) == 0
1245
True
1246
sage: P._line_(Q,-(P+Q)) == 0
1247
True
1248
sage: O._line_(O,Q) == F(1)
1249
True
1250
sage: P._line_(O,Q) == a^4 - a^4 + 1
1251
True
1252
sage: P._line_(13*P,Q) == a^4
1253
True
1254
sage: P._line_(P,Q) == a^4 + a^3 + a^2 + 1
1255
True
1256
1257
See trac #7116::
1258
1259
sage: P._line_ (Q,O)
1260
Traceback (most recent call last):
1261
...
1262
ValueError: Q must be nonzero.
1263
1264
..NOTES:
1265
1266
This function is used in _miller_ algorithm.
1267
1268
AUTHOR:
1269
1270
- David Hansen (2009-01-25)
1271
"""
1272
if Q.is_zero():
1273
raise ValueError, "Q must be nonzero."
1274
1275
if self.is_zero() or R.is_zero():
1276
if self == R:
1277
return self.curve().base_field().one_element()
1278
if self.is_zero():
1279
return Q[0] - R[0]
1280
if R.is_zero():
1281
return Q[0] - self[0]
1282
elif self != R:
1283
if self[0] == R[0]:
1284
return Q[0] - self[0]
1285
else:
1286
l = (R[1] - self[1])/(R[0] - self[0])
1287
return Q[1] - self[1] - l * (Q[0] - self[0])
1288
else:
1289
a1, a2, a3, a4, a6 = self.curve().a_invariants()
1290
numerator = (3*self[0]**2 + 2*a2*self[0] + a4 - a1*self[1])
1291
denominator = (2*self[1] + a1*self[0] + a3)
1292
if denominator == 0:
1293
return Q[0] - self[0]
1294
else:
1295
l = numerator/denominator
1296
return Q[1] - self[1] - l * (Q[0] - self[0])
1297
1298
def _miller_(self,Q,n):
1299
r"""
1300
Return the value at `Q` of the rational function `f_{n,P}`, where the divisor of `f_{n,P}` is `n[P]-[nP]-(n-1)[O]`.
1301
1302
INPUT:
1303
1304
- ``Q`` -- a nonzero point on self.curve().
1305
1306
- ``n`` -- an nonzero integer. If `n<0` then return `Q`
1307
evaluated at `1/(v_{nP}*f_{n,P})` (used in the ate pairing).
1308
1309
OUTPUT:
1310
1311
An element in the base field self.curve().base_field()
1312
A ValueError is raised if `Q` is zero.
1313
1314
EXAMPLES::
1315
1316
sage: F.<a>=GF(2^5)
1317
sage: E=EllipticCurve(F,[0,0,1,1,1])
1318
sage: P = E(a^4 + 1, a^3)
1319
sage: Fx.<b>=GF(2^(4*5))
1320
sage: Ex=EllipticCurve(Fx,[0,0,1,1,1])
1321
sage: phi=Hom(F,Fx)(F.gen().minpoly().roots(Fx)[0][0])
1322
sage: Px=Ex(phi(P.xy()[0]),phi(P.xy()[1]))
1323
sage: Qx = Ex(b^19 + b^18 + b^16 + b^12 + b^10 + b^9 + b^8 + b^5 + b^3 + 1, b^18 + b^13 + b^10 + b^8 + b^5 + b^4 + b^3 + b)
1324
sage: Px._miller_(Qx,41) == b^17 + b^13 + b^12 + b^9 + b^8 + b^6 + b^4 + 1
1325
True
1326
sage: Qx._miller_(Px,41) == b^13 + b^10 + b^8 + b^7 + b^6 + b^5
1327
True
1328
sage: P._miller_(E(0),41)
1329
Traceback (most recent call last):
1330
...
1331
ValueError: Q must be nonzero.
1332
1333
An example of even order::
1334
1335
sage: F.<a> = GF(19^4)
1336
sage: E = EllipticCurve(F,[-1,0])
1337
sage: P = E(15*a^3 + 17*a^2 + 14*a + 13,16*a^3 + 7*a^2 + a + 18)
1338
sage: Q = E(10*a^3 + 16*a^2 + 4*a + 2, 6*a^3 + 4*a^2 + 3*a + 2)
1339
sage: x=P.weil_pairing(Q,360)
1340
sage: x^360 == F(1)
1341
True
1342
1343
You can use the _miller_ function on linearly dependent points, but with the risk of a dividing with zero::
1344
1345
sage: Px._miller_(2*Px,41)
1346
Traceback (most recent call last):
1347
...
1348
ZeroDivisionError: division by zero in finite field.
1349
1350
A small example of embedding degree 6::
1351
1352
sage: q = 401; F = GF(q); a = 146; b = 400; k = 6
1353
sage: E = EllipticCurve([F(a), F(b)])
1354
sage: R.<x> = F[]; K.<a> = GF(q^k, modulus=x^6 + 4*x^4 + 115*x^3 + 81*x^2 + 51*x + 3)
1355
sage: EK = E.base_extend(K)
1356
sage: P = E([F(338), F(227)])
1357
sage: Q_x = 333*a^5 + 391*a^4 + 160*a^3 + 335*a^2 + 71*a + 93
1358
sage: Q_y = 343*a^5 + 273*a^4 + 26*a^3 + 342*a^2 + 340*a + 210
1359
sage: Q = EK([Q_x, Q_y])
1360
sage: P._miller_(Q, 127)
1361
371*a^5 + 39*a^4 + 355*a^3 + 233*a^2 + 20*a + 275
1362
1363
A series of small examples and small torsions. We start with
1364
`n=1`, which is trivial: the function is the constant
1365
1::
1366
1367
sage: E = EllipticCurve([GF(7)(0), 2])
1368
sage: P = E(5, 1); Q = E(0, 3); I = E(0)
1369
sage: I._miller_(P, 1)
1370
1
1371
sage: I._miller_(Q, 1)
1372
1
1373
1374
A two-torsion example. In this case `f_{n,P}(Q) = x_Q - x_P`::
1375
1376
sage: E = EllipticCurve([GF(7)(-1), 0])
1377
sage: P = E(0,0); Q = E(1, 0);
1378
sage: P._miller_(P, 2)
1379
0
1380
sage: Q._miller_(Q, 2)
1381
0
1382
sage: P._miller_(Q, 2)
1383
1
1384
sage: Q._miller_(P, 2)
1385
6
1386
1387
A three-torsion example::
1388
1389
sage: E = EllipticCurve([GF(7)(0), 2])
1390
sage: P = E(5, 1); Q = E(0, 3);
1391
sage: P._miller_(Q, 3)
1392
4
1393
1394
A 4-torsion example::
1395
1396
sage: E = EllipticCurve([GF(7)(-1), 0])
1397
sage: P = E(5, 1); Q = E(4, 2)
1398
sage: P._miller_(Q, 4)
1399
3
1400
1401
A 5-torsion example::
1402
1403
sage: E = EllipticCurve([GF(7)(-1), 4])
1404
sage: P = E(4, 1); Q = E(6, 5)
1405
sage: P._miller_(Q, 5)
1406
1
1407
1408
A 6-torsion example::
1409
1410
sage: E = EllipticCurve([GF(7)(3), 1])
1411
sage: P = E(5, 1); Q = E(3, 3)
1412
sage: P._miller_(Q, 6)
1413
5
1414
1415
An example which is part of an ate pairing calculation. The trace of
1416
the curve is negative, so it should exercise the `n<0` case in the
1417
code::
1418
1419
sage: p = 2017; A = 1; B = 30; r = 29; t = -70; k = 7;
1420
sage: F = GF(p); R.<x> = F[]
1421
sage: E = EllipticCurve([F(A), F(B)]); P = E(369, 716)
1422
sage: K.<a> = GF(p^k, modulus=x^k+2); EK = E.base_extend(K)
1423
sage: Qx = 1226*a^6 + 1778*a^5 + 660*a^4 + 1791*a^3 + 1750*a^2 + 867*a + 770
1424
sage: Qy = 1764*a^6 + 198*a^5 + 1206*a^4 + 406*a^3 + 1200*a^2 + 273*a + 1712
1425
sage: Q = EK(Qx, Qy)
1426
sage: Q._miller_(P, t-1)
1427
1311*a^6 + 1362*a^5 + 1177*a^4 + 807*a^3 + 1331*a^2 + 1530*a + 1931
1428
1429
ALGORITHM:
1430
1431
Double-and-add.
1432
1433
REFERENCES:
1434
1435
- [Mil04] Victor S. Miller, "The Weil pairing, and its efficient calculation", J. Cryptol., 17(4):235-261, 2004
1436
1437
AUTHORS:
1438
1439
- David Hansen (2009-01-25)
1440
- Mariah Lenox (2011-03-07) -- Added more doctests and support for
1441
negative n.
1442
"""
1443
if Q.is_zero():
1444
raise ValueError, "Q must be nonzero."
1445
if n.is_zero():
1446
raise ValueError, "n must be nonzero."
1447
n_is_negative = False
1448
if n < 0:
1449
n = n.abs()
1450
n_is_negative = True
1451
1452
one = self.curve().base_field().one_element()
1453
t = one
1454
V = self
1455
S = 2*V
1456
nbin = n.bits()
1457
i = n.nbits() - 2
1458
while i > -1:
1459
S = 2*V
1460
ell = V._line_(V,Q)
1461
vee = S._line_(-S,Q)
1462
t = (t**2)*(ell/vee)
1463
V = S
1464
if nbin[i] == 1:
1465
S = V+self
1466
ell = V._line_(self,Q)
1467
vee = S._line_(-S,Q)
1468
t=t*(ell/vee)
1469
V = S
1470
i=i-1
1471
if n_is_negative:
1472
vee = V._line_(-V, Q)
1473
t = 1/(t*vee)
1474
return t
1475
1476
def weil_pairing(self, Q, n):
1477
r"""
1478
Compute the Weil pairing of self and `Q` using Miller's algorithm.
1479
1480
INPUT:
1481
1482
- ``Q`` -- a point on self.curve().
1483
1484
- ``n`` -- an integer `n` such that `nP = nQ = (0:1:0)` where `P` = self.
1485
1486
OUTPUT:
1487
1488
An `n`'th root of unity in the base field self.curve().base_field()
1489
1490
EXAMPLES::
1491
1492
sage: F.<a>=GF(2^5)
1493
sage: E=EllipticCurve(F,[0,0,1,1,1])
1494
sage: P = E(a^4 + 1, a^3)
1495
sage: Fx.<b>=GF(2^(4*5))
1496
sage: Ex=EllipticCurve(Fx,[0,0,1,1,1])
1497
sage: phi=Hom(F,Fx)(F.gen().minpoly().roots(Fx)[0][0])
1498
sage: Px=Ex(phi(P.xy()[0]),phi(P.xy()[1]))
1499
sage: O = Ex(0)
1500
sage: Qx = Ex(b^19 + b^18 + b^16 + b^12 + b^10 + b^9 + b^8 + b^5 + b^3 + 1, b^18 + b^13 + b^10 + b^8 + b^5 + b^4 + b^3 + b)
1501
sage: Px.weil_pairing(Qx,41) == b^19 + b^15 + b^9 + b^8 + b^6 + b^4 + b^3 + b^2 + 1
1502
True
1503
sage: Px.weil_pairing(17*Px,41) == Fx(1)
1504
True
1505
sage: Px.weil_pairing(O,41) == Fx(1)
1506
True
1507
1508
An error is raised if either point is not n-torsion::
1509
1510
sage: Px.weil_pairing(O,40)
1511
Traceback (most recent call last):
1512
...
1513
ValueError: points must both be n-torsion
1514
1515
A larger example (see trac \#4964)::
1516
1517
sage: P,Q = EllipticCurve(GF(19^4,'a'),[-1,0]).gens()
1518
sage: P.order(), Q.order()
1519
(360, 360)
1520
sage: z = P.weil_pairing(Q,360)
1521
sage: z.multiplicative_order()
1522
360
1523
1524
An example over a number field::
1525
1526
sage: P,Q = EllipticCurve('11a1').change_ring(CyclotomicField(5)).torsion_subgroup().gens() # long time (10s)
1527
sage: P,Q = (P.element(), Q.element()) # long time
1528
sage: (P.order(),Q.order()) # long time
1529
(5, 5)
1530
sage: P.weil_pairing(Q,5) # long time
1531
zeta5^2
1532
sage: Q.weil_pairing(P,5) # long time
1533
zeta5^3
1534
1535
ALGORITHM:
1536
1537
Implemented using Proposition 8 in [Mil04]. The value 1 is
1538
returned for linearly dependent input points. This condition
1539
is caught via a DivisionByZeroError, since the use of a
1540
discrete logarithm test for linear dependence, is much to slow
1541
for large `n`.
1542
1543
REFERENCES:
1544
1545
[Mil04] Victor S. Miller, "The Weil pairing, and its efficient
1546
calculation", J. Cryptol., 17(4):235-261, 2004
1547
1548
AUTHOR:
1549
1550
- David Hansen (2009-01-25)
1551
"""
1552
P = self
1553
E = P.curve()
1554
1555
if not Q.curve() is E:
1556
raise ValueError, "points must both be on the same curve"
1557
1558
# Test if P, Q are both in E[n]
1559
if not ((n*P).is_zero() and (n*Q).is_zero()):
1560
raise ValueError, "points must both be n-torsion"
1561
1562
one = E.base_field().one_element()
1563
1564
# Case where P = Q
1565
if P == Q:
1566
return one
1567
1568
# Case where P = O or Q = O
1569
if P.is_zero() or Q.is_zero():
1570
return one
1571
1572
# The non-trivial case P != Q
1573
1574
# Note (2010-12-29): The following code block should not be
1575
# used. It attempts to reduce the pairing calculation to order
1576
# d = gcd(|P|,|Q|) instead of order n, but this approach is
1577
# counterproductive, since calculating |P| and |Q| is much
1578
# slower than calculating the pairing [BGN05].
1579
#
1580
# [BGN05] D. Boneh, E. Goh, and K. Nissim, "Evaluating 2-DNF
1581
# Formulas on Ciphertexts", TCC 2005, LNCS 3378, pp. 325-341.
1582
#
1583
# Reduction to order d = gcd(|P|,|Q|); value is a d'th root of unity
1584
# try:
1585
# nP = P.order()
1586
# except AttributeError:
1587
# nP = generic.order_from_multiple(P,n,operation='+')
1588
# try:
1589
# nQ = Q.order()
1590
# except AttributeError:
1591
# nQ = generic.order_from_multiple(Q,n,operation='+')
1592
# d = arith.gcd(nP,nQ)
1593
# if d==1:
1594
# return one
1595
#
1596
# P = (nP//d)*P # order d
1597
# Q = (nQ//d)*Q # order d
1598
# n = d
1599
1600
try:
1601
return ((-1)**n.test_bit(0))*(P._miller_(Q,n)/Q._miller_(P,n))
1602
except ZeroDivisionError, detail:
1603
return one
1604
1605
def tate_pairing(self, Q, n, k, q=None):
1606
r"""
1607
Return Tate pairing of `n`-torsion point `P = self` and point `Q`.
1608
1609
The value returned is `f_{n,P}(Q)^e` where `f_{n,P}` is a function with
1610
divisor `n[P]-n[O].`. This is also known as the "modified Tate
1611
pairing". It is a well-defined bilinear map.
1612
1613
INPUT:
1614
1615
- ``P=self`` -- Elliptic curve point having order n
1616
1617
- ``Q`` -- Elliptic curve point on same curve as P (can be any order)
1618
1619
- ``n`` -- positive integer: order of P
1620
1621
- ``k`` -- positive integer: embedding degree
1622
1623
- ``q`` -- positive integer: size of base field (the "big"
1624
field is `GF(q^k)`. `q` needs to be set only if its value
1625
cannot be deduced.)
1626
1627
OUTPUT:
1628
1629
An `n`'th root of unity in the base field self.curve().base_field()
1630
1631
EXAMPLES:
1632
1633
A simple example, pairing a point with itself, and pairing a point with
1634
another rational point::
1635
1636
sage: p = 103; A = 1; B = 18; E = EllipticCurve(GF(p), [A, B])
1637
sage: P = E(33, 91); n = P.order(); n
1638
19
1639
sage: k = GF(n)(p).multiplicative_order(); k
1640
6
1641
sage: P.tate_pairing(P, n, k)
1642
1
1643
sage: Q = E(87, 51)
1644
sage: P.tate_pairing(Q, n, k)
1645
1
1646
sage: set_random_seed(35)
1647
sage: P.tate_pairing(P,n,k)
1648
1
1649
1650
We now let Q be a point on the same curve as above, but defined over
1651
the pairing extension field, and we also demonstrate the bilinearity of
1652
the pairing::
1653
1654
sage: K.<a> = GF(p^k)
1655
sage: EK = E.base_extend(K); P = EK(P)
1656
sage: Qx = 69*a^5 + 96*a^4 + 22*a^3 + 86*a^2 + 6*a + 35
1657
sage: Qy = 34*a^5 + 24*a^4 + 16*a^3 + 41*a^2 + 4*a + 40
1658
sage: Q = EK(Qx, Qy);
1659
sage: # multiply by cofactor so Q has order n:
1660
sage: h = 551269674; Q = h*Q
1661
sage: P = EK(P); P.tate_pairing(Q, n, k)
1662
24*a^5 + 34*a^4 + 3*a^3 + 69*a^2 + 86*a + 45
1663
sage: s = Integer(randrange(1,n))
1664
sage: ans1 = (s*P).tate_pairing(Q, n, k)
1665
sage: ans2 = P.tate_pairing(s*Q, n, k)
1666
sage: ans3 = P.tate_pairing(Q, n, k)^s
1667
sage: ans1 == ans2 == ans3
1668
True
1669
sage: (ans1 != 1) and (ans1^n == 1)
1670
True
1671
1672
Here is an example of using the Tate pairing to compute the Weil
1673
pairing (using the same data as above)::
1674
1675
sage: e = Integer((p^k-1)/n); e
1676
62844857712
1677
sage: P.weil_pairing(Q, n)^e
1678
94*a^5 + 99*a^4 + 29*a^3 + 45*a^2 + 57*a + 34
1679
sage: P.tate_pairing(Q, n, k) == P._miller_(Q, n)^e
1680
True
1681
sage: Q.tate_pairing(P, n, k) == Q._miller_(P, n)^e
1682
True
1683
sage: P.tate_pairing(Q, n, k)/Q.tate_pairing(P, n, k)
1684
94*a^5 + 99*a^4 + 29*a^3 + 45*a^2 + 57*a + 34
1685
1686
An example where we have to pass the base field size (and we again have
1687
agreement with the Weil pairing)::
1688
1689
sage: F.<a>=GF(2^5)
1690
sage: E=EllipticCurve(F,[0,0,1,1,1])
1691
sage: P = E(a^4 + 1, a^3)
1692
sage: Fx.<b>=GF(2^(4*5))
1693
sage: Ex=EllipticCurve(Fx,[0,0,1,1,1])
1694
sage: phi=Hom(F,Fx)(F.gen().minpoly().roots(Fx)[0][0])
1695
sage: Px=Ex(phi(P.xy()[0]),phi(P.xy()[1]))
1696
sage: Qx = Ex(b^19+b^18+b^16+b^12+b^10+b^9+b^8+b^5+b^3+1, b^18+b^13+b^10+b^8+b^5+b^4+b^3+b)
1697
sage: Px.tate_pairing(Qx, n=41, k=4)
1698
Traceback (most recent call last):
1699
...
1700
ValueError: Unexpected field degree: set keyword argument q equal to the size of the base field (big field is GF(q^4)).
1701
sage: num = Px.tate_pairing(Qx, n=41, k=4, q=32); num
1702
b^19 + b^14 + b^13 + b^12 + b^6 + b^4 + b^3
1703
sage: den = Qx.tate_pairing(Px, n=41, k=4, q=32); den
1704
b^19 + b^17 + b^16 + b^15 + b^14 + b^10 + b^6 + b^2 + 1
1705
sage: e = Integer((32^4-1)/41); e
1706
25575
1707
sage: Px.weil_pairing(Qx, 41)^e == num/den
1708
True
1709
1710
NOTES:
1711
1712
This function uses Miller's algorithm, followed by a naive
1713
exponentiation. It does not do anything fancy. In the case
1714
that there is an issue with `Q` being on one of the lines
1715
generated in the `r*P` calculation, `Q` is offset by a random
1716
point `R` and P.tate_pairing(Q+R,n,k)/P.tate_pairing(R,n,k)
1717
is returned.
1718
1719
AUTHORS:
1720
1721
- Mariah Lenox (2011-03-07)
1722
"""
1723
P = self
1724
E = P.curve()
1725
1726
if not Q.curve() is E:
1727
raise ValueError, "Points must both be on the same curve"
1728
1729
K = E.base_ring()
1730
d = K.degree()
1731
if q == None:
1732
if d == 1:
1733
q = K.order()
1734
elif d == k:
1735
q = K.base_ring().order()
1736
else:
1737
raise ValueError("Unexpected field degree: set keyword argument q equal to the size of the base field (big field is GF(q^%s))."%k)
1738
1739
if n*P != E(0):
1740
raise ValueError('This point is not of order n=%s' %n)
1741
1742
# In small cases, or in the case of pairing an element with
1743
# itself, Q could be on one of the lines in the Miller
1744
# algorithm. If this happens we try again, with an offset of a
1745
# random point.
1746
try:
1747
ret = self._miller_(Q, n)
1748
e = rings.Integer((q**k - 1)/n)
1749
ret = ret**e
1750
except (ZeroDivisionError, ValueError):
1751
R = E.random_point()
1752
ret = self.tate_pairing(Q + R, n, k)/self.tate_pairing(R, n, k)
1753
return ret
1754
1755
def ate_pairing(self, Q, n, k, t, q=None):
1756
r"""
1757
Return ate pairing of `n`-torsion points `P=self` and `Q`.
1758
1759
Also known as the `n`-th modified ate pairing. `P` is `GF(q)`-rational,
1760
and `Q` must be an element of `Ker(\pi-p)`, where `\pi` is the
1761
`q`-frobenius map (and hence `Q` is `GF(q^k)`-rational).
1762
1763
INPUT:
1764
1765
- ``P=self`` -- a point of order `n`, in `ker(\pi-1)`, where
1766
`\pi` is the `q`-Frobenius map (e.g., `P` is `q-rational`).
1767
1768
- ``Q`` -- a point of order `n` in `ker(\pi-q)`
1769
1770
- ``n`` -- the order of `P` and `Q`.
1771
1772
- ``k`` -- the embedding degree.
1773
1774
- ``t`` -- the trace of Frobenius of the curve over `GF(q)`.
1775
1776
- ``q`` -- (default:None) the size of base field (the "big"
1777
field is `GF(q^k)`). `q` needs to be set only if its value
1778
cannot be deduced.
1779
1780
OUTPUT:
1781
1782
FiniteFieldElement in `GF(q^k)` -- the ate pairing of `P` and `Q`.
1783
1784
EXAMPLES:
1785
1786
An example with embedding degree 6::
1787
1788
sage: p = 7549; A = 0; B = 1; n = 157; k = 6; t = 14
1789
sage: F = GF(p); E = EllipticCurve(F, [A, B])
1790
sage: R.<x> = F[]; K.<a> = GF(p^k, modulus=x^k+2)
1791
sage: EK = E.base_extend(K)
1792
sage: P = EK(3050, 5371); Q = EK(6908*a^4, 3231*a^3)
1793
sage: P.ate_pairing(Q, n, k, t)
1794
6708*a^5 + 4230*a^4 + 4350*a^3 + 2064*a^2 + 4022*a + 6733
1795
sage: s = Integer(randrange(1, n))
1796
sage: (s*P).ate_pairing(Q, n, k, t) == P.ate_pairing(s*Q, n, k, t)
1797
True
1798
sage: P.ate_pairing(s*Q, n, k, t) == P.ate_pairing(Q, n, k, t)^s
1799
True
1800
1801
Another example with embedding degree 7 and positive trace::
1802
1803
sage: p = 2213; A = 1; B = 49; n = 1093; k = 7; t = 28
1804
sage: F = GF(p); E = EllipticCurve(F, [A, B])
1805
sage: R.<x> = F[]; K.<a> = GF(p^k, modulus=x^k+2)
1806
sage: EK = E.base_extend(K)
1807
sage: P = EK(1583, 1734)
1808
sage: Qx = 1729*a^6+1767*a^5+245*a^4+980*a^3+1592*a^2+1883*a+722
1809
sage: Qy = 1299*a^6+1877*a^5+1030*a^4+1513*a^3+1457*a^2+309*a+1636
1810
sage: Q = EK(Qx, Qy)
1811
sage: P.ate_pairing(Q, n, k, t)
1812
1665*a^6 + 1538*a^5 + 1979*a^4 + 239*a^3 + 2134*a^2 + 2151*a + 654
1813
sage: s = Integer(randrange(1, n))
1814
sage: (s*P).ate_pairing(Q, n, k, t) == P.ate_pairing(s*Q, n, k, t)
1815
True
1816
sage: P.ate_pairing(s*Q, n, k, t) == P.ate_pairing(Q, n, k, t)^s
1817
True
1818
1819
Another example with embedding degree 7 and negative trace::
1820
1821
sage: p = 2017; A = 1; B = 30; n = 29; k = 7; t = -70
1822
sage: F = GF(p); E = EllipticCurve(F, [A, B])
1823
sage: R.<x> = F[]; K.<a> = GF(p^k, modulus=x^k+2)
1824
sage: EK = E.base_extend(K)
1825
sage: P = EK(369, 716)
1826
sage: Qx = 1226*a^6+1778*a^5+660*a^4+1791*a^3+1750*a^2+867*a+770
1827
sage: Qy = 1764*a^6+198*a^5+1206*a^4+406*a^3+1200*a^2+273*a+1712
1828
sage: Q = EK(Qx, Qy)
1829
sage: P.ate_pairing(Q, n, k, t)
1830
1794*a^6 + 1161*a^5 + 576*a^4 + 488*a^3 + 1950*a^2 + 1905*a + 1315
1831
sage: s = Integer(randrange(1, n))
1832
sage: (s*P).ate_pairing(Q, n, k, t) == P.ate_pairing(s*Q, n, k, t)
1833
True
1834
sage: P.ate_pairing(s*Q, n, k, t) == P.ate_pairing(Q, n, k, t)^s
1835
True
1836
1837
Using the same data, we show that the ate pairing is a power of the
1838
Tate pairing (see [HSV]_ end of section 3.1)::
1839
1840
sage: c = (k*p^(k-1)).mod(n); T = t - 1
1841
sage: N = gcd(T^k - 1, p^k - 1)
1842
sage: s = Integer(N/n)
1843
sage: L = Integer((T^k - 1)/N)
1844
sage: M = (L*s*c.inverse_mod(n)).mod(n)
1845
sage: P.ate_pairing(Q, n, k, t) == Q.tate_pairing(P, n, k)^M
1846
True
1847
1848
An example where we have to pass the base field size (and we again have
1849
agreement with the Tate pairing). Note that though `Px` is not
1850
`F`-rational, (it is the homomorphic image of an `F`-rational point) it
1851
is nonetheless in `ker(\pi-1)`, and so is a legitimate input::
1852
1853
sage: q = 2^5; F.<a>=GF(q)
1854
sage: n = 41; k = 4; t = -8
1855
sage: E=EllipticCurve(F,[0,0,1,1,1])
1856
sage: P = E(a^4 + 1, a^3)
1857
sage: Fx.<b>=GF(q^k)
1858
sage: Ex=EllipticCurve(Fx,[0,0,1,1,1])
1859
sage: phi=Hom(F,Fx)(F.gen().minpoly().roots(Fx)[0][0])
1860
sage: Px=Ex(phi(P.xy()[0]),phi(P.xy()[1]))
1861
sage: Qx = Ex(b^19+b^18+b^16+b^12+b^10+b^9+b^8+b^5+b^3+1, b^18+b^13+b^10+b^8+b^5+b^4+b^3+b)
1862
sage: Qx = Ex(Qx[0]^q, Qx[1]^q) - Qx # ensure Qx is in ker(pi - q)
1863
sage: Px.ate_pairing(Qx, n, k, t)
1864
Traceback (most recent call last):
1865
...
1866
ValueError: Unexpected field degree: set keyword argument q equal to the size of the base field (big field is GF(q^4)).
1867
sage: Px.ate_pairing(Qx, n, k, t, q)
1868
b^19 + b^18 + b^17 + b^16 + b^15 + b^14 + b^13 + b^12 + b^11 + b^9 + b^8 + b^5 + b^4 + b^2 + b + 1
1869
sage: s = Integer(randrange(1, n))
1870
sage: (s*Px).ate_pairing(Qx, n, k, t, q) == Px.ate_pairing(s*Qx, n, k, t, q)
1871
True
1872
sage: Px.ate_pairing(s*Qx, n, k, t, q) == Px.ate_pairing(Qx, n, k, t, q)^s
1873
True
1874
sage: c = (k*q^(k-1)).mod(n); T = t - 1
1875
sage: N = gcd(T^k - 1, q^k - 1)
1876
sage: s = Integer(N/n)
1877
sage: L = Integer((T^k - 1)/N)
1878
sage: M = (L*s*c.inverse_mod(n)).mod(n)
1879
sage: Px.ate_pairing(Qx, n, k, t, q) == Qx.tate_pairing(Px, n, k, q)^M
1880
True
1881
1882
It is an error if `Q` is not in the kernel of `\pi - p`, where `\pi` is
1883
the Frobenius automorphism::
1884
1885
sage: p = 29; A = 1; B = 0; n = 5; k = 2; t = 10
1886
sage: F = GF(p); R.<x> = F[]
1887
sage: E = EllipticCurve(F, [A, B]);
1888
sage: K.<a> = GF(p^k, modulus=x^k+2); EK = E.base_extend(K)
1889
sage: P = EK(13, 8); Q = EK(13, 21)
1890
sage: P.ate_pairing(Q, n, k, t)
1891
Traceback (most recent call last):
1892
...
1893
ValueError: Point (13 : 21 : 1) not in Ker(pi - q)
1894
1895
It is also an error if `P` is not in the kernel os `\pi - 1`::
1896
1897
sage: p = 29; A = 1; B = 0; n = 5; k = 2; t = 10
1898
sage: F = GF(p); R.<x> = F[]
1899
sage: E = EllipticCurve(F, [A, B]);
1900
sage: K.<a> = GF(p^k, modulus=x^k+2); EK = E.base_extend(K)
1901
sage: P = EK(14, 10*a); Q = EK(13, 21)
1902
sage: P.ate_pairing(Q, n, k, t)
1903
Traceback (most recent call last):
1904
...
1905
ValueError: This point (14 : 10*a : 1) is not in Ker(pi - 1)
1906
1907
NOTES:
1908
1909
First defined in the paper of [HSV]_, the ate pairing can be
1910
computationally effective in those cases when the trace of the curve
1911
over the base field is significantly smaller than the expected
1912
value. This implementation is simply Miller's algorithm followed by a
1913
naive exponentiation, and makes no claims towards efficiency.
1914
1915
REFERENCES:
1916
1917
.. [HSV] Hess, Smart, Vercauteren, "The Eta Pairing Revisited",
1918
IEEE Trans. Information Theory, 52(10): 4595-4602, 2006.
1919
1920
AUTHORS:
1921
1922
- Mariah Lenox (2011-03-08)
1923
"""
1924
P = self
1925
# check for same curve
1926
E = P.curve()
1927
O = E(0)
1928
if not Q.curve() is E:
1929
raise ValueError, "Points must both be on the same curve"
1930
1931
# set q to be the order of the base field
1932
K = E.base_ring()
1933
F = K.base_ring()
1934
if q == None:
1935
d = K.degree()
1936
if d == k:
1937
q = K.base_ring().order()
1938
else:
1939
raise ValueError("Unexpected field degree: set keyword argument q equal to the size of the base field (big field is GF(q^%s))." %k)
1940
1941
# check order of P
1942
if n*P != O:
1943
raise ValueError('This point %s is not of order n=%s' % (P,n))
1944
1945
# check for P in kernel pi - 1:
1946
piP = E(P[0]**q, P[1]**q)
1947
if piP - P != O:
1948
raise ValueError('This point %s is not in Ker(pi - 1)' % P)
1949
1950
# check for Q in kernel pi - q:
1951
piQ = E(Q[0]**q, Q[1]**q)
1952
if piQ - q*Q != O:
1953
raise ValueError('Point %s not in Ker(pi - q)' % Q)
1954
1955
T = t-1
1956
ret = Q._miller_(P, T)
1957
e = rings.Integer((q**k - 1)/n)
1958
ret = ret**e
1959
return ret
1960
1961
class EllipticCurvePoint_number_field(EllipticCurvePoint_field):
1962
"""
1963
A point on an elliptic curve over a number field.
1964
1965
Most of the functionality is derived from the parent class
1966
``EllipticCurvePoint_field``. In addition we have support for the
1967
order of a point, and heights (currently only implemented over
1968
`\QQ`).
1969
1970
EXAMPLES::
1971
1972
sage: E = EllipticCurve('37a')
1973
sage: E([0,0])
1974
(0 : 0 : 1)
1975
sage: E(0,0) # brackets are optional
1976
(0 : 0 : 1)
1977
sage: E([GF(5)(0), 0]) # entries are coerced
1978
(0 : 0 : 1)
1979
1980
sage: E(0.000, 0)
1981
(0 : 0 : 1)
1982
1983
sage: E(1,0,0)
1984
Traceback (most recent call last):
1985
...
1986
TypeError: Coordinates [1, 0, 0] do not define a point on
1987
Elliptic Curve defined by y^2 + y = x^3 - x over Rational Field
1988
1989
::
1990
1991
sage: E = EllipticCurve([0,0,1,-1,0])
1992
sage: S = E(QQ); S
1993
Abelian group of points on Elliptic Curve defined by y^2 + y = x^3 - x over Rational Field
1994
1995
TESTS::
1996
1997
sage: loads(S.dumps()) == S
1998
True
1999
sage: P = E(0,0); P
2000
(0 : 0 : 1)
2001
sage: loads(P.dumps()) == P
2002
True
2003
sage: T = 100*P
2004
sage: loads(T.dumps()) == T
2005
True
2006
2007
Test pickling an elliptic curve that has known points on it::
2008
2009
sage: e = EllipticCurve([0, 0, 1, -1, 0]); g = e.gens(); loads(dumps(e)) == e
2010
True
2011
"""
2012
2013
def order(self):
2014
r"""
2015
Return the order of this point on the elliptic curve.
2016
2017
If the point has infinite order, returns +Infinity. For
2018
curves defined over `\QQ`, we call PARI; over other
2019
number fields we implement the function here.
2020
2021
.. note::
2022
2023
:meth:`additive_order` is a synonym for :meth:`order`
2024
2025
EXAMPLES::
2026
2027
sage: E = EllipticCurve([0,0,1,-1,0])
2028
sage: P = E([0,0]); P
2029
(0 : 0 : 1)
2030
sage: P.order()
2031
+Infinity
2032
2033
::
2034
2035
sage: E = EllipticCurve([0,1])
2036
sage: P = E([-1,0])
2037
sage: P.order()
2038
2
2039
sage: P.additive_order()
2040
2
2041
2042
"""
2043
try:
2044
return self._order
2045
except AttributeError:
2046
pass
2047
2048
if self.is_zero():
2049
self._order = rings.Integer(1)
2050
return self._order
2051
2052
E = self.curve()
2053
2054
# Special code for curves over Q, calling PARI
2055
try:
2056
n = int(E.pari_curve().ellorder(self))
2057
if n == 0: n = oo
2058
self._order = n
2059
return n
2060
except PariError:
2061
pass
2062
2063
# Get the torsion order if known, else a bound on (multiple
2064
# of) the order. We do not compute the torsion if it is not
2065
# already known, since computing the bound is faster (and is
2066
# also cached).
2067
2068
try:
2069
N = E._torsion_order
2070
except AttributeError:
2071
N = E._torsion_bound()
2072
2073
# Now self is a torsion point iff it is killed by N:
2074
if not (N*self).is_zero():
2075
self._order = oo
2076
return self._order
2077
2078
# Finally we find the exact order using the generic code:
2079
self._order = generic.order_from_multiple(self,N,operation='+')
2080
return self._order
2081
2082
additive_order = order
2083
2084
def has_finite_order(self):
2085
"""
2086
Return True iff this point has finite order on the elliptic curve.
2087
2088
EXAMPLES::
2089
2090
sage: E = EllipticCurve([0,0,1,-1,0])
2091
sage: P = E([0,0]); P
2092
(0 : 0 : 1)
2093
sage: P.has_finite_order()
2094
False
2095
2096
::
2097
2098
sage: E = EllipticCurve([0,1])
2099
sage: P = E([-1,0])
2100
sage: P.has_finite_order()
2101
True
2102
"""
2103
if self.is_zero(): return True
2104
return self.order() != oo
2105
2106
def has_infinite_order(self):
2107
r"""
2108
Return True iff this point has infinite order on the elliptic curve.
2109
2110
EXAMPLES::
2111
2112
sage: E = EllipticCurve([0,0,1,-1,0])
2113
sage: P = E([0,0]); P
2114
(0 : 0 : 1)
2115
sage: P.has_infinite_order()
2116
True
2117
2118
::
2119
2120
sage: E = EllipticCurve([0,1])
2121
sage: P = E([-1,0])
2122
sage: P.has_infinite_order()
2123
False
2124
"""
2125
if self.is_zero(): return False
2126
return self.order() == oo
2127
2128
def is_on_identity_component(self, embedding=None):
2129
r"""
2130
Returns True iff this point is on the identity component of
2131
its curve with respect to a given (real or complex) embedding.
2132
2133
INPUT:
2134
2135
- ``self`` -- a point on a curve over any ordered field (e.g. `\QQ`)
2136
2137
- ``embedding`` -- an embedding from the base_field of the
2138
point's curve into `\RR` or `\CC`; if None (the
2139
default) it uses the first embedding of the base_field into
2140
`\RR` if any, else the first embedding into `\CC`.
2141
2142
OUTPUT:
2143
2144
(bool) -- True iff the point is on the identity component of
2145
the curve. (If the point is zero then the result is True.)
2146
2147
EXAMPLES:
2148
2149
For `K=\QQ` there is no need to specify an embedding::
2150
2151
sage: E=EllipticCurve('5077a1')
2152
sage: [E.lift_x(x).is_on_identity_component() for x in range(-3,5)]
2153
[False, False, False, False, False, True, True, True]
2154
2155
An example over a field with two real embeddings::
2156
2157
sage: L.<a> = QuadraticField(2)
2158
sage: E=EllipticCurve(L,[0,1,0,a,a])
2159
sage: P=E(-1,0)
2160
sage: [P.is_on_identity_component(e) for e in L.embeddings(RR)]
2161
[False, True]
2162
2163
We can check this as follows::
2164
2165
sage: [e(E.discriminant())>0 for e in L.embeddings(RR)]
2166
[True, False]
2167
sage: e = L.embeddings(RR)[0]
2168
sage: E1 = EllipticCurve(RR,[e(ai) for ai in E.ainvs()])
2169
sage: e1,e2,e3 = E1.two_division_polynomial().roots(RR,multiplicities=False)
2170
sage: e1 < e2 < e3 and e(P[0]) < e3
2171
True
2172
2173
"""
2174
if self.is_zero(): # trivial case
2175
return True
2176
2177
e = embedding
2178
# It is also trivially true if we have a complex embedding
2179
if not e is None:
2180
if not rings.is_RealField(e.codomain()):
2181
return True
2182
2183
# find a suitable embedding if none was supplied:
2184
E = self.curve()
2185
K = E.base_field()
2186
if e is None:
2187
try:
2188
e = K.embeddings(rings.RealField())[0]
2189
except IndexError:
2190
e = K.embeddings(rings.ComplexField())[0]
2191
2192
# If there is only one component, the result is True:
2193
if not rings.is_RealField(e.codomain()): # complex embedding
2194
return True
2195
if e(E.discriminant()) < 0: # only one component
2196
return True
2197
2198
# Now we have a real embedding and two components and have to work:
2199
gx = E.two_division_polynomial()
2200
gxd = gx.derivative()
2201
gxdd = gxd.derivative()
2202
return ( e(gxd(self[0])) > 0 and e(gxdd(self[0])) > 0)
2203
2204
def has_good_reduction(self, P=None):
2205
r"""
2206
Returns True iff this point has good reduction modulo a prime.
2207
2208
INPUT:
2209
2210
- ``P`` -- a prime of the base_field of the point's curve, or
2211
None (default)
2212
2213
OUTPUT:
2214
2215
(bool) If a prime `P` of the base field is specified, returns
2216
True iff the point has good reduction at `P`; otherwise,
2217
return true if the point has god reduction at all primes in
2218
the support of the discriminant of this model.
2219
2220
EXAMPLES::
2221
2222
sage: E = EllipticCurve('990e1')
2223
sage: P = E.gen(0); P
2224
(15 : 51 : 1)
2225
sage: [E.has_good_reduction(p) for p in [2,3,5,7]]
2226
[False, False, False, True]
2227
sage: [P.has_good_reduction(p) for p in [2,3,5,7]]
2228
[True, False, True, True]
2229
sage: [E.tamagawa_exponent(p) for p in [2,3,5,7]]
2230
[2, 2, 1, 1]
2231
sage: [(2*P).has_good_reduction(p) for p in [2,3,5,7]]
2232
[True, True, True, True]
2233
sage: P.has_good_reduction()
2234
False
2235
sage: (2*P).has_good_reduction()
2236
True
2237
sage: (3*P).has_good_reduction()
2238
False
2239
2240
::
2241
2242
sage: K.<i> = NumberField(x^2+1)
2243
sage: E = EllipticCurve(K,[0,1,0,-160,308])
2244
sage: P = E(26,-120)
2245
sage: E.discriminant().support()
2246
[Fractional ideal (i + 1),
2247
Fractional ideal (-i - 2),
2248
Fractional ideal (i - 2),
2249
Fractional ideal (3)]
2250
sage: [E.tamagawa_exponent(p) for p in E.discriminant().support()]
2251
[1, 4, 4, 4]
2252
sage: P.has_good_reduction()
2253
False
2254
sage: (2*P).has_good_reduction()
2255
False
2256
sage: (4*P).has_good_reduction()
2257
True
2258
2259
TESTS:
2260
2261
An example showing that #8498 is fixed::
2262
2263
sage: E = EllipticCurve('11a1')
2264
sage: K.<t> = NumberField(x^2+47)
2265
sage: EK = E.base_extend(K)
2266
sage: T = EK(5,5)
2267
sage: P = EK(-2, -1/2*t - 1/2)
2268
sage: p = K.ideal(11)
2269
sage: T.has_good_reduction(p)
2270
False
2271
sage: P.has_good_reduction(p)
2272
True
2273
"""
2274
if self.is_zero(): # trivial case
2275
return True
2276
2277
E = self.curve()
2278
if P is None:
2279
return all([self.has_good_reduction(P) for P in E.discriminant().support()])
2280
K = E.base_field()
2281
from sage.schemes.elliptic_curves.ell_local_data import check_prime
2282
P = check_prime(K,P)
2283
2284
# If the curve has good reduction at P, the result is True:
2285
t = E.local_data(P).bad_reduction_type()
2286
if t is None:
2287
return True
2288
2289
# Make sure the curve is integral and locally minimal at P:
2290
Emin = E.local_minimal_model(P)
2291
urst = E.isomorphism_to(Emin)
2292
Q = urst(self)
2293
2294
# Scale the homogeneous coordinates of the point to be primitive:
2295
xyz = list(Q)
2296
e = min([c.valuation(P) for c in xyz])
2297
if e !=0:
2298
if K is rings.QQ:
2299
pi = P
2300
else:
2301
pi = K.uniformizer(P)
2302
pie = pi**e
2303
xyz = [c/pie for c in xyz]
2304
2305
# Evaluate the partial derivatives at the point to see if they
2306
# are zero mod P
2307
2308
# See #8498: sometimes evaluating F's derivatives at xyz
2309
# returns a constant polynomial instead of a constant
2310
2311
F = Emin.defining_polynomial()
2312
for v in F.variables():
2313
c = (F.derivative(v))(xyz)
2314
try:
2315
val = c.valuation(P)
2316
except AttributeError:
2317
val = c.constant_coefficient().valuation(P)
2318
if val == 0:
2319
return True
2320
return False
2321
2322
def reduction(self,p):
2323
"""
2324
This finds the reduction of a point `P` on the elliptic curve modulo the prime `p`.
2325
2326
INPUT:
2327
2328
- ``self`` -- A point on an elliptic curve.
2329
2330
- ``p`` -- a prime number
2331
2332
OUTPUT:
2333
2334
The point reduced to be a point on the elliptic curve modulo `p`.
2335
2336
EXAMPLES::
2337
2338
sage: E = EllipticCurve([1,2,3,4,0])
2339
sage: P = E(0,0)
2340
sage: P.reduction(5)
2341
(0 : 0 : 1)
2342
sage: Q = E(98,931)
2343
sage: Q.reduction(5)
2344
(3 : 1 : 1)
2345
sage: Q.reduction(5).curve() == E.reduction(5)
2346
True
2347
2348
::
2349
2350
sage: F.<a> = NumberField(x^2+5)
2351
sage: E = EllipticCurve(F,[1,2,3,4,0])
2352
sage: Q = E(98,931)
2353
sage: Q.reduction(a)
2354
(3 : 1 : 1)
2355
sage: Q.reduction(11)
2356
(10 : 7 : 1)
2357
2358
::
2359
2360
sage: F.<a> = NumberField(x^3+x^2+1)
2361
sage: E = EllipticCurve(F,[a,2])
2362
sage: P = E(a,1)
2363
sage: P.reduction(F.ideal(5))
2364
(abar : 1 : 1)
2365
sage: P.reduction(F.ideal(a^2-4*a-2))
2366
(abar : 1 : 1)
2367
2368
"""
2369
P = self
2370
E = P.curve()
2371
return E.reduction(p)(P)
2372
2373
def height(self, precision=None):
2374
"""
2375
The Neron-Tate canonical height of the point.
2376
2377
INPUT:
2378
2379
- ``self`` -- a point on a curve over a number field
2380
2381
- ``precision`` -- (int or None (default)): the precision in
2382
bits of the result (default real precision if None)
2383
2384
OUTPUT:
2385
2386
The rational number 0, or a nonzero real field element.
2387
2388
The returned height is normalized to be independant of the base field.
2389
Fixing this, there are two normalizations used in the literature,
2390
one of which is double the other. We use the larger of the two,
2391
which is the one appropriate for the BSD conjecture. This is consistant
2392
with [Cre] and double that of [Sil].
2393
2394
REFERENCES:
2395
2396
- [Cre] John Cremona, Algorithms for modular elliptic curves,
2397
Cambridge University Press, 1997.
2398
2399
- [Sil] Silverman, Joseph H. The arithmetic of elliptic curves.
2400
Second edition. Graduate Texts in Mathematics, 106. Springer, 2009.
2401
2402
EXAMPLES::
2403
2404
sage: E = EllipticCurve('11a'); E
2405
Elliptic Curve defined by y^2 + y = x^3 - x^2 - 10*x - 20 over Rational Field
2406
sage: P = E([5,5]); P
2407
(5 : 5 : 1)
2408
sage: P.height()
2409
0
2410
sage: Q = 5*P
2411
sage: Q.height()
2412
0
2413
2414
::
2415
2416
sage: E = EllipticCurve('37a'); E
2417
Elliptic Curve defined by y^2 + y = x^3 - x over Rational Field
2418
sage: P = E([0,0])
2419
sage: P.height()
2420
0.0511114082399688
2421
sage: P.order()
2422
+Infinity
2423
sage: E.regulator()
2424
0.0511114082399688...
2425
2426
sage: def naive_height(P):
2427
... return log(RR(max(abs(P[0].numerator()), abs(P[0].denominator()))))
2428
sage: for n in [1..10]:
2429
... print naive_height(2^n*P)/4^n
2430
0.000000000000000
2431
0.0433216987849966
2432
0.0502949347635656
2433
0.0511006335618645
2434
0.0511007834799612
2435
0.0511013666152466
2436
0.0511034199907743
2437
0.0511106492906471
2438
0.0511114081541082
2439
0.0511114081541180
2440
2441
::
2442
2443
sage: E = EllipticCurve('4602a1'); E
2444
Elliptic Curve defined by y^2 + x*y = x^3 + x^2 - 37746035*x - 89296920339 over Rational Field
2445
sage: x = 77985922458974949246858229195945103471590
2446
sage: y = 19575260230015313702261379022151675961965157108920263594545223
2447
sage: d = 2254020761884782243
2448
sage: E([ x / d^2, y / d^3 ]).height()
2449
86.7406561381275
2450
2451
::
2452
2453
sage: E = EllipticCurve([17, -60, -120, 0, 0]); E
2454
Elliptic Curve defined by y^2 + 17*x*y - 120*y = x^3 - 60*x^2 over Rational Field
2455
sage: E([30, -90]).height()
2456
0
2457
2458
sage: E = EllipticCurve('389a1'); E
2459
Elliptic Curve defined by y^2 + y = x^3 + x^2 - 2*x over Rational Field
2460
sage: [P,Q] = [E(-1,1),E(0,-1)]
2461
sage: P.height(precision=100)
2462
0.68666708330558658572355210295
2463
sage: (3*Q).height(precision=100)/Q.height(precision=100)
2464
9.0000000000000000000000000000
2465
sage: _.parent()
2466
Real Field with 100 bits of precision
2467
2468
Canonical heights over number fields are implemented as well::
2469
2470
sage: R.<x> = QQ[]
2471
sage: K.<a> = NumberField(x^3-2)
2472
sage: E = EllipticCurve([a, 4]); E
2473
Elliptic Curve defined by y^2 = x^3 + a*x + 4 over Number Field in a with defining polynomial x^3 - 2
2474
sage: P = E((0,2))
2475
sage: P.height()
2476
0.810463096585925
2477
sage: P.height(precision=100)
2478
0.81046309658592536863991810577
2479
sage: P.height(precision=200)
2480
0.81046309658592536863991810576865158896130286417155832378086
2481
sage: (2*P).height() / P.height()
2482
4.00000000000000
2483
sage: (100*P).height() / P.height()
2484
10000.0000000000
2485
2486
Some consistency checks::
2487
2488
sage: E = EllipticCurve('5077a1')
2489
sage: P = E([-2,3,1])
2490
sage: P.height()
2491
1.36857250535393
2492
2493
sage: EK = E.change_ring(QuadraticField(-3,'a'))
2494
sage: PK = EK([-2,3,1])
2495
sage: PK.height()
2496
1.36857250535393
2497
2498
sage: K.<i> = NumberField(x^2+1)
2499
sage: E = EllipticCurve(K, [0,0,4,6*i,0])
2500
sage: Q = E.lift_x(-9/4); Q
2501
(-9/4 : 27/8*i - 4 : 1)
2502
sage: Q.height()
2503
2.69518560017909
2504
sage: (15*Q).height() / Q.height()
2505
225.000000000000
2506
2507
sage: E = EllipticCurve('37a')
2508
sage: P = E([0,-1])
2509
sage: P.height()
2510
0.0511114082399688
2511
sage: K.<a> = QuadraticField(-7)
2512
sage: ED = E.quadratic_twist(-7)
2513
sage: Q = E.isomorphism_to(ED.change_ring(K))(P); Q
2514
(0 : -7/2*a - 1/2 : 1)
2515
sage: Q.height()
2516
0.0511114082399688
2517
sage: Q.height(precision=100)
2518
0.051111408239968840235886099755
2519
2520
An example to show that the bug at \#5252 is fixed::
2521
2522
sage: E = EllipticCurve([1, -1, 1, -2063758701246626370773726978, 32838647793306133075103747085833809114881])
2523
sage: P = E([-30987785091199, 258909576181697016447])
2524
sage: P.height()
2525
25.8603170675462
2526
sage: P.height(precision=100)
2527
25.860317067546190743868840741
2528
sage: P.height(precision=250)
2529
25.860317067546190743868840740735110323098872903844416215577171041783572513
2530
sage: P.height(precision=500)
2531
25.8603170675461907438688407407351103230988729038444162155771710417835725129551130570889813281792157278507639909972112856019190236125362914195452321720
2532
2533
An example to show that the bug at \#8319 is fixed (correct height when the curve is not minimal)::
2534
2535
sage: E = EllipticCurve([-5580472329446114952805505804593498080000,-157339733785368110382973689903536054787700497223306368000000])
2536
sage: xP = 204885147732879546487576840131729064308289385547094673627174585676211859152978311600/23625501907057948132262217188983681204856907657753178415430361
2537
sage: P = E.lift_x(xP)
2538
sage: P.height()
2539
157.432598516754
2540
sage: Q = 2*P
2541
sage: Q.height() # long time (4s)
2542
629.730394067016
2543
sage: Q.height()-4*P.height() # long time
2544
0.000000000000000
2545
2546
"""
2547
if self.has_finite_order():
2548
return rings.QQ(0)
2549
2550
if precision is None:
2551
precision = rings.RealField().precision()
2552
2553
try:
2554
height = self.__height
2555
if height.prec() >= precision:
2556
return rings.RealField(precision)(height)
2557
except AttributeError:
2558
pass
2559
2560
if self.curve().base_ring() is rings.QQ:
2561
E = self.curve()
2562
Emin = E.minimal_model()
2563
iso = E.isomorphism_to(Emin)
2564
P = iso(self)
2565
h = Emin.pari_curve(prec=precision).ellheight(P, precision=precision)
2566
height = rings.RealField(precision)(h)
2567
else:
2568
from sage.misc.stopgap import stopgap
2569
stopgap("Computation of heights on elliptic curves over fields other than Q can return very imprecise results.", 12509)
2570
height = (self.nonarchimedian_local_height(prec=precision)
2571
+ self.archimedian_local_height(prec=precision))
2572
2573
self.__height = height
2574
return height
2575
2576
2577
def archimedian_local_height(self, v=None, prec=None):
2578
"""
2579
Computes the local height of self at the archimedian place `v`.
2580
2581
If `v` is None, returns the (weighted) sum of all archimedian
2582
contributions to the height.
2583
2584
The normalization is taken to be independant of the base field,
2585
but twice that in the paper. Note also that this local height depends
2586
on the model of the curve.
2587
2588
INPUT:
2589
2590
- ``v`` -- a real or complex embedding, or None
2591
2592
- ``prec`` -- the precision of the computation. If None, the precision
2593
is deduced from v.
2594
2595
ALGORITHM:
2596
2597
See section 4 of Silverman, J. Computing Heights on Elliptic Curves.
2598
Mathematics of Computation, Vol. 51, No. 183 (Jul., 1988), pp. 339-358
2599
2600
2601
EXAMPLES:
2602
2603
Examples 1, 2, and 3 from the above paper::
2604
2605
sage: K.<a> = QuadraticField(-2)
2606
sage: E = EllipticCurve(K, [0,-1,1,0,0]); E
2607
Elliptic Curve defined by y^2 + y = x^3 + (-1)*x^2 over Number Field in a with defining polynomial x^2 + 2
2608
sage: P = E.lift_x(2+a); P
2609
(a + 2 : 2*a + 1 : 1)
2610
sage: P.archimedian_local_height(K.places(prec=170)[0]) / 2
2611
0.45754773287523276736211210741423654346576029814695
2612
2613
sage: K.<i> = NumberField(x^2+1)
2614
sage: E = EllipticCurve(K, [0,0,4,6*i,0]); E
2615
Elliptic Curve defined by y^2 + 4*y = x^3 + 6*i*x over Number Field in i with defining polynomial x^2 + 1
2616
sage: P = E((0,0))
2617
sage: P.archimedian_local_height(K.places()[0]) / 2
2618
0.510184995162373
2619
2620
sage: Q = E.lift_x(-9/4); Q
2621
(-9/4 : 27/8*i - 4 : 1)
2622
sage: Q.archimedian_local_height(K.places()[0]) / 2
2623
0.654445619529600
2624
"""
2625
if self.has_finite_order():
2626
return QQ(0)
2627
2628
if v is None:
2629
K = self.curve().base_ring()
2630
def local_degree(v):
2631
"""
2632
Computes [ K_v : Q_v ]
2633
"""
2634
return 2 - int(v.im_gens()[0] in rings.RR)
2635
return sum(local_degree(v) * self.archimedian_local_height(v, prec) for v in K.places(prec=prec)) / K.degree()
2636
2637
if prec is None:
2638
prec = v.codomain().prec()
2639
E = self.curve()
2640
b2, b4, b6, b8 = [v(b) for b in E.b_invariants()]
2641
H = max(4, abs(b2), 2*abs(b4), 2*abs(b6), abs(b8))
2642
nterms = int(math.ceil(0.5*prec + 0.5 + 3*math.log(7+4*math.log(H)/3)/4 + math.log(max(1, ~abs(v(E.discriminant()))))/3))
2643
b2p = b2 - 12
2644
b4p = b4 - b2 + 6
2645
b6p = b6 - 2*b4 + b2 - 4
2646
b8p = b8 - 3*b6 + 3*b4 - b2 + 3
2647
2648
t = rings.polygen(v.codomain())
2649
fw = 4*t + b2 * t**2 + 2*b4 * t**3 + b6 * t**4
2650
fz = 1 - b4 * t**2 - 2*b6 * t**3 - b8 * t**4
2651
fwp = 4*t + b2p * t**2 + 2*b4p * t**3 + b6p * t**4
2652
fzp = 1 - b4p * t**2 - 2*b6p * t**3 - b8p * t**4
2653
2654
x = v(self[0])
2655
if abs(x) >= .5:
2656
t = 1/x
2657
beta = True
2658
else:
2659
t = 1/(x+1)
2660
beta = False
2661
lam = -abs(t).log()
2662
mu = 0
2663
four_to_n = rings.QQ(1)
2664
2665
for n in range(nterms):
2666
if beta:
2667
w = fw(t)
2668
z = fz(t)
2669
if abs(w) <= 2 * abs(z):
2670
mu += four_to_n * abs(z).log()
2671
t = w/z
2672
else:
2673
mu += four_to_n * abs(z+w).log()
2674
t = w/(z+w)
2675
beta = not beta
2676
else:
2677
w = fwp(t)
2678
z = fzp(t)
2679
if abs(w) <= 2 * abs(z):
2680
mu += four_to_n * abs(z).log()
2681
t = w/z
2682
else:
2683
mu += four_to_n * abs(z-w).log()
2684
t = w/(z-w)
2685
beta = not beta
2686
four_to_n >>= 2
2687
return lam + mu/4
2688
2689
def nonarchimedian_local_height(self, v=None, prec=None, weighted=False, is_minimal=None):
2690
"""
2691
Computes the local height of self at the non-archimedian place `v`.
2692
2693
If `v` is None, returns the (weighted) sum of all non-archimedian
2694
contributions to the height.
2695
2696
The normalization is taken to be independant of the base field,
2697
but twice that in the paper. Note also that this local height depends
2698
on the model of the curve.
2699
2700
INPUT:
2701
2702
- ``v`` -- a non-archimedean place of the base field of the curve,
2703
or None, in which case the total nonarchimedian contribution
2704
is returned
2705
2706
- ``prec`` -- working precision, or None in which case the height
2707
is returned symbolically
2708
2709
ALGORITHM:
2710
2711
See section 5 of Silverman, J. Computing Heights on Elliptic Curves.
2712
Mathematics of Computation, Vol. 51, No. 183 (Jul., 1988), pp. 339-358
2713
2714
EXAMPLES:
2715
2716
Examples 2 and 3 from the above paper::
2717
2718
sage: K.<i> = NumberField(x^2+1)
2719
sage: E = EllipticCurve(K, [0,0,4,6*i,0]); E
2720
Elliptic Curve defined by y^2 + 4*y = x^3 + 6*i*x over Number Field in i with defining polynomial x^2 + 1
2721
sage: P = E((0,0))
2722
sage: P.nonarchimedian_local_height(K.ideal(i+1))
2723
-1/2*log(2)
2724
sage: P.nonarchimedian_local_height(K.ideal(3))
2725
0
2726
sage: P.nonarchimedian_local_height(K.ideal(1-2*i))
2727
0
2728
2729
sage: Q = E.lift_x(-9/4); Q
2730
(-9/4 : 27/8*i - 4 : 1)
2731
sage: Q.nonarchimedian_local_height(K.ideal(1+i))
2732
2*log(2)
2733
sage: Q.nonarchimedian_local_height(K.ideal(3))
2734
0
2735
sage: Q.nonarchimedian_local_height(K.ideal(1-2*i))
2736
0
2737
sage: Q.nonarchimedian_local_height()
2738
1/2*log(16)
2739
2740
TESTS::
2741
2742
sage: Q.nonarchimedian_local_height(prec=100)
2743
1.3862943611198906188344642429
2744
sage: (3*Q).nonarchimedian_local_height()
2745
1/2*log(75923153929839865104)
2746
2747
sage: F.<a> = NumberField(x^4 + 2*x^3 + 19*x^2 + 18*x + 288)
2748
sage: F.ring_of_integers().gens()
2749
[1, 5/6*a^3 + 1/6*a, 1/6*a^3 + 1/6*a^2, a^3]
2750
sage: F.class_number()
2751
12
2752
sage: E = EllipticCurve('37a').change_ring(F)
2753
sage: P = E((-a^2/6 - a/6 - 1, a)); P
2754
(-1/6*a^2 - 1/6*a - 1 : a : 1)
2755
sage: P[0].is_integral()
2756
True
2757
sage: P.nonarchimedian_local_height()
2758
0
2759
"""
2760
if self.has_finite_order():
2761
return rings.QQ(0)
2762
2763
if prec:
2764
log = lambda x: rings.RealField(prec)(x).log()
2765
else:
2766
from sage.functions.log import log
2767
2768
if v is None:
2769
D = self.curve().discriminant()
2770
K = self.curve().base_ring()
2771
factorD = K.factor(D)
2772
if self[0] == 0:
2773
c = K.ideal(1)
2774
else:
2775
c = K.ideal(self[0]).denominator()
2776
# The last sum is for bad primes that divide c where the model is not minimal.
2777
return (log(c.norm())
2778
+ sum(self.nonarchimedian_local_height(v, prec, True, e < 12) for v,e in factorD if not v.divides(c))
2779
+ sum(self.nonarchimedian_local_height(v, prec, True) - c.valuation(v) * log(v.norm()) for v,e in factorD if e >= 12 and v.divides(c))
2780
) / K.degree()
2781
2782
if is_minimal:
2783
E = self.curve()
2784
P = self
2785
offset = 0
2786
else:
2787
E = self.curve().local_minimal_model(v)
2788
P = self.curve().isomorphism_to(E)(self)
2789
# Silverman's normalization is not invariant under change of model,
2790
# but it all cancels out in the global height.
2791
offset = (self.curve().discriminant()/E.discriminant()).valuation(v)
2792
2793
a1, a2, a3, a4, a6 = E.a_invariants()
2794
b2, b4, b6, b8 = E.b_invariants()
2795
c4 = E.c4()
2796
x,y = P.xy()
2797
D = E.discriminant()
2798
N = D.valuation(v)
2799
A = (3*x**2 + 2*a2*x + a4 - a1*y).valuation(v)
2800
B = (2*y+a1*x+a3).valuation(v)
2801
C = (3*x**4 + b2*x**3 + 3*b4*x**2 + 3*b6*x + b8).valuation(v)
2802
if A <= 0 or B <= 0:
2803
r = max(0, -x.valuation(v))
2804
elif c4.valuation(v) == 0:
2805
n = min(B, N/2)
2806
r = -n*(N-n)/N
2807
elif C > 3*B:
2808
r = -2*B/3
2809
else:
2810
r = -C/4
2811
r -= offset/6
2812
if not r:
2813
return rings.QQ(0)
2814
else:
2815
if E.base_ring() is rings.QQ:
2816
Nv = rings.ZZ(v)
2817
else:
2818
Nv = v.norm()
2819
if not weighted:
2820
r /= v.ramification_index() * v.residue_class_degree()
2821
return r * log(Nv)
2822
2823
def elliptic_logarithm(self, embedding=None, precision=100, algorithm='pari'):
2824
r"""
2825
Returns the elliptic logarithm of this elliptic curve point.
2826
2827
An embedding of the base field into `\RR` or `\CC` (with
2828
arbitrary precision) may be given; otherwise the first real
2829
embedding is used (with the specified precision) if any, else
2830
the first complex embedding.
2831
2832
INPUT:
2833
2834
- ``embedding``: an embedding of the base field into `\RR` or `\CC`
2835
2836
- ``precision``: a positive integer (default 100) setting the
2837
number of bits of precision for the computation
2838
2839
- ``algorithm``: either 'pari' (default for real embeddings)
2840
to use PARI's ``ellpointtoz{}``, or 'sage' for a native
2841
implementation. Ignored for complex embeddings.
2842
2843
ALGORITHM:
2844
2845
See [Co2] Cohen H., A Course in Computational Algebraic Number
2846
Theory GTM 138, Springer 1996 for the case of real embeddings,
2847
and Cremona, J.E. and Thongjunthug , T. 2010 for the complex
2848
case.
2849
2850
AUTHORS:
2851
2852
- Michael Mardaus (2008-07),
2853
- Tobias Nagel (2008-07) -- original version from [Co2].
2854
- John Cremona (2008-07) -- revision following eclib code.
2855
- John Cremona (2010-03) -- implementation for complex embeddings.
2856
2857
EXAMPLES::
2858
2859
sage: E = EllipticCurve('389a')
2860
sage: E.discriminant() > 0
2861
True
2862
sage: P = E([-1,1])
2863
sage: P.is_on_identity_component ()
2864
False
2865
sage: P.elliptic_logarithm (precision=96)
2866
0.4793482501902193161295330101 + 0.985868850775824102211203849...*I
2867
sage: Q=E([3,5])
2868
sage: Q.is_on_identity_component()
2869
True
2870
sage: Q.elliptic_logarithm (precision=96)
2871
1.931128271542559442488585220
2872
2873
An example with negative discriminant, and a torsion point::
2874
2875
sage: E = EllipticCurve('11a1')
2876
sage: E.discriminant() < 0
2877
True
2878
sage: P = E([16,-61])
2879
sage: P.elliptic_logarithm(precision=70)
2880
0.25384186085591068434
2881
sage: E.period_lattice().real_period(prec=70) / P.elliptic_logarithm(precision=70)
2882
5.0000000000000000000
2883
2884
A larger example. The default algorithm uses PARI and makes
2885
sure the result has the requested precision::
2886
2887
sage: E = EllipticCurve([1, 0, 1, -85357462, 303528987048]) #18074g1
2888
sage: P = E([4458713781401/835903744, -64466909836503771/24167649046528, 1])
2889
sage: P.elliptic_logarithm() # 100 bits
2890
0.27656204014107061464076203097
2891
2892
The native algorithm 'sage' used to have trouble with
2893
precision in this example, but no longer::
2894
2895
sage: P.elliptic_logarithm(algorithm='sage') # 100 bits
2896
0.27656204014107061464076203097
2897
2898
This shows that the bug reported at \#4901 has been fixed::
2899
2900
sage: E = EllipticCurve("4390c2")
2901
sage: P = E(683762969925/44944,-565388972095220019/9528128)
2902
sage: P.elliptic_logarithm()
2903
0.00025638725886520225353198932529
2904
sage: P.elliptic_logarithm(precision=64)
2905
0.000256387258865202254
2906
sage: P.elliptic_logarithm(precision=65)
2907
0.0002563872588652022535
2908
sage: P.elliptic_logarithm(precision=128)
2909
0.00025638725886520225353198932528666427412
2910
sage: P.elliptic_logarithm(precision=129)
2911
0.00025638725886520225353198932528666427412
2912
sage: P.elliptic_logarithm(precision=256)
2913
0.0002563872588652022535319893252866642741168388008346370015005142128009610936373
2914
sage: P.elliptic_logarithm(precision=257)
2915
0.00025638725886520225353198932528666427411683880083463700150051421280096109363730
2916
2917
Examples over number fields::
2918
2919
sage: K.<a> = NumberField(x^3-2)
2920
sage: embs = K.embeddings(CC)
2921
sage: E = EllipticCurve([0,1,0,a,a])
2922
sage: Ls = [E.period_lattice(e) for e in embs]
2923
sage: [L.real_flag for L in Ls]
2924
[0, 0, -1]
2925
sage: P = E(-1,0) # order 2
2926
sage: [L.elliptic_logarithm(P) for L in Ls]
2927
[-1.73964256006716 - 1.07861534489191*I, -0.363756518406398 - 1.50699412135253*I, 1.90726488608927]
2928
2929
sage: E = EllipticCurve([-a^2 - a - 1, a^2 + a])
2930
sage: Ls = [E.period_lattice(e) for e in embs]
2931
sage: pts = [E(2*a^2 - a - 1 , -2*a^2 - 2*a + 6 ), E(-2/3*a^2 - 1/3 , -4/3*a - 2/3 ), E(5/4*a^2 - 1/2*a , -a^2 - 1/4*a + 9/4 ), E(2*a^2 + 3*a + 4 , -7*a^2 - 10*a - 12 )]
2932
sage: [[L.elliptic_logarithm(P) for P in pts] for L in Ls]
2933
[[0.250819591818930 - 0.411963479992219*I, -0.290994550611374 - 1.37239400324105*I, -0.693473752205595 - 2.45028458830342*I, -0.151659609775291 - 1.48985406505459*I], [1.33444787667954 - 1.50889756650544*I, 0.792633734249234 - 0.548467043256610*I, 0.390154532655013 + 0.529423541805758*I, 0.931968675085317 - 0.431006981443071*I], [1.14758249500109 + 0.853389664016075*I, 2.59823462472518 + 0.853389664016075*I, 1.75372176444709, 0.303069634723001]]
2934
2935
::
2936
2937
sage: K.<i> = QuadraticField(-1)
2938
sage: E = EllipticCurve([0,0,0,9*i-10,21-i])
2939
sage: emb = K.embeddings(CC)[1]
2940
sage: L = E.period_lattice(emb)
2941
sage: P = E(2-i,4+2*i)
2942
sage: L.elliptic_logarithm(P,prec=100)
2943
0.70448375537782208460499649302 - 0.79246725643650979858266018068*I
2944
2945
"""
2946
from sage.rings.number_field.number_field import refine_embedding
2947
from sage.rings.all import RealField, ComplexField, is_RealField, QQ
2948
2949
# Check the trivial case:
2950
2951
C = ComplexField(precision)
2952
if self.is_zero():
2953
return C.zero()
2954
2955
# find a suitable embedding if none was supplied:
2956
2957
E = self.curve()
2958
K = E.base_field()
2959
rational = (K is QQ)
2960
emb = embedding
2961
2962
if emb is None:
2963
emb = K.embeddings(RealField(precision))
2964
if len(emb)>0:
2965
emb = emb[0]
2966
else:
2967
emb = K.embeddings(ComplexField(precision))[0]
2968
else:
2969
# Get the precision of the supplied embedding
2970
prec = emb.codomain().precision()
2971
# if the precision parameter is greater, refine the embedding:
2972
if precision > prec:
2973
emb = refine_embedding(emb,precision)
2974
2975
L = E.period_lattice(emb)
2976
2977
if algorithm == 'sage' or not is_RealField(emb.codomain):
2978
return L.elliptic_logarithm(self,precision)
2979
2980
if algorithm <> 'pari':
2981
raise ValueError, "algorithm must be either 'pari' or 'sage'"
2982
2983
# From now on emb() is a real embedding of K into
2984
# RealField(precision). We interface with the PARI library.
2985
2986
x, y = self.xy()
2987
if rational: # work with exact coordinates
2988
E_work = E
2989
pt_pari = pari([x,y])
2990
else: # use the embedding to get real coordinates
2991
ai = [emb(a) for a in E.a_invariants()]
2992
E_work = EllipticCurve(ai) # defined over RR
2993
pt_pari = pari([emb(x), emb(y)])
2994
working_prec = precision
2995
E_pari = E_work.pari_curve(prec=working_prec)
2996
log_pari = E_pari.ellpointtoz(pt_pari, precision=working_prec)
2997
2998
while prec_words_to_bits(log_pari.precision()) < precision:
2999
# result is not precise enough, re-compute with double
3000
# precision. if the base field is not QQ, this
3001
# requires modifying the precision of the embedding,
3002
# the curve, and the point
3003
working_prec = 2*working_prec
3004
if not rational:
3005
emb = refine_embedding(emb, working_prec)
3006
ai = [emb(a) for a in E.a_invariants()]
3007
E_work = EllipticCurve(ai) # defined over RR
3008
pt_pari = pari([emb(x), emb(y)])
3009
E_pari = E_work.pari_curve(prec=working_prec)
3010
log_pari = E_pari.ellpointtoz(pt_pari, precision=working_prec)
3011
3012
# normalization step
3013
r, i = C(log_pari)
3014
wR, wI = L.basis(prec=precision)
3015
k = (r/wR).floor()
3016
if k:
3017
r -= k*wR
3018
if self.is_on_identity_component(emb):
3019
return C(r)
3020
# Now there are two components and P is on the non-identity one
3021
return C(r)+C(wI/2)
3022
3023
def padic_elliptic_logarithm(self, p, absprec=20):
3024
r"""
3025
Computes the `p`-adic elliptic logarithm of this point.
3026
3027
INPUT:
3028
3029
``p`` - integer: a prime ``absprec`` - integer (default: 20):
3030
the initial `p`-adic absolute precision of the computation
3031
3032
OUTPUT:
3033
3034
The `p`-adic elliptic logarithm of self, with precision ``absprec``.
3035
3036
AUTHORS:
3037
3038
- Tobias Nagel
3039
- Michael Mardaus
3040
- John Cremona
3041
3042
ALGORITHM:
3043
3044
For points in the formal group (i.e. not integral at `p`) we
3045
take the ``log()`` function from the formal groups module and
3046
evaluate it at `-x/y`. Otherwise we first multiply the point
3047
to get into the formal group, and divide the result
3048
afterwards.
3049
3050
TODO:
3051
3052
See comments at trac \#4805. Currently the absolute precision
3053
of the result may be less than the given value of absprec, and
3054
error-handling is imperfect.
3055
3056
EXAMPLES::
3057
3058
sage: E = EllipticCurve([0,1,1,-2,0])
3059
sage: E(0).padic_elliptic_logarithm(3)
3060
0
3061
sage: P = E(0,0)
3062
sage: P.padic_elliptic_logarithm(3)
3063
2 + 2*3 + 3^3 + 2*3^7 + 3^8 + 3^9 + 3^11 + 3^15 + 2*3^17 + 3^18 + O(3^19)
3064
sage: P.padic_elliptic_logarithm(3).lift()
3065
660257522
3066
sage: P = E(-11/9,28/27)
3067
sage: [(2*P).padic_elliptic_logarithm(p)/P.padic_elliptic_logarithm(p) for p in prime_range(20)] # long time (3s)
3068
[2 + O(2^19), 2 + O(3^20), 2 + O(5^19), 2 + O(7^19), 2 + O(11^19), 2 + O(13^19), 2 + O(17^19), 2 + O(19^19)]
3069
sage: [(3*P).padic_elliptic_logarithm(p)/P.padic_elliptic_logarithm(p) for p in prime_range(12)] # long time (2s)
3070
[1 + 2 + O(2^19), 3 + 3^20 + O(3^21), 3 + O(5^19), 3 + O(7^19), 3 + O(11^19)]
3071
sage: [(5*P).padic_elliptic_logarithm(p)/P.padic_elliptic_logarithm(p) for p in prime_range(12)] # long time (2s)
3072
[1 + 2^2 + O(2^19), 2 + 3 + O(3^20), 5 + O(5^19), 5 + O(7^19), 5 + O(11^19)]
3073
3074
An example which arose during reviewing #4741::
3075
3076
sage: E = EllipticCurve('794a1')
3077
sage: P = E(-1,2)
3078
sage: P.padic_elliptic_logarithm(2) # default precision=20
3079
2^4 + 2^5 + 2^6 + 2^8 + 2^9 + 2^13 + 2^14 + 2^15 + O(2^16)
3080
sage: P.padic_elliptic_logarithm(2, absprec=30)
3081
2^4 + 2^5 + 2^6 + 2^8 + 2^9 + 2^13 + 2^14 + 2^15 + 2^22 + 2^23 + 2^24 + O(2^26)
3082
sage: P.padic_elliptic_logarithm(2, absprec=40)
3083
2^4 + 2^5 + 2^6 + 2^8 + 2^9 + 2^13 + 2^14 + 2^15 + 2^22 + 2^23 + 2^24 + 2^28 + 2^29 + 2^31 + 2^34 + O(2^35)
3084
"""
3085
if not p.is_prime():
3086
raise ValueError,'p must be prime'
3087
debug = False # True
3088
if debug:
3089
print "P=",self,"; p=",p," with precision ",precision
3090
E = self.curve()
3091
Q_p = Qp(p, absprec)
3092
if self.has_finite_order():
3093
return Q_p(0)
3094
while True:
3095
try:
3096
Ep = E.change_ring(Q_p)
3097
P = Ep(self)
3098
x,y = P.xy()
3099
break
3100
except (PrecisionError, ArithmeticError, ZeroDivisionError):
3101
absprec *=2
3102
Q_p = Qp(p, absprec)
3103
if debug:
3104
print "x,y=",(x,y)
3105
f = 1 # f will be such that f*P is in the formal group E^1(Q_p)
3106
if x.valuation() >=0: # P is not in E^1
3107
if not self.has_good_reduction(p): # P is not in E^0
3108
n = E.tamagawa_exponent(p) # n*P has good reduction at p
3109
if debug:
3110
print "Tamagawa exponent = =",n
3111
f = n
3112
P = n*P # lies in E^0
3113
if debug:
3114
print "P=",P
3115
try:
3116
x,y = P.xy()
3117
except ZeroDivisionError:
3118
raise ValueError, "Insufficient precision in p-adic_elliptic_logarithm()"
3119
if debug:
3120
print "x,y=",(x,y)
3121
if x.valuation() >=0: # P is still not in E^1
3122
t = E.local_data(p).bad_reduction_type()
3123
if t is None:
3124
m = E.reduction(p).abelian_group().exponent()
3125
else:
3126
m = p - t
3127
if debug:
3128
print "mod p exponent = =",m
3129
# now m*(n*P) reduces to the identity mod p, so is in E^1(Q_p)
3130
f *= m
3131
P = m*P # lies in E^1
3132
try:
3133
x,y = P.xy()
3134
except ZeroDivisionError:
3135
raise ValueError, "Insufficient precision in p-adic_elliptic_logarithm()"
3136
if debug:
3137
print "f=",f
3138
print "x,y=",(x,y)
3139
vx = x.valuation()
3140
vy = y.valuation()
3141
v = vx-vy
3142
if not ( v>0 and vx==-2*v and vy == -3*v):
3143
raise ValueError, "Insufficient precision in p-adic_elliptic_logarithm()"
3144
try:
3145
t = -x/y
3146
except (ZeroDivisionError, PrecisionError):
3147
raise ValueError, "Insufficient precision in p-adic_elliptic_logarithm()"
3148
if debug:
3149
print "t=",t,", with valuation ",v
3150
phi = Ep.formal().log(prec=1+absprec//v)
3151
return phi(t)/f
3152
3153
3154
class EllipticCurvePoint_finite_field(EllipticCurvePoint_field):
3155
r"""
3156
Class for elliptic curve points over finite fields.
3157
"""
3158
def _magma_init_(self, magma):
3159
"""
3160
Return a string representation of self that ``MAGMA`` can
3161
use for input.
3162
3163
EXAMPLE::
3164
3165
sage: E = EllipticCurve(GF(17), [1,-1])
3166
sage: P = E([13, 4])
3167
sage: P._magma_init_(magma) # optional - magma
3168
'EllipticCurve([_sage_ref...|GF(17)!0,GF(17)!0,GF(17)!0,GF(17)!1,GF(17)!16])![13,4]'
3169
"""
3170
E = self.curve()._magma_init_(magma)
3171
x,y = self.xy()
3172
return "%s![%s,%s]"%(E,x,y)
3173
3174
def discrete_log(self, Q, ord=None):
3175
r"""
3176
Returns discrete log of `Q` with respect to `P` =self.
3177
3178
INPUT:
3179
3180
- ``Q`` (point) -- another point on the same curve as self.
3181
3182
- ``ord`` (integer or None (default)) -- the order of self.
3183
3184
OUTPUT:
3185
3186
(integer) -- The discrete log of `Q` with respect to `P`, which is an
3187
integer `m` with `0\le m<o(P)` such that `mP=Q`, if one
3188
exists. A ValueError is raised if there is no solution.
3189
3190
.. note::
3191
3192
The order of self is computed if not supplied.
3193
3194
AUTHOR:
3195
3196
- John Cremona. Adapted to use generic functions 2008-04-05.
3197
3198
EXAMPLE::
3199
3200
sage: F = GF(3^6,'a')
3201
sage: a = F.gen()
3202
sage: E = EllipticCurve([0,1,1,a,a])
3203
sage: E.cardinality()
3204
762
3205
sage: A = E.abelian_group()
3206
sage: P = A.gen(0).element()
3207
sage: Q=400*P
3208
sage: P.discrete_log(Q)
3209
400
3210
"""
3211
if ord==None: ord=self.order()
3212
try:
3213
return generic.discrete_log(Q,self,ord,operation='+')
3214
except:
3215
raise ValueError, "ECDLog problem has no solution"
3216
3217
3218
def order(self):
3219
r"""
3220
Return the order of this point on the elliptic curve.
3221
3222
ALGORITHM:
3223
3224
Use generic functions from ``sage.groups.generic``. If the
3225
group order is known, use ``order_from_multiple()``, otherwise
3226
use ``order_from_bounds()`` with the Hasse bounds for the base
3227
field. In the latter case, we might find that we have a
3228
generator for the group, in which case it is cached.
3229
3230
We do not cause the group order to be calculated when not
3231
known, since this function is used in determining the group
3232
order via computation of several random points and their
3233
orders.
3234
3235
.. note::
3236
3237
:meth:`additive_order` is a synonym for :meth:`order`
3238
3239
AUTHOR:
3240
3241
- John Cremona, 2008-02-10, adapted 2008-04-05 to use generic functions.
3242
3243
EXAMPLES::
3244
3245
sage: k.<a> = GF(5^5)
3246
sage: E = EllipticCurve(k,[2,4]); E
3247
Elliptic Curve defined by y^2 = x^3 + 2*x + 4 over Finite Field in a of size 5^5
3248
sage: P = E(3*a^4 + 3*a , 2*a + 1 )
3249
sage: P.order()
3250
3227
3251
sage: Q = E(0,2)
3252
sage: Q.order()
3253
7
3254
sage: Q.additive_order()
3255
7
3256
3257
In the next example, the cardinality of E will be computed (using SEA) and cached::
3258
3259
sage: p=next_prime(2^150)
3260
sage: E=EllipticCurve(GF(p),[1,1])
3261
sage: P=E(831623307675610677632782670796608848711856078, 42295786042873366706573292533588638217232964)
3262
sage: P.order()
3263
1427247692705959881058262545272474300628281448
3264
sage: P.order()==E.cardinality()
3265
True
3266
3267
3268
"""
3269
try:
3270
return self._order
3271
except AttributeError:
3272
pass
3273
if self.is_zero():
3274
return rings.Integer(1)
3275
E = self.curve()
3276
K = E.base_ring()
3277
from sage.schemes.plane_curves.projective_curve import Hasse_bounds
3278
bounds = Hasse_bounds(K.order())
3279
3280
try:
3281
M = E._order
3282
try:
3283
plist = E._prime_factors_of_order
3284
except:
3285
plist = M.prime_divisors()
3286
E._prime_factors_of_order = plist
3287
N = generic.order_from_multiple(self,M,plist,operation='+')
3288
except:
3289
if K.is_prime_field():
3290
M = E.cardinality() # computed and cached
3291
plist = M.prime_divisors()
3292
E._prime_factors_of_order = plist
3293
N = generic.order_from_multiple(self,M,plist,operation='+')
3294
else:
3295
N = generic.order_from_bounds(self,bounds,operation='+')
3296
3297
if 2*N>bounds[1]: # then we have a generator, so cache this
3298
try:
3299
dummy = E._order
3300
except:
3301
E._order = N
3302
try:
3303
dummy = E.__abelian_group
3304
except:
3305
E.__abelian_group = AbelianGroup([N]), (self,)
3306
3307
self._order = N
3308
return self._order
3309
3310
additive_order = order
3311
3312
3313
3314