Hide keyboard shortcuts

Hot-keys on this page

r m x p   toggle line displays

j k   next/prev highlighted chunk

0   (zero) top of page

1   (one) first highlighted chunk

1

2

3

4

5

6

7

8

9

10

11

12

13

14

15

16

17

18

19

20

21

22

23

24

25

26

27

28

29

30

31

32

33

34

35

36

37

38

39

40

41

42

43

44

45

46

47

48

49

50

51

52

53

54

55

56

57

58

59

60

61

62

63

64

65

66

67

68

69

70

71

72

73

74

75

76

77

78

79

80

81

82

83

84

85

86

87

88

89

90

91

92

93

94

95

96

97

98

99

100

101

102

103

104

105

106

107

108

109

110

111

112

113

114

115

116

117

118

119

120

121

122

123

124

125

126

127

128

129

130

131

132

133

134

135

136

137

138

139

140

141

142

143

144

145

146

147

148

149

150

151

152

153

154

155

156

157

158

159

160

161

162

163

164

165

166

167

168

169

170

171

172

173

174

175

176

177

178

179

180

181

182

183

184

185

186

187

188

189

190

191

192

193

194

195

196

197

198

199

200

201

202

203

204

205

206

207

208

209

210

211

212

213

214

215

216

217

218

219

220

221

222

223

224

225

226

227

228

229

230

231

232

233

234

235

236

237

238

239

240

241

242

243

244

245

246

247

248

249

250

251

252

253

254

255

256

257

258

259

260

261

262

263

264

265

266

267

268

269

270

271

272

273

274

275

276

277

278

279

280

281

282

283

284

285

286

287

288

289

290

291

292

293

294

295

""" 

Routines for computing special values of L-functions 

 

- :func:`gamma__exact` -- Exact values of the `\Gamma` function at integers and half-integers 

- :func:`zeta__exact` -- Exact values of the Riemann `\zeta` function at critical values 

- :func:`quadratic_L_function__exact` -- Exact values of the Dirichlet L-functions of quadratic characters at critical values 

- :func:`quadratic_L_function__numerical` -- Numerical values of the Dirichlet L-functions of quadratic characters in the domain of convergence 

""" 

# python3 

from __future__ import division, print_function 

 

from sage.combinat.combinat import bernoulli_polynomial 

from sage.misc.functional import denominator 

from sage.rings.all import RealField 

from sage.arith.all import kronecker_symbol, bernoulli, factorial, fundamental_discriminant 

from sage.rings.infinity import infinity 

from sage.rings.integer_ring import ZZ 

from sage.rings.polynomial.polynomial_ring_constructor import PolynomialRing 

from sage.rings.rational_field import QQ 

from sage.rings.real_mpfr import is_RealField 

from sage.symbolic.constants import pi 

from sage.symbolic.all import I 

 

# ---------------- The Gamma Function ------------------ 

 

def gamma__exact(n): 

""" 

Evaluates the exact value of the `\Gamma` function at an integer or 

half-integer argument. 

 

EXAMPLES:: 

 

sage: gamma__exact(4) 

6 

sage: gamma__exact(3) 

2 

sage: gamma__exact(2) 

1 

sage: gamma__exact(1) 

1 

 

sage: gamma__exact(1/2) 

sqrt(pi) 

sage: gamma__exact(3/2) 

1/2*sqrt(pi) 

sage: gamma__exact(5/2) 

3/4*sqrt(pi) 

sage: gamma__exact(7/2) 

15/8*sqrt(pi) 

 

sage: gamma__exact(-1/2) 

-2*sqrt(pi) 

sage: gamma__exact(-3/2) 

4/3*sqrt(pi) 

sage: gamma__exact(-5/2) 

-8/15*sqrt(pi) 

sage: gamma__exact(-7/2) 

16/105*sqrt(pi) 

 

TESTS:: 

 

sage: gamma__exact(1/3) 

Traceback (most recent call last): 

... 

TypeError: you must give an integer or half-integer argument 

""" 

from sage.all import sqrt 

n = QQ(n) 

 

