Book a Demo!
CoCalc Logo Icon
StoreFeaturesDocsShareSupportNewsAboutPoliciesSign UpSign In
sagemath
GitHub Repository: sagemath/sagelib
Path: blob/master/sage/lfunctions/lcalc.py
4066 views
1
r"""
2
Rubinstein's `L`-function Calculator
3
4
This a standard part of Sage. This interface provides complete
5
access to Rubinstein's lcalc calculator with extra PARI
6
functionality compiled in.
7
8
.. note::
9
10
Each call to ``lcalc`` runs a complete
11
``lcalc`` process. On a typical Linux system, this
12
entails about 0.3 seconds overhead.
13
14
AUTHORS:
15
16
- Michael Rubinstein (2005): released under GPL the C++ program lcalc
17
18
- William Stein (2006-03-05): wrote Sage interface to lcalc
19
"""
20
21
########################################################################
22
# Copyright (C) 2006 William Stein <[email protected]>
23
#
24
# Distributed under the terms of the GNU General Public License (GPL)
25
#
26
# http://www.gnu.org/licenses/
27
########################################################################
28
29
import os
30
from sage.structure.sage_object import SageObject
31
from sage.misc.all import pager
32
import sage.rings.all
33
34
prec = 32
35
36
class LCalc(SageObject):
37
r"""
38
Rubinstein's `L`-functions Calculator
39
40
Type ``lcalc.[tab]`` for a list of useful commands that
41
are implemented using the command line interface, but return
42
objects that make sense in Sage. For each command the possible
43
inputs for the L-function are:
44
45
46
- ``"`` - (default) the Riemann zeta function
47
48
- ``'tau'`` - the L function of the Ramanujan delta
49
function
50
51
- elliptic curve E - where E is an elliptic curve over
52
`\mathbb{Q}`; defines `L(E,s)`
53
54
55
You can also use the complete command-line interface of
56
Rubinstein's `L`-functions calculations program via this
57
class. Type ``lcalc.help()`` for a list of commands and
58
how to call them.
59
"""
60
def _repr_(self):
61
return "Rubinsteins L-function Calculator"
62
63
def __call__(self, args):
64
cmd = 'lcalc %s'%args
65
return os.popen(cmd).read().strip()
66
67
def _compute_L(self, L):
68
if isinstance(L, str):
69
if L == 'tau':
70
return '--tau'
71
return L
72
import sage.schemes.all
73
if sage.schemes.all.is_EllipticCurve(L):
74
if L.base_ring() == sage.rings.all.RationalField():
75
L = L.minimal_model()
76
return '-e --a1 %s --a2 %s --a3 %s --a4 %s --a6 %s'%tuple(L.a_invariants())
77
raise TypeError, "$L$-function of %s not known"%L
78
79
def help(self):
80
try:
81
h = self.__help
82
except AttributeError:
83
h = "-"*70 + '\n'
84
h += " Call lcalc with one argument, e.g., \n"
85
h += " sage: lcalc('--tau -z 1000')\n"
86
h += " is translated into the command line\n"
87
h += " $ lcalc --tau -z 1000\n"
88
h += "-"*70 + '\n'
89
h += '\n' + self('--help')
90
self.__help = h
91
pager()(h)
92
93
def zeros(self, n, L=''):
94
"""
95
Return the imaginary parts of the first `n` nontrivial
96
zeros of the `L`-function in the upper half plane, as
97
32-bit reals.
98
99
INPUT:
100
101
102
- ``n`` - integer
103
104
- ``L`` - defines `L`-function (default:
105
Riemann zeta function)
106
107
108
This function also checks the Riemann Hypothesis and makes sure no
109
zeros are missed. This means it looks for several dozen zeros to
110
make sure none have been missed before outputting any zeros at all,
111
so takes longer than
112
``self.zeros_of_zeta_in_interval(...)``.
113
114
EXAMPLES::
115
116
sage: lcalc.zeros(4) # long time
117
[14.1347251, 21.0220396, 25.0108576, 30.4248761]
118
sage: lcalc.zeros(5, L='--tau') # long time
119
[9.22237940, 13.9075499, 17.4427770, 19.6565131, 22.3361036]
120
sage: lcalc.zeros(3, EllipticCurve('37a')) # long time
121
[0.000000000, 5.00317001, 6.87039122]
122
"""
123
L = self._compute_L(L)
124
RR = sage.rings.all.RealField(prec)
125
X = self('-z %s %s'%(int(n), L))
126
return [RR(z) for z in X.split()]
127
128
def zeros_in_interval(self, x, y, stepsize, L=''):
129
r"""
130
Return the imaginary parts of (most of) the nontrivial zeros of the
131
`L`-function on the line `\Re(s)=1/2` with positive
132
imaginary part between `x` and `y`, along with a
133
technical quantity for each.
134
135
INPUT:
136
137
138
- ``x, y, stepsize`` - positive floating point
139
numbers
140
141
- ``L`` - defines `L`-function (default:
142
Riemann zeta function)
143
144
145
OUTPUT: list of pairs (zero, S(T)).
146
147
Rubinstein writes: The first column outputs the imaginary part of
148
the zero, the second column a quantity related to `S(T)`
149
(it increases roughly by 2 whenever a sign change, i.e. pair of
150
zeros, is missed). Higher up the critical strip you should use a
151
smaller stepsize so as not to miss zeros.
152
153
EXAMPLES::
154
155
sage: lcalc.zeros_in_interval(10, 30, 0.1)
156
[(14.1347251, 0.184672916), (21.0220396, -0.0677893290), (25.0108576, -0.0555872781)]
157
"""
158
L = self._compute_L(L)
159
RR = sage.rings.all.RealField(prec)
160
X = self('--zeros-interval -x %s -y %s --stepsize=%s %s'%(
161
float(x), float(y), float(stepsize), L))
162
return [tuple([RR(z) for z in t.split()]) for t in X.split('\n')]
163
164
def value(self, s, L=''):
165
r"""
166
Return `L(s)` for `s` a complex number.
167
168
INPUT:
169
170
171
- ``s`` - complex number
172
173
- ``L`` - defines `L`-function (default:
174
Riemann zeta function)
175
176
177
EXAMPLES::
178
179
sage: I = CC.0
180
sage: lcalc.value(0.5 + 100*I)
181
2.69261989 - 0.0203860296*I
182
183
Note, Sage can also compute zeta at complex numbers (using the PARI
184
C library)::
185
186
sage: (0.5 + 100*I).zeta()
187
2.69261988568132 - 0.0203860296025982*I
188
"""
189
L = self._compute_L(L)
190
CC = sage.rings.all.ComplexField(prec)
191
s = CC(s)
192
x, y = self('-v -x %s -y %s %s'%(s.real(), s.imag(), L)).split()
193
return CC((float(x), float(y)))
194
195
def values_along_line(self, s0, s1, number_samples, L=''):
196
r"""
197
Return values of `L(s)` at ``number_samples``
198
equally-spaced sample points along the line from `s_0` to
199
`s_1` in the complex plane.
200
201
INPUT:
202
203
204
- ``s0, s1`` - complex numbers
205
206
- ``number_samples`` - integer
207
208
- ``L`` - defines `L`-function (default:
209
Riemann zeta function)
210
211
212
OUTPUT:
213
214
215
- ``list`` - list of pairs (s, zeta(s)), where the s
216
are equally spaced sampled points on the line from s0 to s1.
217
218
219
EXAMPLES::
220
221
sage: I = CC.0
222
sage: lcalc.values_along_line(0.5, 0.5+20*I, 5)
223
[(0.500000000, -1.46035451), (0.500000000 + 4.00000000*I, 0.606783764 + 0.0911121400*I), (0.500000000 + 8.00000000*I, 1.24161511 + 0.360047588*I), (0.500000000 + 12.0000000*I, 1.01593665 - 0.745112472*I), (0.500000000 + 16.0000000*I, 0.938545408 + 1.21658782*I)]
224
225
Sometimes warnings are printed (by lcalc) when this command is
226
run::
227
228
sage: E = EllipticCurve('389a')
229
sage: E.lseries().values_along_line(0.5, 3, 5)
230
[(0.000000000, 0.209951303),
231
(0.500000000, -...e-16),
232
(1.00000000, 0.133768433),
233
(1.50000000, 0.360092864),
234
(2.00000000, 0.552975867)]
235
"""
236
L = self._compute_L(L)
237
CC = sage.rings.all.ComplexField(prec)
238
s0 = CC(s0)
239
s1 = CC(s1)
240
v = self('--value-line-segment -x %s -y %s -X %s -Y %s --number-samples %s %s'%(
241
(s0.real(), s0.imag(), s1.real(), s1.imag(), int(number_samples), L)))
242
w = []
243
for a in v.split('\n'):
244
try:
245
x0,y0,x1,y1 = a.split()
246
w.append((CC(x0,y0), CC(x1,y1)))
247
except ValueError:
248
print 'lcalc: ', a
249
return w
250
251
def twist_values(self, s, dmin, dmax, L=''):
252
r"""
253
Return values of `L(s, \chi_k)` for each quadratic
254
character `\chi_k` whose discriminant `d` satisfies
255
`d_{\min} \leq d \leq d_{\max}`.
256
257
INPUT:
258
259
260
- ``s`` - complex numbers
261
262
- ``dmin`` - integer
263
264
- ``dmax`` - integer
265
266
- ``L`` - defines `L`-function (default:
267
Riemann zeta function)
268
269
270
OUTPUT:
271
272
273
- ``list`` - list of pairs (d, L(s,chi_d))
274
275
276
EXAMPLES::
277
278
sage: lcalc.twist_values(0.5, -10, 10)
279
[(-8, 1.10042141), (-7, 1.14658567), (-4, 0.667691457), (-3, 0.480867558), (5, 0.231750947), (8, 0.373691713)]
280
"""
281
L = self._compute_L(L)
282
CC = sage.rings.all.ComplexField(prec)
283
Z = sage.rings.all.Integer
284
s = CC(s)
285
typ = '--twist-quadratic'
286
dmin = int(dmin)
287
dmax = int(dmax)
288
v = self('-v -x %s -y %s %s --start %s --finish %s %s'%(
289
(s.real(), s.imag(), typ, dmin, dmax, L)))
290
w = []
291
if len(v) == 0:
292
return w
293
if len(v) == 0:
294
return w
295
for a in v.split('\n'):
296
d,x,y = a.split()
297
w.append((Z(d), CC(x,y)))
298
return w
299
300
def twist_zeros(self, n, dmin, dmax, L=''):
301
r"""
302
Return first `n` real parts of nontrivial zeros for each
303
quadratic character `\chi_k` whose discriminant `d` satisfies
304
`d_{\min} \leq d \leq d_{\max}`.
305
306
INPUT:
307
308
309
- ``n`` - integer
310
311
- ``dmin`` - integer
312
313
- ``dmax`` - integer
314
315
- ``L`` - defines `L`-function (default:
316
Riemann zeta function)
317
318
319
OUTPUT:
320
321
322
- ``dict`` - keys are the discriminants `d`,
323
and values are list of corresponding zeros.
324
325
326
EXAMPLES::
327
328
sage: lcalc.twist_zeros(3, -3, 6)
329
{-3: [8.03973716, 11.2492062, 15.7046192], 5: [6.64845335, 9.83144443, 11.9588456]}
330
"""
331
L = self._compute_L(L)
332
RR = sage.rings.all.RealField(prec)
333
Z = sage.rings.all.Integer
334
typ = '--twist-quadratic'
335
n = int(n)
336
v = self('-z %s %s --start %s --finish %s %s'%(
337
(n, typ, dmin, dmax, L)))
338
w = {}
339
if len(v) == 0:
340
return w
341
for a in v.split('\n'):
342
d, x = a.split()
343
x = RR(x)
344
d = Z(d)
345
if w.has_key(d):
346
w[d].append(x)
347
else:
348
w[d] = [x]
349
return w
350
351
def analytic_rank(self, L=''):
352
r"""
353
Return the analytic rank of the `L`-function at the central
354
critical point.
355
356
INPUT:
357
358
359
- ``L`` - defines `L`-function (default:
360
Riemann zeta function)
361
362
363
OUTPUT: integer
364
365
.. note::
366
367
Of course this is not provably correct in general, since it
368
is an open problem to compute analytic ranks provably
369
correctly in general.
370
371
EXAMPLES::
372
373
sage: E = EllipticCurve('37a')
374
sage: lcalc.analytic_rank(E)
375
1
376
"""
377
L = self._compute_L(L)
378
Z = sage.rings.all.Integer
379
s = self('--rank-compute %s'%L)
380
i = s.find('equals')
381
return Z(s[i+6:])
382
383
384
385
# An instance
386
lcalc = LCalc()
387
388
389
390