Book a Demo!
CoCalc Logo Icon
StoreFeaturesDocsShareSupportNewsAboutPoliciesSign UpSign In
sagemath
GitHub Repository: sagemath/sagesmc
Path: blob/master/src/sage/libs/mwrank/interface.py
8815 views
1
r"""
2
Sage interface to Cremona's eclib library (also known as mwrank)
3
4
This is the Sage interface to John Cremona's ``eclib`` C++ library for
5
arithmetic on elliptic curves. The classes defined in this module
6
give Sage interpreter-level access to some of the functionality of
7
``eclib``. For most purposes, it is not necessary to directly use these
8
classes. Instead, one can create an
9
:class:`EllipticCurve <sage.schemes.elliptic_curves.constructor.EllipticCurve>`
10
and call methods that are implemented using this module.
11
12
.. note::
13
14
This interface is a direct library-level interface to ``eclib``,
15
including the 2-descent program ``mwrank``.
16
"""
17
18
from sage.structure.sage_object import SageObject
19
from sage.rings.integer_ring import IntegerRing
20
21
# Need to permit tabs in order to doctest verbose output.
22
"""
23
SAGE_DOCTEST_ALLOW_TABS
24
"""
25
26
def get_precision():
27
r"""
28
Return the global NTL real number precision.
29
30
See also :meth:`set_precision`.
31
32
.. warning::
33
34
The internal precision is binary. This function multiplies the
35
binary precision by 0.3 (`=\log_2(10)` approximately) and
36
truncates.
37
38
OUTPUT:
39
40
(int) The current decimal precision.
41
42
EXAMPLES::
43
44
sage: mwrank_get_precision()
45
50
46
"""
47
# don't want to load mwrank every time Sage starts up, so we do
48
# the import here.
49
from sage.libs.mwrank.mwrank import get_precision
50
return get_precision()
51
52
def set_precision(n):
53
r"""
54
Set the global NTL real number precision. This has a massive
55
effect on the speed of mwrank calculations. The default (used if
56
this function is not called) is ``n=50``, but it might have to be
57
increased if a computation fails. See also :meth:`get_precision`.
58
59
INPUT:
60
61
- ``n`` (long) -- real precision used for floating point
62
computations in the library, in decimal digits.
63
64
.. warning::
65
66
This change is global and affects *all* future calls of eclib
67
functions by Sage.
68
69
EXAMPLES::
70
71
sage: mwrank_set_precision(20)
72
"""
73
# don't want to load mwrank every time Sage starts up, so we do
74
# the import here.
75
from sage.libs.mwrank.mwrank import set_precision
76
set_precision(n)
77
78
class mwrank_EllipticCurve(SageObject):
79
r"""
80
The :class:`mwrank_EllipticCurve` class represents an elliptic
81
curve using the ``Curvedata`` class from ``eclib``, called here an 'mwrank
82
elliptic curve'.
83
84
Create the mwrank elliptic curve with invariants
85
``ainvs``, which is a list of 5 or less *integers* `a_1`,
86
`a_2`, `a_3`, `a_4`, and `a_5`.
87
88
If strictly less than 5 invariants are given, then the *first*
89
ones are set to 0, so, e.g., ``[3,4]`` means `a_1=a_2=a_3=0` and
90
`a_4=3`, `a_5=4`.
91
92
INPUT:
93
94
- ``ainvs`` (list or tuple) -- a list of 5 or less integers, the
95
coefficients of a nonsingular Weierstrass equation.
96
97
- ``verbose`` (bool, default ``False``) -- verbosity flag. If ``True``,
98
then all Selmer group computations will be verbose.
99
100
EXAMPLES:
101
102
We create the elliptic curve `y^2 + y = x^3 + x^2 - 2x`::
103
104
sage: e = mwrank_EllipticCurve([0, 1, 1, -2, 0])
105
sage: e.ainvs()
106
[0, 1, 1, -2, 0]
107
108
This example illustrates that omitted `a`-invariants default to `0`::
109
110
sage: e = mwrank_EllipticCurve([3, -4])
111
sage: e
112
y^2 = x^3 + 3*x - 4
113
sage: e.ainvs()
114
[0, 0, 0, 3, -4]
115
116
The entries of the input list are coerced to :class:`int`.
117
If this is impossible, then an error is raised::
118
119
sage: e = mwrank_EllipticCurve([3, -4.8]); e
120
Traceback (most recent call last):
121
...
122
TypeError: ainvs must be a list or tuple of integers.
123
124
When you enter a singular model you get an exception::
125
126
sage: e = mwrank_EllipticCurve([0, 0])
127
Traceback (most recent call last):
128
...
129
ArithmeticError: Invariants (= [0, 0, 0, 0, 0]) do not describe an elliptic curve.
130
"""
131
132
def __init__(self, ainvs, verbose=False):
133
r"""
134
Create the mwrank elliptic curve with invariants
135
``ainvs``, which is a list of 5 or less *integers* `a_1`,
136
`a_2`, `a_3`, `a_4`, and `a_5`.
137
138
See the docstring of this class for full documentation.
139
140
EXAMPLES:
141
142
We create the elliptic curve `y^2 + y = x^3 + x^2 - 2x`::
143
144
sage: e = mwrank_EllipticCurve([0, 1, 1, -2, 0])
145
sage: e.ainvs()
146
[0, 1, 1, -2, 0]
147
"""
148
# import here to save time during startup (mwrank takes a while to init)
149
150
from sage.libs.mwrank.mwrank import _Curvedata
151
152
# if not isinstance(ainvs, list) and len(ainvs) <= 5:
153
if not isinstance(ainvs, (list,tuple)) or not len(ainvs) <= 5:
154
raise TypeError, "ainvs must be a list or tuple of length at most 5."
155
156
# Pad ainvs on the beginning by 0's, so e.g.
157
# [a4,a5] works.
158
ainvs = [0]*(5-len(ainvs)) + ainvs
159
160
# Convert each entry to an int
161
try:
162
a_int = [IntegerRing()(x) for x in ainvs]
163
except (TypeError, ValueError):
164
raise TypeError, "ainvs must be a list or tuple of integers."
165
self.__ainvs = a_int
166
self.__curve = _Curvedata(a_int[0], a_int[1], a_int[2],
167
a_int[3], a_int[4])
168
169
if verbose:
170
self.__verbose = True
171
else:
172
self.__verbose = False
173
174
# place holders
175
self.__saturate = -2 # not yet saturated
176
177
def __reduce__(self):
178
r"""
179
Standard Python function used in pickling.
180
181
EXAMPLES::
182
183
sage: E = mwrank_EllipticCurve([0,0,1,-7,6])
184
sage: E.__reduce__()
185
(<class 'sage.libs.mwrank.interface.mwrank_EllipticCurve'>, ([0, 0, 1, -7, 6], False))
186
187
188
"""
189
return mwrank_EllipticCurve, (self.__ainvs, self.__verbose)
190
191
192
def set_verbose(self, verbose):
193
"""
194
Set the verbosity of printing of output by the :meth:`two_descent()` and
195
other functions.
196
197
INPUT:
198
199
- ``verbose`` (int) -- if positive, print lots of output when
200
doing 2-descent.
201
202
EXAMPLES::
203
204
sage: E = mwrank_EllipticCurve([0, 0, 1, -1, 0])
205
sage: E.saturate() # no output
206
sage: E.gens()
207
[[0, -1, 1]]
208
209
sage: E = mwrank_EllipticCurve([0, 0, 1, -1, 0])
210
sage: E.set_verbose(1)
211
sage: E.saturate() # tol 1e-14
212
Basic pair: I=48, J=-432
213
disc=255744
214
2-adic index bound = 2
215
By Lemma 5.1(a), 2-adic index = 1
216
2-adic index = 1
217
One (I,J) pair
218
Looking for quartics with I = 48, J = -432
219
Looking for Type 2 quartics:
220
Trying positive a from 1 up to 1 (square a first...)
221
(1,0,-6,4,1) --trivial
222
Trying positive a from 1 up to 1 (...then non-square a)
223
Finished looking for Type 2 quartics.
224
Looking for Type 1 quartics:
225
Trying positive a from 1 up to 2 (square a first...)
226
(1,0,0,4,4) --nontrivial...(x:y:z) = (1 : 1 : 0)
227
Point = [0:0:1]
228
height = 0.0511114082399688402358
229
Rank of B=im(eps) increases to 1 (The previous point is on the egg)
230
Exiting search for Type 1 quartics after finding one which is globally soluble.
231
Mordell rank contribution from B=im(eps) = 1
232
Selmer rank contribution from B=im(eps) = 1
233
Sha rank contribution from B=im(eps) = 0
234
Mordell rank contribution from A=ker(eps) = 0
235
Selmer rank contribution from A=ker(eps) = 0
236
Sha rank contribution from A=ker(eps) = 0
237
Searching for points (bound = 8)...done:
238
found points which generate a subgroup of rank 1
239
and regulator 0.0511114082399688402358
240
Processing points found during 2-descent...done:
241
now regulator = 0.0511114082399688402358
242
Saturating (with bound = -1)...done:
243
points were already saturated.
244
"""
245
self.__verbose = verbose
246
247
248
def _curve_data(self):
249
r"""
250
Returns the underlying :class:`_Curvedata` class for this mwrank elliptic curve.
251
252
EXAMPLES::
253
254
sage: E = mwrank_EllipticCurve([0,0,1,-1,0])
255
sage: E._curve_data()
256
[0,0,1,-1,0]
257
b2 = 0 b4 = -2 b6 = 1 b8 = -1
258
c4 = 48 c6 = -216
259
disc = 37 (# real components = 2)
260
#torsion not yet computed
261
"""
262
return self.__curve
263
264
def ainvs(self):
265
r"""
266
Returns the `a`-invariants of this mwrank elliptic curve.
267
268
EXAMPLES::
269
270
sage: E = mwrank_EllipticCurve([0,0,1,-1,0])
271
sage: E.ainvs()
272
[0, 0, 1, -1, 0]
273
"""
274
return self.__ainvs
275
276
def isogeny_class(self, verbose=False):
277
r"""
278
Returns the isogeny class of this mwrank elliptic curve.
279
280
EXAMPLES::
281
282
sage: E = mwrank_EllipticCurve([0,-1,1,0,0])
283
sage: E.isogeny_class()
284
([[0, -1, 1, 0, 0], [0, -1, 1, -10, -20], [0, -1, 1, -7820, -263580]], [[0, 5, 0], [5, 0, 5], [0, 5, 0]])
285
"""
286
return self.__curve.isogeny_class(verbose)
287
288
def __repr__(self):
289
r"""
290
Returns the string representation of this mwrank elliptic curve.
291
292
EXAMPLES::
293
294
sage: E = mwrank_EllipticCurve([0,-1,1,0,0])
295
sage: E.__repr__()
296
'y^2+ y = x^3 - x^2 '
297
"""
298
# TODO: Is the use (or omission) of spaces here intentional?
299
a = self.ainvs()
300
s = "y^2"
301
if a[0] == -1:
302
s += "- x*y "
303
elif a[0] == 1:
304
s += "+ x*y "
305
elif a[0] != 0:
306
s += "+ %s*x*y "%a[0]
307
if a[2] == -1:
308
s += " - y"
309
elif a[2] == 1:
310
s += "+ y"
311
elif a[2] != 0:
312
s += "+ %s*y"%a[2]
313
s += " = x^3 "
314
if a[1] == -1:
315
s += "- x^2 "
316
elif a[1] == 1:
317
s += "+ x^2 "
318
elif a[1] != 0:
319
s += "+ %s*x^2 "%a[1]
320
if a[3] == -1:
321
s += "- x "
322
elif a[3] == 1:
323
s += "+ x "
324
elif a[3] != 0:
325
s += "+ %s*x "%a[3]
326
if a[4] == -1:
327
s += "-1"
328
elif a[4] == 1:
329
s += "+1"
330
elif a[4] != 0:
331
s += "+ %s"%a[4]
332
s = s.replace("+ -","- ")
333
return s
334
335
336
def two_descent(self,
337
verbose = True,
338
selmer_only = False,
339
first_limit = 20,
340
second_limit = 8,
341
n_aux = -1,
342
second_descent = True):
343
"""
344
Compute 2-descent data for this curve.
345
346
INPUT:
347
348
- ``verbose`` (bool, default ``True``) -- print what mwrank is doing.
349
350
- ``selmer_only`` (bool, default ``False``) -- ``selmer_only`` switch.
351
352
- ``first_limit`` (int, default 20) -- bound on `|x|+|z|` in
353
quartic point search.
354
355
- ``second_limit`` (int, default 8) -- bound on
356
`\log \max(|x|,|z|)`, i.e. logarithmic.
357
358
- ``n_aux`` (int, default -1) -- (only relevant for general
359
2-descent when 2-torsion trivial) number of primes used for
360
quartic search. ``n_aux=-1`` causes default (8) to be used.
361
Increase for curves of higher rank.
362
363
- ``second_descent`` (bool, default ``True``) -- (only relevant
364
for curves with 2-torsion, where mwrank uses descent via
365
2-isogeny) flag determining whether or not to do second
366
descent. *Default strongly recommended.*
367
368
369
OUTPUT:
370
371
Nothing -- nothing is returned.
372
373
TESTS (see #7992)::
374
375
sage: EllipticCurve([0, prod(prime_range(10))]).mwrank_curve().two_descent()
376
Basic pair: I=0, J=-5670
377
disc=-32148900
378
2-adic index bound = 2
379
2-adic index = 2
380
Two (I,J) pairs
381
Looking for quartics with I = 0, J = -5670
382
Looking for Type 3 quartics:
383
Trying positive a from 1 up to 5 (square a first...)
384
Trying positive a from 1 up to 5 (...then non-square a)
385
(2,0,-12,19,-6) --nontrivial...(x:y:z) = (2 : 4 : 1)
386
Point = [-2488:-4997:512]
387
height = 6.46767239...
388
Rank of B=im(eps) increases to 1
389
Trying negative a from -1 down to -3
390
Finished looking for Type 3 quartics.
391
Looking for quartics with I = 0, J = -362880
392
Looking for Type 3 quartics:
393
Trying positive a from 1 up to 20 (square a first...)
394
Trying positive a from 1 up to 20 (...then non-square a)
395
Trying negative a from -1 down to -13
396
Finished looking for Type 3 quartics.
397
Mordell rank contribution from B=im(eps) = 1
398
Selmer rank contribution from B=im(eps) = 1
399
Sha rank contribution from B=im(eps) = 0
400
Mordell rank contribution from A=ker(eps) = 0
401
Selmer rank contribution from A=ker(eps) = 0
402
Sha rank contribution from A=ker(eps) = 0
403
sage: EllipticCurve([0, prod(prime_range(100))]).mwrank_curve().two_descent()
404
Traceback (most recent call last):
405
...
406
RuntimeError: Aborted
407
408
409
"""
410
from sage.libs.mwrank.mwrank import _two_descent # import here to save time
411
first_limit = int(first_limit)
412
second_limit = int(second_limit)
413
n_aux = int(n_aux)
414
second_descent = int(second_descent) # convert from bool to (int) 0 or 1
415
# TODO: Don't allow limits above some value...???
416
# (since otherwise mwrank just sets limit tiny)
417
self.__descent = _two_descent()
418
self.__descent.do_descent(self.__curve,
419
verbose,
420
selmer_only,
421
first_limit,
422
second_limit,
423
n_aux,
424
second_descent)
425
if not self.__two_descent_data().ok():
426
raise RuntimeError, "A 2-descent did not complete successfully."
427
428
def __two_descent_data(self):
429
r"""
430
Returns the 2-descent data for this elliptic curve.
431
432
EXAMPLES::
433
434
sage: E = mwrank_EllipticCurve([0,-1,1,0,0])
435
sage: E._mwrank_EllipticCurve__two_descent_data()
436
<sage.libs.mwrank.mwrank._two_descent object at ...>
437
"""
438
try:
439
return self.__descent
440
except AttributeError:
441
self.two_descent(self.__verbose)
442
return self.__descent
443
444
def conductor(self):
445
"""
446
Return the conductor of this curve, computed using Cremona's
447
implementation of Tate's algorithm.
448
449
.. note::
450
451
This is independent of PARI's.
452
453
EXAMPLES::
454
455
sage: E = mwrank_EllipticCurve([1, 1, 0, -6958, -224588])
456
sage: E.conductor()
457
2310
458
"""
459
return self.__curve.conductor()
460
461
def rank(self):
462
"""
463
Returns the rank of this curve, computed using :meth:`two_descent()`.
464
465
In general this may only be a lower bound for the rank; an
466
upper bound may be obtained using the function :meth:`rank_bound()`.
467
To test whether the value has been proved to be correct, use
468
the method :meth:`certain()`.
469
470
EXAMPLES::
471
472
sage: E = mwrank_EllipticCurve([0, -1, 0, -900, -10098])
473
sage: E.rank()
474
0
475
sage: E.certain()
476
True
477
478
::
479
480
sage: E = mwrank_EllipticCurve([0, -1, 1, -929, -10595])
481
sage: E.rank()
482
0
483
sage: E.certain()
484
False
485
486
"""
487
return self.__two_descent_data().getrank()
488
489
def rank_bound(self):
490
"""
491
Returns an upper bound for the rank of this curve, computed
492
using :meth:`two_descent()`.
493
494
If the curve has no 2-torsion, this is equal to the 2-Selmer
495
rank. If the curve has 2-torsion, the upper bound may be
496
smaller than the bound obtained from the 2-Selmer rank minus
497
the 2-rank of the torsion, since more information is gained
498
from the 2-isogenous curve or curves.
499
500
EXAMPLES:
501
502
The following is the curve 960D1, which has rank 0,
503
but Sha of order 4::
504
505
sage: E = mwrank_EllipticCurve([0, -1, 0, -900, -10098])
506
sage: E.rank_bound()
507
0
508
sage: E.rank()
509
0
510
511
In this case the rank was computed using a second descent,
512
which is able to determine (by considering a 2-isogenous
513
curve) that Sha is nontrivial. If we deliberately stop the
514
second descent, the rank bound is larger::
515
516
sage: E = mwrank_EllipticCurve([0, -1, 0, -900, -10098])
517
sage: E.two_descent(second_descent = False, verbose=False)
518
sage: E.rank_bound()
519
2
520
521
In contrast, for the curve 571A, also with rank 0 and Sha
522
of order 4, we only obtain an upper bound of 2::
523
524
sage: E = mwrank_EllipticCurve([0, -1, 1, -929, -10595])
525
sage: E.rank_bound()
526
2
527
528
In this case the value returned by :meth:`rank()` is only a
529
lower bound in general (though this is correct)::
530
531
sage: E.rank()
532
0
533
sage: E.certain()
534
False
535
536
"""
537
return self.__two_descent_data().getrankbound()
538
539
def selmer_rank(self):
540
r"""
541
Returns the rank of the 2-Selmer group of the curve.
542
543
EXAMPLES:
544
545
The following is the curve 960D1, which has rank 0, but Sha of
546
order 4. The 2-torsion has rank 2, and the Selmer rank is 3::
547
548
sage: E = mwrank_EllipticCurve([0, -1, 0, -900, -10098])
549
sage: E.selmer_rank()
550
3
551
552
Nevertheless, we can obtain a tight upper bound on the rank
553
since a second descent is performed which establishes the
554
2-rank of Sha::
555
556
sage: E.rank_bound()
557
0
558
559
To show that this was resolved using a second descent, we do
560
the computation again but turn off ``second_descent``::
561
562
sage: E = mwrank_EllipticCurve([0, -1, 0, -900, -10098])
563
sage: E.two_descent(second_descent = False, verbose=False)
564
sage: E.rank_bound()
565
2
566
567
For the curve 571A, also with rank 0 and Sha of order 4,
568
but with no 2-torsion, the Selmer rank is strictly greater
569
than the rank::
570
571
sage: E = mwrank_EllipticCurve([0, -1, 1, -929, -10595])
572
sage: E.selmer_rank()
573
2
574
sage: E.rank_bound()
575
2
576
577
In cases like this with no 2-torsion, the rank upper bound is
578
always equal to the 2-Selmer rank. If we ask for the rank,
579
all we get is a lower bound::
580
581
sage: E.rank()
582
0
583
sage: E.certain()
584
False
585
586
"""
587
return self.__two_descent_data().getselmer()
588
589
def regulator(self):
590
r"""
591
Return the regulator of the saturated Mordell-Weil group.
592
593
EXAMPLES::
594
595
sage: E = mwrank_EllipticCurve([0, 0, 1, -1, 0])
596
sage: E.regulator()
597
0.05111140823996884
598
"""
599
self.saturate()
600
if not self.certain():
601
raise RuntimeError, "Unable to saturate Mordell-Weil group."
602
R = self.__two_descent_data().regulator()
603
return float(R)
604
605
def saturate(self, bound=-1):
606
"""
607
Compute the saturation of the Mordell-Weil group at all
608
primes up to ``bound``.
609
610
INPUT:
611
612
- ``bound`` (int, default -1) -- Use `-1` (the default) to
613
saturate at *all* primes, `0` for no saturation, or `n` (a
614
positive integer) to saturate at all primes up to `n`.
615
616
EXAMPLES:
617
618
Since the 2-descent automatically saturates at primes up to
619
20, it is not easy to come up with an example where saturation
620
has any effect::
621
622
sage: E = mwrank_EllipticCurve([0, 0, 0, -1002231243161, 0])
623
sage: E.gens()
624
[[-1001107, -4004428, 1]]
625
sage: E.saturate()
626
sage: E.gens()
627
[[-1001107, -4004428, 1]]
628
629
"""
630
bound = int(bound)
631
if self.__saturate < bound:
632
self.__two_descent_data().saturate(bound)
633
self.__saturate = bound
634
635
def gens(self):
636
"""
637
Return a list of the generators for the Mordell-Weil group.
638
639
EXAMPLES::
640
641
sage: E = mwrank_EllipticCurve([0, 0, 1, -1, 0])
642
sage: E.gens()
643
[[0, -1, 1]]
644
"""
645
self.saturate()
646
from sage.misc.all import preparse
647
from sage.rings.all import Integer
648
return eval(preparse(self.__two_descent_data().getbasis().replace(":",",")))
649
650
def certain(self):
651
r"""
652
Returns ``True`` if the last :meth:`two_descent()` call provably correctly
653
computed the rank. If :meth:`two_descent()` hasn't been
654
called, then it is first called by :meth:`certain()`
655
using the default parameters.
656
657
The result is ``True`` if and only if the results of the methods
658
:meth:`rank()` and :meth:`rank_bound()` are equal.
659
660
EXAMPLES:
661
662
A 2-descent does not determine `E(\QQ)` with certainty
663
for the curve `y^2 + y = x^3 - x^2 - 120x - 2183`::
664
665
sage: E = mwrank_EllipticCurve([0, -1, 1, -120, -2183])
666
sage: E.two_descent(False)
667
...
668
sage: E.certain()
669
False
670
sage: E.rank()
671
0
672
673
The previous value is only a lower bound; the upper bound is greater::
674
675
sage: E.rank_bound()
676
2
677
678
In fact the rank of `E` is actually 0 (as one could see by
679
computing the `L`-function), but Sha has order 4 and the
680
2-torsion is trivial, so mwrank cannot conclusively
681
determine the rank in this case.
682
"""
683
return bool(self.__two_descent_data().getcertain())
684
685
#def fullmw(self):
686
# return self.__two_descent_data().getfullmw()
687
688
def CPS_height_bound(self):
689
r"""
690
Return the Cremona-Prickett-Siksek height bound. This is a
691
floating point number `B` such that if `P` is a point on the
692
curve, then the naive logarithmic height `h(P)` is less than
693
`B+\hat{h}(P)`, where `\hat{h}(P)` is the canonical height of
694
`P`.
695
696
.. warning::
697
698
We assume the model is minimal!
699
700
EXAMPLES::
701
702
sage: E = mwrank_EllipticCurve([0, 0, 0, -1002231243161, 0])
703
sage: E.CPS_height_bound()
704
14.163198527061496
705
sage: E = mwrank_EllipticCurve([0,0,1,-7,6])
706
sage: E.CPS_height_bound()
707
0.0
708
"""
709
return self.__curve.cps_bound()
710
711
def silverman_bound(self):
712
r"""
713
Return the Silverman height bound. This is a floating point
714
number `B` such that if `P` is a point on the curve, then the
715
naive logarithmic height `h(P)` is less than `B+\hat{h}(P)`,
716
where `\hat{h}(P)` is the canonical height of `P`.
717
718
.. warning::
719
720
We assume the model is minimal!
721
722
EXAMPLES::
723
724
sage: E = mwrank_EllipticCurve([0, 0, 0, -1002231243161, 0])
725
sage: E.silverman_bound()
726
18.29545210468247
727
sage: E = mwrank_EllipticCurve([0,0,1,-7,6])
728
sage: E.silverman_bound()
729
6.284833369972403
730
"""
731
return self.__curve.silverman_bound()
732
733
734
class mwrank_MordellWeil(SageObject):
735
r"""
736
The :class:`mwrank_MordellWeil` class represents a subgroup of a
737
Mordell-Weil group. Use this class to saturate a specified list
738
of points on an :class:`mwrank_EllipticCurve`, or to search for
739
points up to some bound.
740
741
INPUT:
742
743
- ``curve`` (:class:`mwrank_EllipticCurve`) -- the underlying
744
elliptic curve.
745
746
- ``verbose`` (bool, default ``False``) -- verbosity flag (controls
747
amount of output produced in point searches).
748
749
- ``pp`` (int, default 1) -- process points flag (if nonzero,
750
the points found are processed, so that at all times only a
751
`\ZZ`-basis for the subgroup generated by the points found
752
so far is stored; if zero, no processing is done and all
753
points found are stored).
754
755
- ``maxr`` (int, default 999) -- maximum rank (quit point
756
searching once the points found generate a subgroup of this
757
rank; useful if an upper bound for the rank is already
758
known).
759
760
EXAMPLE::
761
762
sage: E = mwrank_EllipticCurve([1,0,1,4,-6])
763
sage: EQ = mwrank_MordellWeil(E)
764
sage: EQ
765
Subgroup of Mordell-Weil group: []
766
sage: EQ.search(2)
767
P1 = [0:1:0] is torsion point, order 1
768
P1 = [1:-1:1] is torsion point, order 2
769
P1 = [2:2:1] is torsion point, order 3
770
P1 = [9:23:1] is torsion point, order 6
771
772
sage: E = mwrank_EllipticCurve([0,0,1,-7,6])
773
sage: EQ = mwrank_MordellWeil(E)
774
sage: EQ.search(2)
775
P1 = [0:1:0] is torsion point, order 1
776
P1 = [-3:0:1] is generator number 1
777
...
778
P4 = [-91:804:343] = -2*P1 + 2*P2 + 1*P3 (mod torsion)
779
sage: EQ
780
Subgroup of Mordell-Weil group: [[1:-1:1], [-2:3:1], [-14:25:8]]
781
782
Example to illustrate the verbose parameter::
783
784
sage: E = mwrank_EllipticCurve([0,0,1,-7,6])
785
sage: EQ = mwrank_MordellWeil(E, verbose=False)
786
sage: EQ.search(1)
787
sage: EQ
788
Subgroup of Mordell-Weil group: [[1:-1:1], [-2:3:1], [-14:25:8]]
789
790
sage: EQ = mwrank_MordellWeil(E, verbose=True)
791
sage: EQ.search(1)
792
P1 = [0:1:0] is torsion point, order 1
793
P1 = [-3:0:1] is generator number 1
794
saturating up to 20...Checking 2-saturation
795
Points have successfully been 2-saturated (max q used = 7)
796
Checking 3-saturation
797
Points have successfully been 3-saturated (max q used = 7)
798
Checking 5-saturation
799
Points have successfully been 5-saturated (max q used = 23)
800
Checking 7-saturation
801
Points have successfully been 7-saturated (max q used = 41)
802
Checking 11-saturation
803
Points have successfully been 11-saturated (max q used = 17)
804
Checking 13-saturation
805
Points have successfully been 13-saturated (max q used = 43)
806
Checking 17-saturation
807
Points have successfully been 17-saturated (max q used = 31)
808
Checking 19-saturation
809
Points have successfully been 19-saturated (max q used = 37)
810
done
811
P2 = [-2:3:1] is generator number 2
812
saturating up to 20...Checking 2-saturation
813
possible kernel vector = [1,1]
814
This point may be in 2E(Q): [14:-52:1]
815
...and it is!
816
Replacing old generator #1 with new generator [1:-1:1]
817
Points have successfully been 2-saturated (max q used = 7)
818
Index gain = 2^1
819
Checking 3-saturation
820
Points have successfully been 3-saturated (max q used = 13)
821
Checking 5-saturation
822
Points have successfully been 5-saturated (max q used = 67)
823
Checking 7-saturation
824
Points have successfully been 7-saturated (max q used = 53)
825
Checking 11-saturation
826
Points have successfully been 11-saturated (max q used = 73)
827
Checking 13-saturation
828
Points have successfully been 13-saturated (max q used = 103)
829
Checking 17-saturation
830
Points have successfully been 17-saturated (max q used = 113)
831
Checking 19-saturation
832
Points have successfully been 19-saturated (max q used = 47)
833
done (index = 2).
834
Gained index 2, new generators = [ [1:-1:1] [-2:3:1] ]
835
P3 = [-14:25:8] is generator number 3
836
saturating up to 20...Checking 2-saturation
837
Points have successfully been 2-saturated (max q used = 11)
838
Checking 3-saturation
839
Points have successfully been 3-saturated (max q used = 13)
840
Checking 5-saturation
841
Points have successfully been 5-saturated (max q used = 71)
842
Checking 7-saturation
843
Points have successfully been 7-saturated (max q used = 101)
844
Checking 11-saturation
845
Points have successfully been 11-saturated (max q used = 127)
846
Checking 13-saturation
847
Points have successfully been 13-saturated (max q used = 151)
848
Checking 17-saturation
849
Points have successfully been 17-saturated (max q used = 139)
850
Checking 19-saturation
851
Points have successfully been 19-saturated (max q used = 179)
852
done (index = 1).
853
P4 = [-1:3:1] = -1*P1 + -1*P2 + -1*P3 (mod torsion)
854
P4 = [0:2:1] = 2*P1 + 0*P2 + 1*P3 (mod torsion)
855
P4 = [2:13:8] = -3*P1 + 1*P2 + -1*P3 (mod torsion)
856
P4 = [1:0:1] = -1*P1 + 0*P2 + 0*P3 (mod torsion)
857
P4 = [2:0:1] = -1*P1 + 1*P2 + 0*P3 (mod torsion)
858
P4 = [18:7:8] = -2*P1 + -1*P2 + -1*P3 (mod torsion)
859
P4 = [3:3:1] = 1*P1 + 0*P2 + 1*P3 (mod torsion)
860
P4 = [4:6:1] = 0*P1 + -1*P2 + -1*P3 (mod torsion)
861
P4 = [36:69:64] = 1*P1 + -2*P2 + 0*P3 (mod torsion)
862
P4 = [68:-25:64] = -2*P1 + -1*P2 + -2*P3 (mod torsion)
863
P4 = [12:35:27] = 1*P1 + -1*P2 + -1*P3 (mod torsion)
864
sage: EQ
865
Subgroup of Mordell-Weil group: [[1:-1:1], [-2:3:1], [-14:25:8]]
866
867
Example to illustrate the process points (``pp``) parameter::
868
869
sage: E = mwrank_EllipticCurve([0,0,1,-7,6])
870
sage: EQ = mwrank_MordellWeil(E, verbose=False, pp=1)
871
sage: EQ.search(1); EQ # generators only
872
Subgroup of Mordell-Weil group: [[1:-1:1], [-2:3:1], [-14:25:8]]
873
sage: EQ = mwrank_MordellWeil(E, verbose=False, pp=0)
874
sage: EQ.search(1); EQ # all points found
875
Subgroup of Mordell-Weil group: [[-3:0:1], [-2:3:1], [-14:25:8], [-1:3:1], [0:2:1], [2:13:8], [1:0:1], [2:0:1], [18:7:8], [3:3:1], [4:6:1], [36:69:64], [68:-25:64], [12:35:27]]
876
"""
877
878
def __init__(self, curve, verbose=True, pp=1, maxr=999):
879
r"""
880
Constructor for the :class:`mwrank_MordellWeil` class.
881
882
See the docstring of this class for full documentation.
883
884
EXAMPLE::
885
886
sage: E = mwrank_EllipticCurve([1,0,1,4,-6])
887
sage: EQ = mwrank_MordellWeil(E)
888
sage: EQ
889
Subgroup of Mordell-Weil group: []
890
"""
891
if not isinstance(curve, mwrank_EllipticCurve):
892
raise TypeError, "curve (=%s) must be an mwrank_EllipticCurve"%curve
893
self.__curve = curve
894
self.__verbose = verbose
895
self.__pp = pp
896
self.__maxr = maxr
897
if verbose:
898
verb = 1
899
else:
900
verb = 0
901
from sage.libs.mwrank.mwrank import _mw # import here to save time
902
self.__mw = _mw(curve._curve_data(), verb, pp, maxr)
903
904
def __reduce__(self):
905
r"""
906
Standard Python function used in pickling.
907
908
EXAMPLES::
909
910
sage: E = mwrank_EllipticCurve([0,0,1,-7,6])
911
sage: EQ = mwrank_MordellWeil(E)
912
sage: EQ.__reduce__()
913
(<class 'sage.libs.mwrank.interface.mwrank_MordellWeil'>, (y^2+ y = x^3 - 7*x + 6, True, 1, 999))
914
"""
915
return mwrank_MordellWeil, (self.__curve, self.__verbose, self.__pp, self.__maxr)
916
917
def __repr__(self):
918
r"""
919
String representation of this Mordell-Weil subgroup.
920
921
OUTPUT:
922
923
(string) String representation of this Mordell-Weil subgroup.
924
925
EXAMPLE::
926
927
sage: E = mwrank_EllipticCurve([0,0,1,-7,6])
928
sage: EQ = mwrank_MordellWeil(E, verbose=False)
929
sage: EQ.__repr__()
930
'Subgroup of Mordell-Weil group: []'
931
sage: EQ.search(1)
932
sage: EQ.__repr__()
933
'Subgroup of Mordell-Weil group: [[1:-1:1], [-2:3:1], [-14:25:8]]'
934
"""
935
# FIXED: Changing "Mordell Weil" to "Mordell-Weil" did break doctests...
936
return "Subgroup of Mordell-Weil group: %s"%self.__mw
937
938
def process(self, v, sat=0):
939
"""
940
This function allows one to add points to a :class:`mwrank_MordellWeil` object.
941
942
Process points in the list ``v``, with saturation at primes up to
943
``sat``. If ``sat`` is zero (the default), do no saturation.
944
945
INPUT:
946
947
- ``v`` (list of 3-tuples or lists of ints or Integers) -- a
948
list of triples of integers, which define points on the
949
curve.
950
951
- ``sat`` (int, default 0) -- saturate at primes up to ``sat``, or at
952
*all* primes if ``sat`` is zero.
953
954
OUTPUT:
955
956
None. But note that if the ``verbose`` flag is set, then there
957
will be some output as a side-effect.
958
959
EXAMPLES::
960
961
sage: E = mwrank_EllipticCurve([0,0,1,-7,6])
962
sage: E.gens()
963
[[1, -1, 1], [-2, 3, 1], [-14, 25, 8]]
964
sage: EQ = mwrank_MordellWeil(E)
965
sage: EQ.process([[1, -1, 1], [-2, 3, 1], [-14, 25, 8]])
966
P1 = [1:-1:1] is generator number 1
967
P2 = [-2:3:1] is generator number 2
968
P3 = [-14:25:8] is generator number 3
969
970
::
971
972
sage: EQ.points()
973
[[1, -1, 1], [-2, 3, 1], [-14, 25, 8]]
974
975
Example to illustrate the saturation parameter ``sat``::
976
977
sage: E = mwrank_EllipticCurve([0,0,1,-7,6])
978
sage: EQ = mwrank_MordellWeil(E)
979
sage: EQ.process([[1547, -2967, 343], [2707496766203306, 864581029138191, 2969715140223272], [-13422227300, -49322830557, 12167000000]], sat=20)
980
P1 = [1547:-2967:343] is generator number 1
981
...
982
Gained index 5, new generators = [ [-2:3:1] [-14:25:8] [1:-1:1] ]
983
984
sage: EQ.points()
985
[[-2, 3, 1], [-14, 25, 8], [1, -1, 1]]
986
987
Here the processing was followed by saturation at primes up to
988
20. Now we prevent this initial saturation::
989
990
sage: E = mwrank_EllipticCurve([0,0,1,-7,6])
991
sage: EQ = mwrank_MordellWeil(E)
992
sage: EQ.process([[1547, -2967, 343], [2707496766203306, 864581029138191, 2969715140223272], [-13422227300, -49322830557, 12167000000]], sat=0)
993
P1 = [1547:-2967:343] is generator number 1
994
P2 = [2707496766203306:864581029138191:2969715140223272] is generator number 2
995
P3 = [-13422227300:-49322830557:12167000000] is generator number 3
996
sage: EQ.points()
997
[[1547, -2967, 343], [2707496766203306, 864581029138191, 2969715140223272], [-13422227300, -49322830557, 12167000000]]
998
sage: EQ.regulator()
999
375.42919921875
1000
sage: EQ.saturate(2) # points were not 2-saturated
1001
saturating basis...Saturation index bound = 93
1002
WARNING: saturation at primes p > 2 will not be done;
1003
...
1004
Gained index 2
1005
New regulator = 93.857300720636393209
1006
(False, 2, '[ ]')
1007
sage: EQ.points()
1008
[[-2, 3, 1], [2707496766203306, 864581029138191, 2969715140223272], [-13422227300, -49322830557, 12167000000]]
1009
sage: EQ.regulator()
1010
93.8572998046875
1011
sage: EQ.saturate(3) # points were not 3-saturated
1012
saturating basis...Saturation index bound = 46
1013
WARNING: saturation at primes p > 3 will not be done;
1014
...
1015
Gained index 3
1016
New regulator = 10.4285889689595992455
1017
(False, 3, '[ ]')
1018
sage: EQ.points()
1019
[[-2, 3, 1], [-14, 25, 8], [-13422227300, -49322830557, 12167000000]]
1020
sage: EQ.regulator()
1021
10.4285888671875
1022
sage: EQ.saturate(5) # points were not 5-saturated
1023
saturating basis...Saturation index bound = 15
1024
WARNING: saturation at primes p > 5 will not be done;
1025
...
1026
Gained index 5
1027
New regulator = 0.417143558758383969818
1028
(False, 5, '[ ]')
1029
sage: EQ.points()
1030
[[-2, 3, 1], [-14, 25, 8], [1, -1, 1]]
1031
sage: EQ.regulator()
1032
0.4171435534954071
1033
sage: EQ.saturate() # points are now saturated
1034
saturating basis...Saturation index bound = 3
1035
Checking saturation at [ 2 3 ]
1036
Checking 2-saturation
1037
Points were proved 2-saturated (max q used = 11)
1038
Checking 3-saturation
1039
Points were proved 3-saturated (max q used = 13)
1040
done
1041
(True, 1, '[ ]')
1042
"""
1043
if not isinstance(v, list):
1044
raise TypeError, "v (=%s) must be a list"%v
1045
sat = int(sat)
1046
for P in v:
1047
if not isinstance(P, (list,tuple)) or len(P) != 3:
1048
raise TypeError, "v (=%s) must be a list of 3-tuples (or 3-element lists) of ints"%v
1049
self.__mw.process(P, sat)
1050
1051
def regulator(self):
1052
"""
1053
Return the regulator of the points in this subgroup of
1054
the Mordell-Weil group.
1055
1056
.. note::
1057
1058
``eclib`` can compute the regulator to arbitrary precision,
1059
but the interface currently returns the output as a ``float``.
1060
1061
OUTPUT:
1062
1063
(float) The regulator of the points in this subgroup.
1064
1065
EXAMPLES::
1066
1067
sage: E = mwrank_EllipticCurve([0,-1,1,0,0])
1068
sage: E.regulator()
1069
1.0
1070
1071
sage: E = mwrank_EllipticCurve([0,0,1,-7,6])
1072
sage: E.regulator()
1073
0.417143558758384
1074
"""
1075
return self.__mw.regulator()
1076
1077
def rank(self):
1078
"""
1079
Return the rank of this subgroup of the Mordell-Weil group.
1080
1081
OUTPUT:
1082
1083
(int) The rank of this subgroup of the Mordell-Weil group.
1084
1085
EXAMPLES::
1086
1087
sage: E = mwrank_EllipticCurve([0,-1,1,0,0])
1088
sage: E.rank()
1089
0
1090
1091
A rank 3 example::
1092
1093
sage: E = mwrank_EllipticCurve([0,0,1,-7,6])
1094
sage: EQ = mwrank_MordellWeil(E)
1095
sage: EQ.rank()
1096
0
1097
sage: EQ.regulator()
1098
1.0
1099
1100
The preceding output is correct, since we have not yet tried
1101
to find any points on the curve either by searching or
1102
2-descent::
1103
1104
sage: EQ
1105
Subgroup of Mordell-Weil group: []
1106
1107
Now we do a very small search::
1108
1109
sage: EQ.search(1)
1110
P1 = [0:1:0] is torsion point, order 1
1111
P1 = [-3:0:1] is generator number 1
1112
saturating up to 20...Checking 2-saturation
1113
...
1114
P4 = [12:35:27] = 1*P1 + -1*P2 + -1*P3 (mod torsion)
1115
sage: EQ
1116
Subgroup of Mordell-Weil group: [[1:-1:1], [-2:3:1], [-14:25:8]]
1117
sage: EQ.rank()
1118
3
1119
sage: EQ.regulator()
1120
0.4171435534954071
1121
1122
We do in fact now have a full Mordell-Weil basis.
1123
1124
"""
1125
return self.__mw.rank()
1126
1127
def saturate(self, max_prime=-1, odd_primes_only=False):
1128
r"""
1129
Saturate this subgroup of the Mordell-Weil group.
1130
1131
INPUT:
1132
1133
- ``max_prime`` (int, default -1) -- saturation is performed for
1134
all primes up to ``max_prime``. If `-1` (the default), an
1135
upper bound is computed for the primes at which the subgroup
1136
may not be saturated, and this is used; however, if the
1137
computed bound is greater than a value set by the ``eclib``
1138
library (currently 97) then no saturation will be attempted
1139
at primes above this.
1140
1141
- ``odd_primes_only`` (bool, default ``False``) -- only do
1142
saturation at odd primes. (If the points have been found
1143
via :meth:``two_descent()`` they should already be 2-saturated.)
1144
1145
OUTPUT:
1146
1147
(3-tuple) (``ok``, ``index``, ``unsatlist``) where:
1148
1149
- ``ok`` (bool) -- ``True`` if and only if the saturation was
1150
provably successful at all primes attempted. If the default
1151
was used for ``max_prime`` and no warning was output about
1152
the computed saturation bound being too high, then ``True``
1153
indicates that the subgroup is saturated at *all*
1154
primes.
1155
1156
- ``index`` (int) -- the index of the group generated by the
1157
original points in their saturation.
1158
1159
- ``unsatlist`` (list of ints) -- list of primes at which
1160
saturation could not be proved or achieved. Increasing the
1161
decimal precision should correct this, since it happens when
1162
a linear combination of the points appears to be a multiple
1163
of `p` but cannot be divided by `p`. (Note that ``eclib``
1164
uses floating point methods based on elliptic logarithms to
1165
divide points.)
1166
1167
.. note::
1168
1169
We emphasize that if this function returns ``True`` as the
1170
first return argument (``ok``), and if the default was used for the
1171
parameter ``max_prime``, then the points in the basis after
1172
calling this function are saturated at *all* primes,
1173
i.e., saturating at the primes up to ``max_prime`` are
1174
sufficient to saturate at all primes. Note that the
1175
function might not have needed to saturate at all primes up
1176
to ``max_prime``. It has worked out what prime you need to
1177
saturate up to, and that prime might be smaller than ``max_prime``.
1178
1179
.. note::
1180
1181
Currently (May 2010), this does not remember the result of
1182
calling :meth:`search()`. So calling :meth:`search()` up
1183
to height 20 then calling :meth:`saturate()` results in
1184
another search up to height 18.
1185
1186
EXAMPLES::
1187
1188
sage: E = mwrank_EllipticCurve([0,0,1,-7,6])
1189
sage: EQ = mwrank_MordellWeil(E)
1190
1191
We initialise with three points which happen to be 2, 3 and 5
1192
times the generators of this rank 3 curve. To prevent
1193
automatic saturation at this stage we set the parameter
1194
``sat`` to 0 (which is in fact the default)::
1195
1196
sage: EQ.process([[1547, -2967, 343], [2707496766203306, 864581029138191, 2969715140223272], [-13422227300, -49322830557, 12167000000]], sat=0)
1197
P1 = [1547:-2967:343] is generator number 1
1198
P2 = [2707496766203306:864581029138191:2969715140223272] is generator number 2
1199
P3 = [-13422227300:-49322830557:12167000000] is generator number 3
1200
sage: EQ
1201
Subgroup of Mordell-Weil group: [[1547:-2967:343], [2707496766203306:864581029138191:2969715140223272], [-13422227300:-49322830557:12167000000]]
1202
sage: EQ.regulator()
1203
375.42919921875
1204
1205
Now we saturate at `p=2`, and gain index 2::
1206
1207
sage: EQ.saturate(2) # points were not 2-saturated
1208
saturating basis...Saturation index bound = 93
1209
WARNING: saturation at primes p > 2 will not be done;
1210
...
1211
Gained index 2
1212
New regulator = 93.857300720636393209
1213
(False, 2, '[ ]')
1214
sage: EQ
1215
Subgroup of Mordell-Weil group: [[-2:3:1], [2707496766203306:864581029138191:2969715140223272], [-13422227300:-49322830557:12167000000]]
1216
sage: EQ.regulator()
1217
93.8572998046875
1218
1219
Now we saturate at `p=3`, and gain index 3::
1220
1221
sage: EQ.saturate(3) # points were not 3-saturated
1222
saturating basis...Saturation index bound = 46
1223
WARNING: saturation at primes p > 3 will not be done;
1224
...
1225
Gained index 3
1226
New regulator = 10.4285889689595992455
1227
(False, 3, '[ ]')
1228
sage: EQ
1229
Subgroup of Mordell-Weil group: [[-2:3:1], [-14:25:8], [-13422227300:-49322830557:12167000000]]
1230
sage: EQ.regulator()
1231
10.4285888671875
1232
1233
Now we saturate at `p=5`, and gain index 5::
1234
1235
sage: EQ.saturate(5) # points were not 5-saturated
1236
saturating basis...Saturation index bound = 15
1237
WARNING: saturation at primes p > 5 will not be done;
1238
...
1239
Gained index 5
1240
New regulator = 0.417143558758383969818
1241
(False, 5, '[ ]')
1242
sage: EQ
1243
Subgroup of Mordell-Weil group: [[-2:3:1], [-14:25:8], [1:-1:1]]
1244
sage: EQ.regulator()
1245
0.4171435534954071
1246
1247
Finally we finish the saturation. The output here shows that
1248
the points are now provably saturated at all primes::
1249
1250
sage: EQ.saturate() # points are now saturated
1251
saturating basis...Saturation index bound = 3
1252
Checking saturation at [ 2 3 ]
1253
Checking 2-saturation
1254
Points were proved 2-saturated (max q used = 11)
1255
Checking 3-saturation
1256
Points were proved 3-saturated (max q used = 13)
1257
done
1258
(True, 1, '[ ]')
1259
1260
Of course, the :meth:`process()` function would have done all this
1261
automatically for us::
1262
1263
sage: E = mwrank_EllipticCurve([0,0,1,-7,6])
1264
sage: EQ = mwrank_MordellWeil(E)
1265
sage: EQ.process([[1547, -2967, 343], [2707496766203306, 864581029138191, 2969715140223272], [-13422227300, -49322830557, 12167000000]], sat=5)
1266
P1 = [1547:-2967:343] is generator number 1
1267
...
1268
Gained index 5, new generators = [ [-2:3:1] [-14:25:8] [1:-1:1] ]
1269
sage: EQ
1270
Subgroup of Mordell-Weil group: [[-2:3:1], [-14:25:8], [1:-1:1]]
1271
sage: EQ.regulator()
1272
0.4171435534954071
1273
1274
But we would still need to use the :meth:`saturate()` function to
1275
verify that full saturation has been done::
1276
1277
sage: EQ.saturate()
1278
saturating basis...Saturation index bound = 3
1279
Checking saturation at [ 2 3 ]
1280
Checking 2-saturation
1281
Points were proved 2-saturated (max q used = 11)
1282
Checking 3-saturation
1283
Points were proved 3-saturated (max q used = 13)
1284
done
1285
(True, 1, '[ ]')
1286
1287
Note the output of the preceding command: it proves that the
1288
index of the points in their saturation is at most 3, then
1289
proves saturation at 2 and at 3, by reducing the points modulo
1290
all primes of good reduction up to 11, respectively 13.
1291
"""
1292
ok, index, unsat = self.__mw.saturate(int(max_prime), odd_primes_only)
1293
return bool(ok), int(str(index)), unsat
1294
1295
def search(self, height_limit=18, verbose=False):
1296
r"""
1297
Search for new points, and add them to this subgroup of the
1298
Mordell-Weil group.
1299
1300
INPUT:
1301
1302
- ``height_limit`` (float, default: 18) -- search up to this
1303
logarithmic height.
1304
1305
.. note::
1306
1307
On 32-bit machines, this *must* be < 21.48 else
1308
`\exp(h_{\text{lim}}) > 2^{31}` and overflows. On 64-bit machines, it
1309
must be *at most* 43.668. However, this bound is a logarithmic
1310
bound and increasing it by just 1 increases the running time
1311
by (roughly) `\exp(1.5)=4.5`, so searching up to even 20
1312
takes a very long time.
1313
1314
.. note::
1315
1316
The search is carried out with a quadratic sieve, using
1317
code adapted from a version of Michael Stoll's
1318
``ratpoints`` program. It would be preferable to use a
1319
newer version of ``ratpoints``.
1320
1321
- ``verbose`` (bool, default ``False``) -- turn verbose operation on
1322
or off.
1323
1324
EXAMPLES:
1325
1326
A rank 3 example, where a very small search is sufficient to
1327
find a Mordell-Weil basis::
1328
1329
sage: E = mwrank_EllipticCurve([0,0,1,-7,6])
1330
sage: EQ = mwrank_MordellWeil(E)
1331
sage: EQ.search(1)
1332
P1 = [0:1:0] is torsion point, order 1
1333
P1 = [-3:0:1] is generator number 1
1334
...
1335
P4 = [12:35:27] = 1*P1 + -1*P2 + -1*P3 (mod torsion)
1336
sage: EQ
1337
Subgroup of Mordell-Weil group: [[1:-1:1], [-2:3:1], [-14:25:8]]
1338
1339
In the next example, a search bound of 12 is needed to find a
1340
non-torsion point::
1341
1342
sage: E = mwrank_EllipticCurve([0, -1, 0, -18392, -1186248]) #1056g4
1343
sage: EQ = mwrank_MordellWeil(E)
1344
sage: EQ.search(11); EQ
1345
P1 = [0:1:0] is torsion point, order 1
1346
P1 = [161:0:1] is torsion point, order 2
1347
Subgroup of Mordell-Weil group: []
1348
sage: EQ.search(12); EQ
1349
P1 = [0:1:0] is torsion point, order 1
1350
P1 = [161:0:1] is torsion point, order 2
1351
P1 = [4413270:10381877:27000] is generator number 1
1352
...
1353
Subgroup of Mordell-Weil group: [[4413270:10381877:27000]]
1354
"""
1355
height_limit = float(height_limit)
1356
if height_limit >= 21.4: # TODO: docstring says 21.48 (for 32-bit machines; what about 64-bit...?)
1357
raise ValueError, "The height limit must be < 21.4."
1358
1359
moduli_option = 0 # Use Stoll's sieving program... see strategies in ratpoints-1.4.c
1360
1361
## moduli_option -- int (default: 0); if > 0; a flag used to determine
1362
## the moduli that are used in sieving
1363
## 1 -- first 10 odd primes; the first one you
1364
## would think of.
1365
## 2 -- three composites; $2^6\cdot 3^4$, ... TODO
1366
## (from German mathematician J. Gebel;
1367
## personal conversation about SIMATH)
1368
## 3 -- nine prime powers; $2^5, \ldots$;
1369
## like 1 but includes powers of small primes
1370
## TODO: Extract the meaning from mwprocs.cc; line 776 etc.
1371
1372
verbose = bool(verbose)
1373
self.__mw.search(height_limit, moduli_option, verbose)
1374
1375
def points(self):
1376
"""
1377
Return a list of the generating points in this Mordell-Weil
1378
group.
1379
1380
OUTPUT:
1381
1382
(list) A list of lists of length 3, each holding the
1383
primitive integer coordinates `[x,y,z]` of a generating
1384
point.
1385
1386
EXAMPLES::
1387
1388
sage: E = mwrank_EllipticCurve([0,0,1,-7,6])
1389
sage: EQ = mwrank_MordellWeil(E)
1390
sage: EQ.search(1)
1391
P1 = [0:1:0] is torsion point, order 1
1392
P1 = [-3:0:1] is generator number 1
1393
...
1394
P4 = [12:35:27] = 1*P1 + -1*P2 + -1*P3 (mod torsion)
1395
sage: EQ.points()
1396
[[1, -1, 1], [-2, 3, 1], [-14, 25, 8]]
1397
1398
"""
1399
from sage.misc.all import preparse
1400
from sage.rings.all import Integer
1401
return eval(preparse(self.__mw.getbasis().replace(":",",")))
1402
1403
1404
1405