if denominator(n) == 1: 

if n <= 0: 

return infinity 

if n > 0: 

return factorial(n-1) 

elif denominator(n) == 2: 

ans = QQ.one() 

while n != QQ((1, 2)): 

if n < 0: 

ans /= n 

n += 1 

elif n > 0: 

n += -1 

ans *= n 

 

ans *= sqrt(pi) 

return ans 

else: 

raise TypeError("you must give an integer or half-integer argument") 

 

# ------------- The Riemann Zeta Function -------------- 

 

def zeta__exact(n): 

r""" 

Returns the exact value of the Riemann Zeta function 

 

The argument must be a critical value, namely either positive even 

or negative odd. 

 

See for example [Iwa1972]_, p13, Special value of `\zeta(2k)` 

 

EXAMPLES: 

 

Let us test the accuracy for negative special values:: 

 

sage: RR = RealField(100) 

sage: for i in range(1,10): 

....: print("zeta({}): {}".format(1-2*i, RR(zeta__exact(1-2*i)) - zeta(RR(1-2*i)))) 

zeta(-1): 0.00000000000000000000000000000 

zeta(-3): 0.00000000000000000000000000000 

zeta(-5): 0.00000000000000000000000000000 

zeta(-7): 0.00000000000000000000000000000 

zeta(-9): 0.00000000000000000000000000000 

zeta(-11): 0.00000000000000000000000000000 

zeta(-13): 0.00000000000000000000000000000 

zeta(-15): 0.00000000000000000000000000000 

zeta(-17): 0.00000000000000000000000000000 

 

Let us test the accuracy for positive special values:: 

 

sage: all(abs(RR(zeta__exact(2*i))-zeta(RR(2*i))) < 10**(-28) for i in range(1,10)) 

True 

 

TESTS:: 

 

sage: zeta__exact(4) 

1/90*pi^4 

sage: zeta__exact(-3) 

1/120 

sage: zeta__exact(0) 

-1/2 

sage: zeta__exact(5) 

Traceback (most recent call last): 

... 

TypeError: n must be a critical value (i.e. even > 0 or odd < 0) 

 

REFERENCES: 

 

- [Iwa1972]_ 

- [IR1990]_ 

- [Was1997]_ 

""" 

if n < 0: 

return bernoulli(1-n)/(n-1) 

elif n > 1: 

if (n % 2 == 0): 

