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