return ZZ(-1)**(n//2 + 1) * ZZ(2)**(n-1) * pi**n * bernoulli(n) / factorial(n) 

else: 

raise TypeError("n must be a critical value (i.e. even > 0 or odd < 0)") 

elif n == 1: 

return infinity 

elif n == 0: 

return QQ((-1, 2)) 

 

# ---------- Dirichlet L-functions with quadratic characters ---------- 

 

def QuadraticBernoulliNumber(k, d): 

r""" 

Compute `k`-th Bernoulli number for the primitive 

quadratic character associated to `\chi(x) = \left(\frac{d}{x}\right)`. 

 

EXAMPLES: 

 

Let us create a list of some odd negative fundamental discriminants:: 

 

sage: test_set = [d for d in range(-163, -3, 4) if is_fundamental_discriminant(d)] 

 

In general, we have `B_{1, \chi_d} = -2 h/w` for odd negative fundamental 

discriminants:: 

 

sage: all(QuadraticBernoulliNumber(1, d) == -len(BinaryQF_reduced_representatives(d)) for d in test_set) 

True 

 

REFERENCES: 

 

- [Iwa1972]_, pp 7-16. 

""" 

# Ensure the character is primitive 

d1 = fundamental_discriminant(d) 

f = abs(d1) 

 

# Make the (usual) k-th Bernoulli polynomial 

x = PolynomialRing(QQ, 'x').gen() 

bp = bernoulli_polynomial(x, k) 

 

# Make the k-th quadratic Bernoulli number 

total = sum([kronecker_symbol(d1, i) * bp(i/f) for i in range(f)]) 

total *= (f ** (k-1)) 

 

return total 

 

def quadratic_L_function__exact(n, d): 

r""" 

Returns the exact value of a quadratic twist of the Riemann Zeta function 

by `\chi_d(x) = \left(\frac{d}{x}\right)`. 

 

The input `n` must be a critical value. 

 

EXAMPLES:: 

 

sage: quadratic_L_function__exact(1, -4) 

1/4*pi 

sage: quadratic_L_function__exact(-4, -4) 

5/2 

sage: quadratic_L_function__exact(2, 1) 

1/6*pi^2 

 

TESTS:: 

 

sage: quadratic_L_function__exact(2, -4) 

Traceback (most recent call last): 

... 

TypeError: n must be a critical value (i.e. odd > 0 or even <= 0) 

 

REFERENCES: 

 

- [Iwa1972]_, pp 16-17, Special values of `L(1-n, \chi)` and `L(n, \chi)` 

- [IR1990]_ 

- [Was1997]_ 

""" 

from sage.all import SR, sqrt 

if n <= 0: 

return QuadraticBernoulliNumber(1-n,d)/(n-1) 

elif n >= 1: 

# Compute the kind of critical values (p10) 

if kronecker_symbol(fundamental_discriminant(d), -1) == 1: 

delta = 0 

else: 

delta = 1 

 

# Compute the positive special values (p17) 

if ((n - delta) % 2 == 0): 

f = abs(fundamental_discriminant(d)) 

if delta == 0: 

GS = sqrt(f) 

else: 

GS = I * sqrt(f) 

ans = SR(ZZ(-1)**(1+(n-delta)/2)) 

ans *= (2*pi/f)**n 

ans *= GS # Evaluate the Gauss sum here! =0 

ans *= QQ.one()/(2 * I**delta) 

ans *= QuadraticBernoulliNumber(n,d)/factorial(n) 

return ans 

else: 

if delta == 0: 

raise TypeError("n must be a critical value (i.e. even > 0 or odd < 0)") 

if delta == 1: 

raise TypeError("n must be a critical value (i.e. odd > 0 or even <= 0)") 

 

def quadratic_L_function__numerical(n, d, num_terms=1000): 

""" 

Evaluate the Dirichlet L-function (for quadratic character) numerically 

(in a very naive way). 

 

EXAMPLES: 

 

First, let us test several values for a given character:: 

 

sage: RR = RealField(100) 

sage: for i in range(5): 

....: print("L({}, (-4/.)): {}".format(1+2*i, RR(quadratic_L_function__exact(1+2*i, -4)) - quadratic_L_function__numerical(RR(1+2*i),-4, 10000))) 

L(1, (-4/.)): 0.000049999999500000024999996962707 

L(3, (-4/.)): 4.99999970000003...e-13 

L(5, (-4/.)): 4.99999922759382...e-21 

L(7, (-4/.)): ...e-29 

L(9, (-4/.)): ...e-29 

 

This procedure fails for negative special values, as the Dirichlet 

series does not converge here:: 

 

sage: quadratic_L_function__numerical(-3,-4, 10000) 

Traceback (most recent call last): 

... 

ValueError: the Dirichlet series does not converge here 

 

Test for several characters that the result agrees with the exact 

value, to a given accuracy :: 

 

sage: for d in range(-20,0): # long time (2s on sage.math 2014) 

....: if abs(RR(quadratic_L_function__numerical(1, d, 10000) - quadratic_L_function__exact(1, d))) > 0.001: 

....: print("Oops! We have a problem at d = {}: exact = {}, numerical = {}".format(d, RR(quadratic_L_function__exact(1, d)), RR(quadratic_L_function__numerical(1, d)))) 

""" 

# Set the correct precision if it is given (for n). 

if is_RealField(n.parent()): 

R = n.parent() 

else: 

R = RealField() 

 

if n < 0: 

raise ValueError('the Dirichlet series does not converge here') 

 

d1 = fundamental_discriminant(d) 

ans = R.zero() 

for i in range(1,num_terms): 

ans += R(kronecker_symbol(d1,i) / R(i)**n) 

return ans