Book a Demo!
CoCalc Logo Icon
StoreFeaturesDocsShareSupportNewsAboutPoliciesSign UpSign In
sagemath
GitHub Repository: sagemath/sagelib
Path: blob/master/sage/quadratic_forms/genera/genus.py
4077 views
1
#*****************************************************************************
2
# Copyright (C) 2007 David Kohel <[email protected]>
3
# Gabriele Nebe <[email protected]>
4
#
5
# Distributed under the terms of the GNU General Public License (GPL)
6
#
7
# http://www.gnu.org/licenses/
8
#*****************************************************************************
9
10
import sage.misc.misc as misc
11
from sage.rings.arith import LCM
12
from sage.matrix.matrix_space import MatrixSpace
13
from sage.rings.integer_ring import IntegerRing
14
from sage.rings.rational_field import RationalField
15
from sage.rings.integer import Integer
16
from sage.rings.finite_rings.constructor import FiniteField
17
18
def Genus(A):
19
"""
20
Given a nonsingular symmetric matrix A, return the genus of A.
21
22
INPUT:
23
A -- a symmetric matrix with coefficients in ZZ
24
25
OUTPUT:
26
A GenusSymbol_global_ring object, encoding the Conway-Sloane
27
genus symbol of the quadratic form whose Gram matrix is A.
28
29
EXAMPLES:
30
sage: from sage.quadratic_forms.genera.genus import GenusSymbol_global_ring
31
sage: from sage.quadratic_forms.genera.genus import Genus
32
33
sage: A = Matrix(ZZ, 2, 2, [1,1,1,2])
34
sage: Genus(A)
35
Genus of [1 1]
36
[1 2]
37
"""
38
return GenusSymbol_global_ring(A)
39
40
41
42
def LocalGenusSymbol(A,p):
43
"""
44
Given a nonsingular symmetric matrix A, return the local symbol of A at the prime p.
45
46
INPUT:
47
A -- a symmetric matrix with coefficients in ZZ
48
p -- an integer prime p > 0
49
50
OUTPUT:
51
A Genus_Symbol_p_adic_ring object, encoding the Conway-Sloane
52
genus symbol at p of the quadratic form whose Gram matrix is A.
53
54
EXAMPLES:
55
sage: from sage.quadratic_forms.genera.genus import LocalGenusSymbol
56
57
sage: A = Matrix(ZZ, 2, 2, [1,1,1,2])
58
sage: LocalGenusSymbol(A, 2)
59
Genus symbol at 2 : [[0, 2, 1, 1, 2]]
60
sage: LocalGenusSymbol(A, 3)
61
Genus symbol at 3 : [[0, 2, 1]]
62
63
sage: A = Matrix(ZZ, 2, 2, [1,0,0,2])
64
sage: LocalGenusSymbol(A, 2)
65
Genus symbol at 2 : [[0, 1, 1, 1, 1], [1, 1, 1, 1, 1]]
66
sage: LocalGenusSymbol(A, 3)
67
Genus symbol at 3 : [[0, 2, -1]]
68
"""
69
val = A.determinant().valuation(p)
70
symbol = p_adic_symbol(A, p, val = val)
71
return Genus_Symbol_p_adic_ring(p, symbol)
72
73
74
75
def is_GlobalGenus(G):
76
"""
77
Given a genus symbol G (specified by a collection of local symbols), return
78
True in G represents the genus of a global quadratic form or lattice.
79
80
INPUT:
81
G -- GenusSymbol_global_ring object
82
83
OUTPUT:
84
boolean
85
86
EXAMPLES:
87
sage: from sage.quadratic_forms.genera.genus import GenusSymbol_global_ring
88
sage: from sage.quadratic_forms.genera.genus import Genus, is_GlobalGenus
89
90
sage: A = Matrix(ZZ, 2, 2, [1,1,1,2])
91
sage: G = Genus(A)
92
sage: is_GlobalGenus(G)
93
True
94
"""
95
D = G.determinant()
96
r, s = G.signature_pair_of_matrix()
97
oddity = r - s
98
for loc in G._local_symbols:
99
p = loc._prime
100
sym = loc._symbol
101
v = sum([ s[0]*s[1] for s in sym ])
102
a = D // (p**v)
103
b = Integer(misc.prod([ s[2] for s in sym ]))
104
if p == 2:
105
if not is_2_adic_genus(sym):
106
# print "False in is_2_adic_genus(sym)"
107
return False
108
if (a*b).kronecker(p) != 1:
109
# print "False in (%s*%s).kronecker(%s)"%(a,b,p)
110
return False
111
oddity -= loc.excess()
112
else:
113
if a.kronecker(p) != b:
114
# print "False in %s.kronecker(%s) != *%s"%(a,p,b)
115
return False
116
oddity += loc.excess()
117
if oddity%8 != 0:
118
# print "False in oddity"
119
return False
120
return True
121
122
123
124
def is_2_adic_genus(genus_symbol_quintuple_list):
125
"""
126
Given a 2-adic local symbol (as the underlying list of quintuples)
127
check whether it is the 2-adic symbol of a 2-adic form.
128
129
INPUT:
130
genus_symbol_quintuple_list -- a quintuple of integers (with certain
131
restrictions).
132
133
OUTPUT:
134
boolean
135
136
EXAMPLES:
137
sage: from sage.quadratic_forms.genera.genus import LocalGenusSymbol
138
139
sage: A = Matrix(ZZ, 2, 2, [1,1,1,2])
140
sage: G2 = LocalGenusSymbol(A, 2)
141
sage: is_2_adic_genus(G2.symbol_tuple_list())
142
True
143
144
sage: A = Matrix(ZZ, 2, 2, [1,1,1,2])
145
sage: G3 = LocalGenusSymbol(A, 3)
146
sage: is_2_adic_genus(G3.symbol_tuple_list()) ## This raises an error
147
Traceback (most recent call last):
148
...
149
TypeError: The genus symbols are not quintuples, so it's not a genus symbol for the prime p=2.
150
151
sage: A = Matrix(ZZ, 2, 2, [1,0,0,2])
152
sage: G2 = LocalGenusSymbol(A, 2)
153
sage: is_2_adic_genus(G2.symbol_tuple_list())
154
True
155
"""
156
## TO DO: Add explicit checking for the prime p here to ensure it's p=2... not just the quintuple checking below
157
158
for s in genus_symbol_quintuple_list:
159
160
## Check that we have a quintuple (i.e. that p=2 and not p >2)
161
if len(s) != 5:
162
raise TypeError, "The genus symbols are not quintuples, so it's not a genus symbol for the prime p=2."
163
164
## Check the Conway-Sloane conditions
165
if s[1] == 1:
166
if s[3] == 0 or s[2] != s[4]:
167
return False
168
if s[1] == 2 and s[3] == 1:
169
if s[2] in (1,-1):
170
if not s[4] in (0,2,6):
171
return False
172
if s[2] in (3,-3):
173
if not s[4] in (2,4,6):
174
return False
175
if (s[1] - s[4])% 2 == 1:
176
return False
177
if s[3] == 0 and s[4] != 0:
178
return False
179
return True
180
181
182
183
def canonical_2_adic_compartments(genus_symbol_quintuple_list):
184
"""
185
Given a 2-adic local symbol (as the underlying list of quintuples)
186
this returns a list of lists of indices of the
187
genus_symbol_quintuple_list which are in the same compartment. A
188
compartment is defined to be a maximal interval of Jordan
189
components all (scaled) of type I (i.e. odd).
190
191
INPUT:
192
genus_symbol_quintuple_list -- a quintuple of integers (with certain
193
restrictions).
194
195
OUTPUT:
196
a list of lists of integers.
197
198
EXAMPLES:
199
sage: from sage.quadratic_forms.genera.genus import LocalGenusSymbol
200
sage: from sage.quadratic_forms.genera.genus import canonical_2_adic_compartments
201
202
sage: A = Matrix(ZZ, 2, 2, [1,1,1,2])
203
sage: G2 = LocalGenusSymbol(A, 2); G2
204
Genus symbol at 2 : [[0, 2, 1, 1, 2]]
205
sage: canonical_2_adic_compartments(G2.symbol_tuple_list())
206
[[0]]
207
208
sage: A = Matrix(ZZ, 2, 2, [1,0,0,2])
209
sage: G2 = LocalGenusSymbol(A, 2); G2
210
Genus symbol at 2 : [[0, 1, 1, 1, 1], [1, 1, 1, 1, 1]]
211
sage: canonical_2_adic_compartments(G2.symbol_tuple_list())
212
[[0, 1]]
213
214
sage: A = DiagonalQuadraticForm(ZZ, [1,2,3,4]).Hessian_matrix()
215
sage: G2 = LocalGenusSymbol(A, 2); G2
216
Genus symbol at 2 : [[1, 2, 3, 1, 4], [2, 1, 1, 1, 1], [3, 1, 1, 1, 1]]
217
sage: canonical_2_adic_compartments(G2.symbol_tuple_list())
218
[[0, 1, 2]]
219
220
sage: A = Matrix(ZZ, 2, 2, [2,1,1,2])
221
sage: G2 = LocalGenusSymbol(A, 2); G2
222
Genus symbol at 2 : [[0, 2, 3, 0, 0]]
223
sage: canonical_2_adic_compartments(G2.symbol_tuple_list()) ## No compartments here!
224
[]
225
226
NOTES:
227
See Conway-Sloane 3rd edition, pp. 381-382 for definitions and examples.
228
"""
229
symbol = genus_symbol_quintuple_list
230
compartments = []
231
i = 0
232
r = len(symbol)
233
while i < r:
234
s = symbol[i]
235
if s[3] == 1:
236
v = s[0]
237
c = []
238
while i < r and symbol[i][3] == 1 and symbol[i][0] == v:
239
c.append(i)
240
i += 1
241
v += 1
242
compartments.append(c)
243
else:
244
i += 1
245
return compartments
246
247
248
249
250
251
def canonical_2_adic_trains(genus_symbol_quintuple_list, compartments=None):
252
"""
253
Given a 2-adic local symbol (as the underlying list of quintuples)
254
this returns a list of lists of indices of the
255
genus_symbol_quintuple_list which are in the same train. A train
256
is defined to be a maximal interval of Jordan components so that
257
at least one of each adjacent pair (allowing zero-dimensional
258
Jordan components) is (scaled) of type I (i.e. odd).
259
260
INPUT:
261
genus_symbol_quintuple_list -- a quintuple of integers (with certain
262
restrictions).
263
compartments -- a list of lists of distinct integers (optional)
264
265
OUTPUT:
266
a list of lists of distinct integers.
267
268
EXAMPLES:
269
sage: from sage.quadratic_forms.genera.genus import LocalGenusSymbol
270
sage: from sage.quadratic_forms.genera.genus import canonical_2_adic_compartments
271
sage: from sage.quadratic_forms.genera.genus import canonical_2_adic_trains
272
273
sage: A = Matrix(ZZ, 2, 2, [1,1,1,2])
274
sage: G2 = LocalGenusSymbol(A, 2); G2
275
Genus symbol at 2 : [[0, 2, 1, 1, 2]]
276
sage: c = canonical_2_adic_compartments(G2.symbol_tuple_list()); c
277
[[0]]
278
sage: canonical_2_adic_trains(G2.symbol_tuple_list(), c)
279
[[0]]
280
281
sage: A = Matrix(ZZ, 2, 2, [1,0,0,2])
282
sage: G2 = LocalGenusSymbol(A, 2); G2
283
Genus symbol at 2 : [[0, 1, 1, 1, 1], [1, 1, 1, 1, 1]]
284
sage: c = canonical_2_adic_compartments(G2.symbol_tuple_list()); c
285
[[0, 1]]
286
sage: canonical_2_adic_trains(G2.symbol_tuple_list(), c)
287
[[0, 1]]
288
sage: canonical_2_adic_trains(G2.symbol_tuple_list())
289
[[0, 1]]
290
291
sage: A = DiagonalQuadraticForm(ZZ, [1,2,3,4]).Hessian_matrix()
292
sage: G2 = LocalGenusSymbol(A, 2); G2
293
Genus symbol at 2 : [[1, 2, 3, 1, 4], [2, 1, 1, 1, 1], [3, 1, 1, 1, 1]]
294
sage: c = canonical_2_adic_compartments(G2.symbol_tuple_list()); c
295
[[0, 1, 2]]
296
sage: canonical_2_adic_trains(G2.symbol_tuple_list())
297
[[0, 1, 2]]
298
299
sage: A = Matrix(ZZ, 2, 2, [2,1,1,2])
300
sage: G2 = LocalGenusSymbol(A, 2); G2
301
Genus symbol at 2 : [[0, 2, 3, 0, 0]]
302
sage: c = canonical_2_adic_compartments(G2.symbol_tuple_list()); c ## No compartments here!
303
[]
304
sage: canonical_2_adic_trains(G2.symbol_tuple_list())
305
[]
306
307
NOTES:
308
See Conway-Sloane 3rd edition, pp. 381-382 for definitions and examples.
309
310
TO DO:
311
- Add a non-trivial example in the doctest here!
312
"""
313
## Recompute compartments if none are passed.
314
if compartments == None:
315
compartments = canonical_2_adic_compartments(genus_symbol_quintuple_list)
316
317
symbol = genus_symbol_quintuple_list
318
trains = []
319
i = 0
320
while i < len(compartments):
321
flag = True
322
train = [ ]
323
while flag:
324
ci = compartments[i]
325
j = ci[0]
326
if j == 0 or symbol[j-1][0] != symbol[j][0] - 1:
327
train += ci
328
else:
329
train += [j-1] + ci
330
act = ci[len(ci)-1]+1
331
if i+1 < len(compartments):
332
ci_plus = compartments[i+1]
333
else:
334
if act < len(symbol) and symbol[act][0] == symbol[act-1][0] +1:
335
train += [act]
336
flag = False
337
if flag and symbol[ci[len(ci)-1]][0]+2 != symbol[ci_plus[0]][0]:
338
if act != ci_plus[0] and symbol[act][0] == symbol[act-1][0] +1:
339
train += [act]
340
flag = False
341
i += 1
342
trains.append(train)
343
return trains
344
345
346
347
348
def canonical_2_adic_reduction(genus_symbol_quintuple_list):
349
"""
350
Given a 2-adic local symbol (as the underlying list of quintuples)
351
this returns a canonical 2-adic symbol (again as a raw list of
352
quintuples of integers) which has at most one minus sign per train
353
and this sign appears on the smallest dimensional Jordan component
354
in each train. This results from applying the "sign-walking" and
355
"oddity fusion" equivalences.
356
357
INPUT:
358
genus_symbol_quintuple_list -- a quintuple of integers (with certain
359
restrictions).
360
compartments -- a list of lists of distinct integers (optional)
361
362
OUTPUT:
363
a list of lists of distinct integers.
364
365
EXAMPLES:
366
sage: from sage.quadratic_forms.genera.genus import LocalGenusSymbol
367
sage: from sage.quadratic_forms.genera.genus import canonical_2_adic_reduction
368
369
sage: A = Matrix(ZZ, 2, 2, [1,1,1,2])
370
sage: G2 = LocalGenusSymbol(A, 2); G2
371
Genus symbol at 2 : [[0, 2, 1, 1, 2]]
372
sage: canonical_2_adic_reduction(G2.symbol_tuple_list())
373
[[0, 2, 1, 1, 2]]
374
375
sage: A = Matrix(ZZ, 2, 2, [1,0,0,2])
376
sage: G2 = LocalGenusSymbol(A, 2); G2
377
Genus symbol at 2 : [[0, 1, 1, 1, 1], [1, 1, 1, 1, 1]]
378
sage: canonical_2_adic_reduction(G2.symbol_tuple_list()) ## Oddity fusion occurred here!
379
[[0, 1, 1, 1, 2], [1, 1, 1, 1, 0]]
380
381
sage: A = DiagonalQuadraticForm(ZZ, [1,2,3,4]).Hessian_matrix()
382
sage: G2 = LocalGenusSymbol(A, 2); G2
383
Genus symbol at 2 : [[1, 2, 3, 1, 4], [2, 1, 1, 1, 1], [3, 1, 1, 1, 1]]
384
sage: canonical_2_adic_reduction(G2.symbol_tuple_list()) ## Oddity fusion occurred here!
385
[[1, 2, -1, 1, 6], [2, 1, 1, 1, 0], [3, 1, 1, 1, 0]]
386
387
sage: A = Matrix(ZZ, 2, 2, [2,1,1,2])
388
sage: G2 = LocalGenusSymbol(A, 2); G2
389
Genus symbol at 2 : [[0, 2, 3, 0, 0]]
390
sage: canonical_2_adic_reduction(G2.symbol_tuple_list())
391
[[0, 2, -1, 0, 0]]
392
393
NOTES:
394
See Conway-Sloane 3rd edition, pp. 381-382 for definitions and examples.
395
396
TO DO:
397
- Add an example where sign walking occurs!
398
"""
399
canonical_symbol = genus_symbol_quintuple_list
400
# Canonical determinants:
401
for i in range(len(genus_symbol_quintuple_list)):
402
d = genus_symbol_quintuple_list[i][2]
403
if d in (1,7):
404
canonical_symbol[i][2] = 1
405
else:
406
canonical_symbol[i][2] = -1
407
# Oddity fusion:
408
compartments = canonical_2_adic_compartments(genus_symbol_quintuple_list)
409
for compart in compartments:
410
oddity = sum([ genus_symbol_quintuple_list[i][4] for i in compart ]) % 8
411
for i in compart:
412
genus_symbol_quintuple_list[i][4] = 0
413
genus_symbol_quintuple_list[compart[0]][4] = oddity
414
#print "End oddity fusion:", canonical_symbol
415
# Sign walking:
416
trains = canonical_2_adic_trains(genus_symbol_quintuple_list, compartments)
417
for train in trains:
418
t = len(train)
419
for i in range(t-1):
420
t1 = train[t-i-1]
421
if canonical_symbol[t1][2] == -1:
422
canonical_symbol[t1][2] = 1
423
canonical_symbol[t1-1][2] *= -1
424
for compart in compartments:
425
if t1-1 in compart or t1 in compart:
426
o = canonical_symbol[compart[0]][4]
427
canonical_symbol[compart[0]][4] = (o+4) % 8
428
#print "End sign walking:", canonical_symbol
429
return canonical_symbol
430
431
432
433
434
435
def basis_complement(B):
436
"""
437
Given an echelonized basis matrix (over a field), calculate a
438
matrix whose rows form a basis complement (to the rows of B).
439
440
INPUT:
441
B -- matrix over a field in row echelon form
442
443
OUTPUT:
444
a rectangular matrix over a field
445
446
EXAMPLES:
447
sage: from sage.quadratic_forms.genera.genus import basis_complement
448
449
sage: A = Matrix(ZZ, 2, 2, [1,1,1,1])
450
sage: B = A.kernel().echelonized_basis_matrix(); B
451
[ 1 -1]
452
sage: basis_complement(B)
453
[0 1]
454
"""
455
F = B.parent().base_ring()
456
m = B.nrows()
457
n = B.ncols()
458
C = MatrixSpace(F,n-m,n,sparse=True)(0)
459
k = 0
460
l = 0
461
for i in range(m):
462
for j in range(k,n):
463
if B[i,j] == 0:
464
C[l,j] = 1
465
l += 1
466
else:
467
k = j+1
468
break
469
for j in range(k,n):
470
C[l+j-k,j] = 1
471
return C
472
473
474
475
def signature_pair_of_matrix(A):
476
"""
477
Computes the signature pair (p, n) of a non-degenerate symmetric
478
matrix, where
479
p = number of positive eigenvalues of A
480
n = number of negative eigenvalues of A
481
482
INPUT:
483
A -- symmetric matrix (assumed to be non-degenerate)
484
485
OUTPUT:
486
a pair (tuple) of integers.
487
488
EXAMPLES:
489
sage: from sage.quadratic_forms.genera.genus import signature_pair_of_matrix
490
491
sage: A = Matrix(ZZ, 2, 2, [-1,0,0,3])
492
sage: signature_pair_of_matrix(A)
493
(1, 1)
494
495
sage: A = Matrix(ZZ, 2, 2, [-1,1,1,7])
496
sage: signature_pair_of_matrix(A)
497
(1, 1)
498
499
sage: A = Matrix(ZZ, 2, 2, [3,1,1,7])
500
sage: signature_pair_of_matrix(A)
501
(2, 0)
502
503
sage: A = Matrix(ZZ, 2, 2, [-3,1,1,-11])
504
sage: signature_pair_of_matrix(A)
505
(0, 2)
506
507
508
sage: A = Matrix(ZZ, 2, 2, [1,1,1,1])
509
sage: signature_pair_of_matrix(A) ## Raises an error -- degenerate matrix
510
Traceback (most recent call last):
511
...
512
TypeError: A is assumed to be non-degenerate, but it's det = 0.
513
514
"""
515
from sage.quadratic_forms.quadratic_form import QuadraticForm
516
s_vec = QuadraticForm(A.base_extend(A.base_ring().fraction_field())).signature_vector()
517
518
## Check that the matrix is non-degenerate (i.e. no zero eigenvalues)
519
if s_vec[2] != 0:
520
raise TypeError, "A is assumed to be non-degenerate, but it's det = 0."
521
522
## Return the pair (p,n)
523
return s_vec[:2]
524
525
526
527
def p_adic_symbol(A, p, val):
528
"""
529
Given a symmetric matrix A and prime p, return the genus symbol at p.
530
531
val = valuation of the maximal elementary divisor of A
532
needed to obtain enough precision
533
calculation is modulo p to the val+3
534
TODO: Some description of the definition of the genus symbol.
535
536
INPUT:
537
A -- symmetric matrix with integer coefficients
538
p -- prime number > 0
539
val -- integer >= 0
540
541
OUTPUT:
542
a list of lists of integers
543
544
EXAMPLES:
545
sage: from sage.quadratic_forms.genera.genus import p_adic_symbol
546
547
sage: A = DiagonalQuadraticForm(ZZ, [1,2,3,4]).Hessian_matrix()
548
sage: p_adic_symbol(A, 2, 2)
549
[[1, 2, 3, 1, 4], [2, 1, 1, 1, 1], [3, 1, 1, 1, 1]]
550
551
sage: p_adic_symbol(A, 3, 1)
552
[[0, 3, 1], [1, 1, -1]]
553
554
"""
555
if p % 2 == 0:
556
return two_adic_symbol(A, val)
557
m0 = min([ c.valuation(p) for c in A.list() ])
558
q = p**m0
559
n = A.nrows()
560
A = MatrixSpace(IntegerRing(),n,n)([ c // q for c in A.list() ])
561
A_p = MatrixSpace(FiniteField(p),n,n)(A)
562
B_p = A_p.kernel().echelonized_basis_matrix()
563
if B_p.nrows() == 0:
564
e0 = Integer(A_p.det()).kronecker(p)
565
n0 = A.nrows()
566
return [ [m0,n0,e0] ]
567
else:
568
C_p = basis_complement(B_p)
569
e0 = Integer((C_p*A_p*C_p.transpose()).det()).kronecker(p)
570
n0 = C_p.nrows()
571
sym = [ [0,n0,e0] ]
572
r = B_p.nrows()
573
B = MatrixSpace(IntegerRing(),r,n)(B_p)
574
C = MatrixSpace(IntegerRing(),n-r,n)(C_p)
575
# Construct the blocks for the Jordan decomposition [F,X;X,A_new]
576
F = MatrixSpace(RationalField(),n-r,n-r)(C*A*C.transpose())
577
U = F**-1
578
d = LCM([ c.denominator() for c in U.list() ])
579
R = IntegerRing().quotient_ring(Integer(p)**(val+3))
580
u = R(d)**-1
581
MatR = MatrixSpace(R,n-r,n-r)
582
MatZ = MatrixSpace(IntegerRing(),n-r,n-r)
583
U = MatZ(MatR(MatZ(U*d))*u)
584
# X = C*A*B.transpose()
585
# A = B*A*B.transpose() - X.transpose()*U*X
586
X = C*A
587
A = B*(A - X.transpose()*U*X)*B.transpose()
588
return [ [s[0]+m0] + s[1:] for s in sym + p_adic_symbol(A, p, val) ]
589
590
591
592
def is_even_matrix(A):
593
"""
594
Determines if the integral symmetric matrix A is even
595
(i.e. represents only even numbers). If not, then it returns the
596
index of an odd diagonal entry. If it is even, then we return the
597
index -1.
598
599
INPUT:
600
A -- symmetric integer matrix
601
602
OUTPUT:
603
a pair of the form (boolean, integer)
604
605
EXAMPLES:
606
sage: from sage.quadratic_forms.genera.genus import is_even_matrix
607
608
sage: A = Matrix(ZZ, 2, 2, [1,1,1,1])
609
sage: is_even_matrix(A)
610
(False, 0)
611
612
sage: A = Matrix(ZZ, 2, 2, [2,1,1,2])
613
sage: is_even_matrix(A)
614
(True, -1)
615
"""
616
for i in range(A.nrows()):
617
if A[i,i]%2 == 1:
618
return False, i
619
return True, -1
620
621
622
623
def split_odd(A):
624
"""
625
Given a non-degenerate Gram matrix A (mod 8), return a splitting [u] + B
626
such that u is odd and B is not even.
627
628
INPUT:
629
A -- an odd symmetric matrix with integer coefficients (which
630
admits a splitting as above).
631
632
OUTPUT:
633
a pair (u, B) consisting of an odd integer u and an odd
634
integral symmetric matrix B.
635
636
EXAMPLES:
637
sage: from sage.quadratic_forms.genera.genus import is_even_matrix
638
sage: from sage.quadratic_forms.genera.genus import split_odd
639
640
sage: A = Matrix(ZZ, 2, 2, [1,2,2,3])
641
sage: is_even_matrix(A)
642
(False, 0)
643
sage: split_odd(A)
644
(1, [-1])
645
646
sage: A = Matrix(ZZ, 2, 2, [1,2,2,5])
647
sage: split_odd(A)
648
(1, [1])
649
650
sage: A = Matrix(ZZ, 2, 2, [1,1,1,1])
651
sage: is_even_matrix(A)
652
(False, 0)
653
sage: split_odd(A) ## This fails because no such splitting exists. =(
654
Traceback (most recent call last):
655
...
656
RuntimeError: The matrix A does not admit a non-even splitting.
657
658
sage: A = Matrix(ZZ, 2, 2, [1,2,2,6])
659
sage: split_odd(A) ## This fails because no such splitting exists. =(
660
Traceback (most recent call last):
661
...
662
RuntimeError: The matrix A does not admit a non-even splitting.
663
664
"""
665
n0 = A.nrows()
666
if n0 == 1:
667
return A[0,0], MatrixSpace(IntegerRing(),0,A.ncols())([])
668
even, i = is_even_matrix(A)
669
R = A.parent().base_ring()
670
C = MatrixSpace(R,n0-1,n0)(0)
671
u = A[i,i]
672
for j in range(n0-1):
673
if j < i:
674
C[j,j] = 1
675
C[j,i] = -A[j,i]*u
676
else:
677
C[j,j+1] = 1
678
C[j,i] = -A[j+1,i]*u
679
B = C*A*C.transpose()
680
even, j = is_even_matrix(B)
681
if even:
682
I = A.parent()(1)
683
# TODO: we could manually (re)construct the kernel here...
684
if i == 0:
685
I[1,0] = 1 - A[1,0]*u
686
i = 1
687
else:
688
I[0,i] = 1 - A[0,i]*u
689
i = 0
690
A = I*A*I.transpose()
691
u = A[i,i]
692
C = MatrixSpace(R,n0-1,n0)(0)
693
for j in range(n0-1):
694
if j < i:
695
C[j,j] = 1
696
C[j,i] = -A[j,i]*u
697
else:
698
C[j,j+1] = 1
699
C[j,i] = -A[j+1,i]*u
700
B = C*A*C.transpose()
701
even, j = is_even_matrix(B)
702
if even:
703
print "B:"
704
print B
705
raise RuntimeError, "The matrix A does not admit a non-even splitting."
706
return u, B
707
708
709
710
def trace_diag_mod_8(A):
711
"""
712
Return the trace of the diagonalised form of A of an integral
713
symmetric matrix which is diagonalizable mod 8. (Note that since
714
the Jordan decomposition into blocks of size <= 2 is not unique
715
here, this is not the same as saying that A is always diagonal in
716
any 2-adic Jordan decomposition!)
717
718
INPUT:
719
A -- symmetric matrix with coefficients in Z which is odd in
720
Z/2Z and has determinant not divisible by 8.
721
722
OUTPUT:
723
an integer
724
725
EXAMPLES:
726
sage: from sage.quadratic_forms.genera.genus import is_even_matrix
727
sage: from sage.quadratic_forms.genera.genus import split_odd
728
sage: from sage.quadratic_forms.genera.genus import trace_diag_mod_8
729
730
sage: A = Matrix(ZZ, 2, 2, [1,2,2,3])
731
sage: is_even_matrix(A)
732
(False, 0)
733
sage: split_odd(A)
734
(1, [-1])
735
sage: trace_diag_mod_8(A)
736
0
737
738
sage: A = Matrix(ZZ, 2, 2, [1,2,2,5])
739
sage: split_odd(A)
740
(1, [1])
741
sage: trace_diag_mod_8(A)
742
2
743
"""
744
tr = 0
745
while A.nrows() > 0:
746
u, A = split_odd(A)
747
tr += u
748
return IntegerRing()(tr)
749
750
751
752
def two_adic_symbol(A, val):
753
"""
754
Given a symmetric matrix A and prime p, return the genus symbol at p.
755
756
val = valuation of maximal 2-elementary divisor
757
758
The genus symbol of a component 2^m*f is of the form (m,n,s,d[,o]),
759
where
760
m = valuation of the component
761
n = dimension of f
762
d = det(f) in {1,3,5,7}
763
s = 0 (or 1) if even (or odd)
764
o = oddity of f (= 0 if s = 0) in Z/8Z
765
766
INPUT:
767
A -- symmetric matrix with integer coefficients
768
val -- integer >=0
769
770
OUTPUT:
771
a list of lists of integers (representing a Conway-Sloane 2-adic symbol)
772
773
EXAMPLES:
774
sage: from sage.quadratic_forms.genera.genus import two_adic_symbol
775
776
sage: A = diagonal_matrix(ZZ, [1,2,3,4])
777
sage: two_adic_symbol(A, 2)
778
[[0, 2, 3, 1, 4], [1, 1, 1, 1, 1], [2, 1, 1, 1, 1]]
779
780
"""
781
m0 = min([ c.valuation(2) for c in A.list() ])
782
q = 2**m0
783
A = A.parent()([ c // q for c in A.list() ])
784
ZZ = IntegerRing()
785
n = A.nrows()
786
A_2 = MatrixSpace(FiniteField(2),n,n)(A)
787
K_2 = A_2.kernel()
788
R_8 = ZZ.quotient_ring(Integer(8))
789
790
## Deal with the matrix being non-degenerate mod 2.
791
if K_2.dimension() == 0:
792
A_8 = MatrixSpace(R_8,n)(A)
793
n0 = A.nrows()
794
# d0 = ZZ(A_8.determinant()) # no determinant over Z/8Z
795
d0 = ZZ(R_8(MatrixSpace(ZZ,n)(A_8).determinant()))
796
if d0 == 0: ## SANITY CHECK: The mod 8 determinant shouldn't be zero.
797
print "A:"
798
print A
799
assert False
800
even, i = is_even_matrix(A_2) ## Determine whether the matrix is even or odd.
801
if even:
802
return [ [m0,n0,d0,0,0] ]
803
else:
804
tr8 = trace_diag_mod_8(A_8) ## Here we already know that A_8 is odd and diagonalizable mod 8.
805
return [ [m0,n0,d0,1,tr8] ]
806
807
## Deal with the matrix being degenerate mod 2.
808
else:
809
B_2 = K_2.echelonized_basis_matrix()
810
C_2 = basis_complement(B_2)
811
n0 = C_2.nrows()
812
C = MatrixSpace(ZZ,n0,n)(C_2)
813
A_new = C*A*C.transpose()
814
# compute oddity modulo 8:
815
A_8 = MatrixSpace(R_8,n0,n0)(A_new)
816
# d0 = A_8.det() # no determinant over Z/8Z
817
d0 = ZZ(R_8(MatrixSpace(ZZ,n0,n0)(A_8).determinant()))
818
if d0 == 0:
819
print "A:"
820
print A_new
821
assert False
822
even, i = is_even_matrix(A_new)
823
if even:
824
sym = [ [0,n0,d0,0,0] ]
825
else:
826
tr8 = trace_diag_mod_8(A_8)
827
sym = [ [0,n0,d0,1,tr8] ]
828
r = B_2.nrows()
829
B = MatrixSpace(ZZ,r,n)(B_2)
830
C = MatrixSpace(IntegerRing(),n-r,n)(C_2)
831
F = MatrixSpace(RationalField(),n-r,n-r)(C*A*C.transpose())
832
U = F**-1
833
d = LCM([ c.denominator() for c in U.list() ])
834
R = IntegerRing().quotient_ring(Integer(2)**(val+3))
835
u = R(d)**-1
836
MatR = MatrixSpace(R,n-r,n-r)
837
MatZ = MatrixSpace(IntegerRing(),n-r,n-r)
838
U = MatZ(MatR(MatZ(U*d))*u)
839
X = C*A
840
A = B*(A - X.transpose()*U*X)*B.transpose()
841
return [ [s[0]+m0] + s[1:] for s in sym + two_adic_symbol(A, val) ]
842
843
844
845
846
847
## Removed because it was unused and undocumented!
848
#
849
#def is_trivial_symbol(p, sym):
850
# """
851
# """
852
# if len(sym) != 1:
853
# return False
854
# if sym[0] != 0 or sym[2] != 1:
855
# return False
856
# if p != 2:
857
# return True
858
# return sym[3] == 1 and sym[1] % 8 == sym[4]
859
860
861
862
863
864
865
866
867
class Genus_Symbol_p_adic_ring(object):
868
"""
869
Local genus symbol over a p-adic ring.
870
"""
871
def __init__(self, prime, symbol, check = True):
872
"""
873
Create the local genus symbol of given prime and local invariants.
874
875
The genus symbol of a component p^m*A for odd prime = p is of the
876
form (m,n,d), where
877
878
m = valuation of the component
879
n = rank of A
880
d = det(A) in {1,u} for normalized quadratic non-residue u.
881
882
The genus symbol of a component 2^m*A is of the form (m,n,s,d,o),
883
where
884
885
m = valuation of the component
886
n = rank of A
887
d = det(A) in {1,3,5,7}
888
s = 0 (or 1) if even (or odd)
889
o = oddity of A (= 0 if s = 0) in Z/8Z
890
= the trace of the diagonalization of A
891
892
The genus symbol is a list of such symbols (ordered by m) for each
893
of the Jordan blocks A_1,...,A_t.
894
895
Reference: Conway and Sloane 3rd edition, Chapter 15, Section 7.
896
897
898
WARNING/NOTE: This normalization seems non-standard, and we
899
should review this entire class to make sure that we have our
900
doubling conventions straight throughout! This is especially
901
noticeable in the determinant and excess methods!!
902
903
904
INPUT:
905
prime -- a prime integer > 0
906
symbol -- the list of invariants for Jordan blocks
907
A_t,...,A_t given as a list of lists of integers
908
909
OUTPUT:
910
None
911
912
EXAMPLES:
913
sage: from sage.quadratic_forms.genera.genus import p_adic_symbol
914
sage: from sage.quadratic_forms.genera.genus import Genus_Symbol_p_adic_ring
915
916
sage: A = diagonal_matrix(ZZ, [1,2,3,4])
917
sage: p = 2
918
sage: s2 = p_adic_symbol(A, p, 2); s2
919
[[0, 2, 3, 1, 4], [1, 1, 1, 1, 1], [2, 1, 1, 1, 1]]
920
sage: G = Genus_Symbol_p_adic_ring(p,s2);G
921
Genus symbol at 2 : [[0, 2, 3, 1, 4], [1, 1, 1, 1, 1], [2, 1, 1, 1, 1]]
922
sage: G == loads(dumps(G))
923
True
924
925
sage: A = diagonal_matrix(ZZ, [1,2,3,4])
926
sage: p = 3
927
sage: s3 = p_adic_symbol(A, p, 1); s3
928
[[0, 3, -1], [1, 1, 1]]
929
sage: G = Genus_Symbol_p_adic_ring(p,s3);G
930
Genus symbol at 3 : [[0, 3, -1], [1, 1, 1]]
931
sage: G == loads(dumps(G))
932
True
933
934
935
"""
936
if check:
937
pass
938
self._prime = prime
939
self._symbol = symbol
940
self._canonical_symbol = None
941
942
def __repr__(self):
943
"""
944
Gives a string representation for the p-adic genus symbol
945
946
INPUT:
947
None
948
949
OUTPUT:
950
a string
951
952
EXAMPLES:
953
sage: from sage.quadratic_forms.genera.genus import two_adic_symbol
954
sage: from sage.quadratic_forms.genera.genus import Genus_Symbol_p_adic_ring
955
956
sage: A = diagonal_matrix(ZZ, [1,2,3,4])
957
sage: s2 = two_adic_symbol(A, 2); s2
958
[[0, 2, 3, 1, 4], [1, 1, 1, 1, 1], [2, 1, 1, 1, 1]]
959
sage: G = Genus_Symbol_p_adic_ring(2, s2)
960
sage: G.__repr__()
961
'Genus symbol at 2 : [[0, 2, 3, 1, 4], [1, 1, 1, 1, 1], [2, 1, 1, 1, 1]]'
962
963
"""
964
return "Genus symbol at %s : %s"%(self._prime, self._symbol)
965
966
967
def __eq__(self, other):
968
"""
969
Determines if two genus symbols are equal (not just equivalent!).
970
971
INPUT:
972
a Genus_Symbol_p_adic_ring object
973
974
OUTPUT:
975
boolean
976
977
EXAMPLES:
978
sage: from sage.quadratic_forms.genera.genus import p_adic_symbol
979
sage: from sage.quadratic_forms.genera.genus import Genus_Symbol_p_adic_ring
980
981
sage: A = diagonal_matrix(ZZ, [1,2,3,4])
982
sage: p = 2
983
sage: G2 = Genus_Symbol_p_adic_ring(p, p_adic_symbol(A, p, 2))
984
sage: p = 3
985
sage: G3 = Genus_Symbol_p_adic_ring(p, p_adic_symbol(A, p, 1))
986
987
sage: G2 == G3
988
False
989
sage: G3 == G2
990
False
991
sage: G2 == G2
992
True
993
sage: G3 == G3
994
True
995
996
"""
997
p = self._prime
998
if p != other._prime:
999
return False
1000
return self.canonical_symbol() == other.canonical_symbol()
1001
1002
1003
def __ne__(self, other):
1004
"""
1005
Determines if two genus symbols are unequal (not just inequivalent!).
1006
1007
INPUT:
1008
a Genus_Symbol_p_adic_ring object
1009
1010
OUTPUT:
1011
boolean
1012
1013
sage: from sage.quadratic_forms.genera.genus import p_adic_symbol
1014
sage: from sage.quadratic_forms.genera.genus import Genus_Symbol_p_adic_ring
1015
1016
sage: A = diagonal_matrix(ZZ, [1,2,3,4])
1017
sage: p = 2
1018
sage: G2 = Genus_Symbol_p_adic_ring(p, p_adic_symbol(A, p, 2))
1019
sage: p = 3
1020
sage: G3 = Genus_Symbol_p_adic_ring(p, p_adic_symbol(A, p, 1))
1021
1022
sage: G2 != G3
1023
True
1024
sage: G3 != G2
1025
True
1026
sage: G2 != G2
1027
False
1028
sage: G3 != G3
1029
False
1030
1031
"""
1032
return not self.__eq__(other)
1033
1034
1035
## Added these two methods to make this class iterable...
1036
#def __getitem__(self, i):
1037
# return self._symbol[i]
1038
#
1039
#def len(self):
1040
# return len(self._symbol)
1041
## ------------------------------------------------------
1042
1043
1044
def canonical_symbol(self):
1045
"""
1046
Return (and cache) the canonical p-adic genus symbol. This is
1047
only really affects the 2-adic symbol, since when p > 2 the
1048
symbol is already canonical.
1049
1050
INPUT:
1051
None
1052
1053
OUTPUT:
1054
a list of lists of integers
1055
1056
EXAMPLES:
1057
sage: from sage.quadratic_forms.genera.genus import p_adic_symbol
1058
sage: from sage.quadratic_forms.genera.genus import Genus_Symbol_p_adic_ring
1059
1060
sage: A = Matrix(ZZ, 2, 2, [1,1,1,2])
1061
sage: p = 2
1062
sage: G2 = Genus_Symbol_p_adic_ring(p, p_adic_symbol(A, p, 2)); G2
1063
Genus symbol at 2 : [[0, 2, 1, 1, 2]]
1064
sage: G2.canonical_symbol()
1065
[[0, 2, 1, 1, 2]]
1066
1067
sage: A = Matrix(ZZ, 2, 2, [1,0,0,2])
1068
sage: p = 2
1069
sage: G2 = Genus_Symbol_p_adic_ring(p, p_adic_symbol(A, p, 2)); G2
1070
Genus symbol at 2 : [[0, 1, 1, 1, 1], [1, 1, 1, 1, 1]]
1071
sage: G2.canonical_symbol() ## Oddity fusion occurred here!
1072
[[0, 1, 1, 1, 2], [1, 1, 1, 1, 0]]
1073
1074
sage: A = DiagonalQuadraticForm(ZZ, [1,2,3,4]).Hessian_matrix()
1075
sage: p = 2
1076
sage: G2 = Genus_Symbol_p_adic_ring(p, p_adic_symbol(A, p, 2)); G2
1077
Genus symbol at 2 : [[1, 2, 3, 1, 4], [2, 1, 1, 1, 1], [3, 1, 1, 1, 1]]
1078
sage: G2.canonical_symbol() ## Oddity fusion occurred here!
1079
[[1, 2, -1, 1, 6], [2, 1, 1, 1, 0], [3, 1, 1, 1, 0]]
1080
1081
sage: A = Matrix(ZZ, 2, 2, [2,1,1,2])
1082
sage: p = 2
1083
sage: G2 = Genus_Symbol_p_adic_ring(p, p_adic_symbol(A, p, 2)); G2
1084
Genus symbol at 2 : [[0, 2, 3, 0, 0]]
1085
sage: G2.canonical_symbol()
1086
[[0, 2, -1, 0, 0]]
1087
1088
1089
sage: A = DiagonalQuadraticForm(ZZ, [1,2,3,4]).Hessian_matrix()
1090
sage: p = 3
1091
sage: G3 = Genus_Symbol_p_adic_ring(p, p_adic_symbol(A, p, 2)); G3
1092
Genus symbol at 3 : [[0, 3, 1], [1, 1, -1]]
1093
sage: G3.canonical_symbol()
1094
[[0, 3, 1], [1, 1, -1]]
1095
1096
1097
NOTES:
1098
See Conway-Sloane 3rd edition, pp. 381-382 for definitions and examples.
1099
1100
TO DO:
1101
- Add an example where sign walking occurs!
1102
"""
1103
symbol = self._symbol
1104
if self._prime == 2:
1105
if self._canonical_symbol is None:
1106
self._canonical_symbol = canonical_2_adic_reduction(symbol)
1107
return self._canonical_symbol
1108
else:
1109
return self._symbol
1110
1111
1112
1113
def symbol_tuple_list(self):
1114
"""
1115
Returns the underlying list of lists of integers defining the genus symbol.
1116
1117
INPUT:
1118
None
1119
1120
OUTPUT:
1121
list of lists of integers
1122
1123
EXAMPLES:
1124
sage: from sage.quadratic_forms.genera.genus import p_adic_symbol
1125
sage: from sage.quadratic_forms.genera.genus import Genus_Symbol_p_adic_ring
1126
1127
sage: A = DiagonalQuadraticForm(ZZ, [1,2,3,4]).Hessian_matrix()
1128
sage: p = 3
1129
sage: G3 = Genus_Symbol_p_adic_ring(p, p_adic_symbol(A, p, 2)); G3
1130
Genus symbol at 3 : [[0, 3, 1], [1, 1, -1]]
1131
sage: G3.symbol_tuple_list()
1132
[[0, 3, 1], [1, 1, -1]]
1133
sage: type(G3.symbol_tuple_list())
1134
<type 'list'>
1135
1136
sage: A = DiagonalQuadraticForm(ZZ, [1,2,3,4]).Hessian_matrix()
1137
sage: p = 2
1138
sage: G2 = Genus_Symbol_p_adic_ring(p, p_adic_symbol(A, p, 2)); G2
1139
Genus symbol at 2 : [[1, 2, 3, 1, 4], [2, 1, 1, 1, 1], [3, 1, 1, 1, 1]]
1140
sage: G2.symbol_tuple_list()
1141
[[1, 2, 3, 1, 4], [2, 1, 1, 1, 1], [3, 1, 1, 1, 1]]
1142
sage: type(G2.symbol_tuple_list())
1143
<type 'list'>
1144
1145
"""
1146
return self._symbol
1147
1148
1149
1150
def number_of_blocks(self):
1151
"""
1152
Returns the number of positive dimensional symbols/Jordan blocks
1153
1154
INPUT:
1155
None
1156
1157
OUTPUT:
1158
integer >= 0
1159
1160
EXAMPLES:
1161
sage: from sage.quadratic_forms.genera.genus import p_adic_symbol
1162
sage: from sage.quadratic_forms.genera.genus import Genus_Symbol_p_adic_ring
1163
1164
sage: A = DiagonalQuadraticForm(ZZ, [1,2,3,4]).Hessian_matrix()
1165
sage: p = 2
1166
sage: G2 = Genus_Symbol_p_adic_ring(p, p_adic_symbol(A, p, 2)); G2
1167
Genus symbol at 2 : [[1, 2, 3, 1, 4], [2, 1, 1, 1, 1], [3, 1, 1, 1, 1]]
1168
sage: G2.number_of_blocks()
1169
3
1170
1171
sage: A = DiagonalQuadraticForm(ZZ, [1,2,3,4]).Hessian_matrix()
1172
sage: p = 3
1173
sage: G3 = Genus_Symbol_p_adic_ring(p, p_adic_symbol(A, p, 2)); G3
1174
Genus symbol at 3 : [[0, 3, 1], [1, 1, -1]]
1175
sage: G3.number_of_blocks()
1176
2
1177
1178
"""
1179
return len(self._symbol)
1180
1181
1182
def determinant(self):
1183
"""
1184
Returns the (p-part of the) determinant (square-class) of the
1185
Hessian matrix of the quadratic form (given by regarding the
1186
integral symmetric matrix which generated this genus symbol as
1187
the Gram matrix of Q) associated to this local genus symbol.
1188
1189
INPUT:
1190
None
1191
1192
OUTPUT:
1193
an integer
1194
1195
EXAMPLES:
1196
sage: from sage.quadratic_forms.genera.genus import p_adic_symbol
1197
sage: from sage.quadratic_forms.genera.genus import Genus_Symbol_p_adic_ring
1198
1199
sage: A = DiagonalQuadraticForm(ZZ, [1,2,3,4]).Hessian_matrix()
1200
sage: p = 2
1201
sage: G2 = Genus_Symbol_p_adic_ring(p, p_adic_symbol(A, p, 2)); G2
1202
Genus symbol at 2 : [[1, 2, 3, 1, 4], [2, 1, 1, 1, 1], [3, 1, 1, 1, 1]]
1203
sage: G2.determinant()
1204
128
1205
1206
sage: A = DiagonalQuadraticForm(ZZ, [1,2,3,4]).Hessian_matrix()
1207
sage: p = 3
1208
sage: G3 = Genus_Symbol_p_adic_ring(p, p_adic_symbol(A, p, 2)); G3
1209
Genus symbol at 3 : [[0, 3, 1], [1, 1, -1]]
1210
sage: G3.determinant()
1211
3
1212
"""
1213
p = self._prime
1214
return misc.prod([ p**(s[0]*s[1]) for s in self._symbol ])
1215
1216
1217
def rank(self):
1218
"""
1219
Returns the dimension of a quadratic form associated to this genus symbol.
1220
1221
TO DO: DELETE THIS METHOD IN FAVOR OF THE dimension() METHOD BELOW!
1222
1223
1224
INPUT:
1225
None
1226
1227
OUTPUT:
1228
an integer >= 0
1229
1230
EXAMPLES:
1231
sage: from sage.quadratic_forms.genera.genus import p_adic_symbol
1232
sage: from sage.quadratic_forms.genera.genus import Genus_Symbol_p_adic_ring
1233
1234
sage: A = DiagonalQuadraticForm(ZZ, [1,2,3,4]).Hessian_matrix()
1235
sage: p = 2
1236
sage: G2 = Genus_Symbol_p_adic_ring(p, p_adic_symbol(A, p, 2)); G2
1237
Genus symbol at 2 : [[1, 2, 3, 1, 4], [2, 1, 1, 1, 1], [3, 1, 1, 1, 1]]
1238
sage: G2.rank()
1239
4
1240
1241
sage: A = DiagonalQuadraticForm(ZZ, [1,2,3,4]).Hessian_matrix()
1242
sage: p = 3
1243
sage: G3 = Genus_Symbol_p_adic_ring(p, p_adic_symbol(A, p, 2)); G3
1244
Genus symbol at 3 : [[0, 3, 1], [1, 1, -1]]
1245
sage: G3.rank()
1246
4
1247
1248
"""
1249
return sum([ s[1] for s in self._symbol ])
1250
1251
1252
def dimension(self):
1253
"""
1254
Returns the dimension of a quadratic form associated to this genus symbol.
1255
1256
INPUT:
1257
None
1258
1259
OUTPUT:
1260
an integer >= 0
1261
1262
EXAMPLES:
1263
sage: from sage.quadratic_forms.genera.genus import p_adic_symbol
1264
sage: from sage.quadratic_forms.genera.genus import Genus_Symbol_p_adic_ring
1265
1266
sage: A = DiagonalQuadraticForm(ZZ, [1,2,3,4]).Hessian_matrix()
1267
sage: p = 2
1268
sage: G2 = Genus_Symbol_p_adic_ring(p, p_adic_symbol(A, p, 2)); G2
1269
Genus symbol at 2 : [[1, 2, 3, 1, 4], [2, 1, 1, 1, 1], [3, 1, 1, 1, 1]]
1270
sage: G2.dimension()
1271
4
1272
1273
sage: A = DiagonalQuadraticForm(ZZ, [1,2,3,4]).Hessian_matrix()
1274
sage: p = 3
1275
sage: G3 = Genus_Symbol_p_adic_ring(p, p_adic_symbol(A, p, 2)); G3
1276
Genus symbol at 3 : [[0, 3, 1], [1, 1, -1]]
1277
sage: G3.dimension()
1278
4
1279
1280
"""
1281
return self.rank()
1282
1283
1284
def excess(self):
1285
"""
1286
Returns the p-excess of the quadratic form whose Hessian
1287
matrix is the symmetric matrix A. When p = 2 the p-excess is
1288
called the oddity.
1289
1290
WARNING/NOTE: This normalization seems non-standard, and we
1291
should review this entire class to make sure that we have our
1292
doubling conventions straight throughout!
1293
1294
REFERENCE:
1295
Conway and Sloane Book, 3rd edition, pp 370-371.
1296
1297
INPUT:
1298
None
1299
1300
OUTPUT:
1301
an integer
1302
1303
EXAMPLES:
1304
sage: from sage.quadratic_forms.genera.genus import p_adic_symbol
1305
sage: from sage.quadratic_forms.genera.genus import Genus_Symbol_p_adic_ring
1306
1307
sage: AC = diagonal_matrix(ZZ, [1,3,-3])
1308
sage: p=2; Genus_Symbol_p_adic_ring(p, p_adic_symbol(AC, p, 2)).excess()
1309
1
1310
sage: p=3; Genus_Symbol_p_adic_ring(p, p_adic_symbol(AC, p, 2)).excess()
1311
0
1312
sage: p=5; Genus_Symbol_p_adic_ring(p, p_adic_symbol(AC, p, 2)).excess()
1313
0
1314
sage: p=7; Genus_Symbol_p_adic_ring(p, p_adic_symbol(AC, p, 2)).excess()
1315
0
1316
sage: p=11; Genus_Symbol_p_adic_ring(p, p_adic_symbol(AC, p, 2)).excess()
1317
0
1318
1319
sage: AC = 2 * diagonal_matrix(ZZ, [1,3,-3])
1320
sage: p=2; Genus_Symbol_p_adic_ring(p, p_adic_symbol(AC, p, 2)).excess()
1321
1
1322
sage: p=3; Genus_Symbol_p_adic_ring(p, p_adic_symbol(AC, p, 2)).excess()
1323
0
1324
sage: p=5; Genus_Symbol_p_adic_ring(p, p_adic_symbol(AC, p, 2)).excess()
1325
0
1326
sage: p=7; Genus_Symbol_p_adic_ring(p, p_adic_symbol(AC, p, 2)).excess()
1327
0
1328
sage: p=11; Genus_Symbol_p_adic_ring(p, p_adic_symbol(AC, p, 2)).excess()
1329
0
1330
1331
sage: A = 2*diagonal_matrix(ZZ, [1,2,3,4])
1332
sage: p=2; Genus_Symbol_p_adic_ring(p, p_adic_symbol(A, p, 2)).excess()
1333
2
1334
sage: p=3; Genus_Symbol_p_adic_ring(p, p_adic_symbol(A, p, 2)).excess()
1335
6
1336
sage: p=5; Genus_Symbol_p_adic_ring(p, p_adic_symbol(A, p, 2)).excess()
1337
0
1338
sage: p=7; Genus_Symbol_p_adic_ring(p, p_adic_symbol(A, p, 2)).excess()
1339
0
1340
sage: p=11; Genus_Symbol_p_adic_ring(p, p_adic_symbol(A, p, 2)).excess()
1341
0
1342
1343
"""
1344
p = self._prime
1345
if self._prime == 2:
1346
k = 0
1347
for s in self._symbol:
1348
if s[0]%2 == 1 and s[2] in (3,5):
1349
k += 1
1350
return Integer(sum([ s[4] for s in self._symbol ]) + 4*k).mod(8)
1351
else:
1352
k = 0
1353
for s in self._symbol:
1354
if s[0]%2 == 1 and s[2] == -1:
1355
k += 1
1356
return Integer(sum([ s[1]*(p**s[0]-1) for s in self._symbol ]) + 4*k).mod(8)
1357
1358
1359
1360
def trains(self):
1361
"""
1362
Compute the indices for each of the trains in this local genus
1363
symbol if it is associated to the prime p=2 (and raise an
1364
error for all other primes).
1365
1366
INPUT:
1367
None
1368
1369
OUTPUT:
1370
a list of integers >= 0
1371
1372
EXAMPLES:
1373
sage: from sage.quadratic_forms.genera.genus import p_adic_symbol
1374
sage: from sage.quadratic_forms.genera.genus import Genus_Symbol_p_adic_ring
1375
1376
sage: A = DiagonalQuadraticForm(ZZ, [1,2,3,4]).Hessian_matrix()
1377
sage: p = 2
1378
sage: G2 = Genus_Symbol_p_adic_ring(p, p_adic_symbol(A, p, 2)); G2
1379
Genus symbol at 2 : [[1, 2, 3, 1, 4], [2, 1, 1, 1, 1], [3, 1, 1, 1, 1]]
1380
sage: G2.trains()
1381
[[0, 1, 2]]
1382
1383
"""
1384
## Check that p = 2
1385
if self._prime != 2:
1386
raise TypeError, "trains() only makes sense when the prime of the p_adic_Genus_Symbol is p=2"
1387
symbol = self._symbol
1388
compartments = canonical_2_adic_compartments(symbol)
1389
return canonical_2_adic_trains(symbol, compartments)
1390
1391
1392
def compartments(self):
1393
"""
1394
Compute the indices for each of the compartments in this local genus
1395
symbol if it is associated to the prime p=2 (and raise an
1396
error for all other primes).
1397
1398
INPUT:
1399
None
1400
1401
OUTPUT:
1402
a list of integers >= 0
1403
1404
EXAMPLES:
1405
sage: from sage.quadratic_forms.genera.genus import p_adic_symbol
1406
sage: from sage.quadratic_forms.genera.genus import Genus_Symbol_p_adic_ring
1407
1408
sage: A = DiagonalQuadraticForm(ZZ, [1,2,3,4]).Hessian_matrix()
1409
sage: p = 2
1410
sage: G2 = Genus_Symbol_p_adic_ring(p, p_adic_symbol(A, p, 2)); G2
1411
Genus symbol at 2 : [[1, 2, 3, 1, 4], [2, 1, 1, 1, 1], [3, 1, 1, 1, 1]]
1412
sage: G2.compartments()
1413
[[0, 1, 2]]
1414
1415
"""
1416
## Check that p = 2
1417
if self._prime != 2:
1418
raise TypeError, "compartments() only makes sense when the prime of the p_adic_Genus_Symbol is p=2"
1419
symbol = self._symbol
1420
return canonical_2_adic_compartments(symbol)
1421
1422
1423
1424
1425
1426
1427
class GenusSymbol_global_ring(object):
1428
"""
1429
This represents a collection of local genus symbols (at primes)
1430
and signature information which represent the genus of a
1431
non-degenerate integral lattice.
1432
"""
1433
1434
def __init__(self, A, max_elem_divisors=None):
1435
"""
1436
Initialize a global genus symbol from a non-degenerate
1437
integral gram matrix (and possibly information about its
1438
largest elementary divisors).
1439
1440
INPUT:
1441
A -- a symmetric matrix with integer coefficients
1442
max_elem_divisors -- the input precision for valuation of
1443
maximal p-elementary divisor. (OPTIONAL)
1444
1445
OUTPUT:
1446
None
1447
1448
EXAMPLES:
1449
sage: from sage.quadratic_forms.genera.genus import GenusSymbol_global_ring
1450
1451
sage: A = DiagonalQuadraticForm(ZZ, [1,2,3,4]).Hessian_matrix()
1452
sage: G = GenusSymbol_global_ring(A);G
1453
Genus of [2 0 0 0]
1454
[0 4 0 0]
1455
[0 0 6 0]
1456
[0 0 0 8]
1457
sage: G == loads(dumps(G))
1458
True
1459
1460
"""
1461
D = A.determinant()
1462
D = 2*D
1463
prms = [ p[0] for p in D.factor() ]
1464
self._representative = A
1465
self._signature = signature_pair_of_matrix(A)
1466
self._local_symbols = []
1467
for p in prms:
1468
if max_elem_divisors is None:
1469
val = D.valuation(p)
1470
symbol = p_adic_symbol(A, p, val = val)
1471
G = Genus_Symbol_p_adic_ring(p, symbol)
1472
self._local_symbols.append(G)
1473
1474
1475
def __repr__(self):
1476
"""
1477
Returns a string representing the global genus symbol.
1478
1479
INPUT:
1480
None
1481
1482
OUTPUT:
1483
a string
1484
1485
EXAMPLES:
1486
sage: from sage.quadratic_forms.genera.genus import GenusSymbol_global_ring
1487
1488
sage: A = DiagonalQuadraticForm(ZZ, [1,2,3,4]).Hessian_matrix()
1489
sage: GS = GenusSymbol_global_ring(A)
1490
sage: GS.__repr__()
1491
'Genus of [2 0 0 0]\n[0 4 0 0]\n[0 0 6 0]\n[0 0 0 8]'
1492
1493
"""
1494
return "Genus of %s"%self._representative
1495
1496
1497
1498
def __eq__(self, other):
1499
"""
1500
Determines if two global genus symbols are equal (not just equivalent!).
1501
1502
INPUT:
1503
a GenusSymbol_global_ring object
1504
1505
OUTPUT:
1506
boolean
1507
1508
EXAMPLES:
1509
sage: from sage.quadratic_forms.genera.genus import GenusSymbol_global_ring
1510
1511
sage: A1 = DiagonalQuadraticForm(ZZ, [1,2,3,4]).Hessian_matrix()
1512
sage: GS1 = GenusSymbol_global_ring(A1)
1513
sage: A2 = DiagonalQuadraticForm(ZZ, [1,2,3,5]).Hessian_matrix()
1514
sage: GS2 = GenusSymbol_global_ring(A2)
1515
1516
sage: GS1 == GS2
1517
False
1518
1519
sage: GS2 == GS1
1520
False
1521
1522
sage: GS1 == GS1
1523
True
1524
1525
sage: GS2 == GS2
1526
True
1527
1528
"""
1529
if self is other:
1530
return True
1531
t = len(self._local_symbols)
1532
if t != len(other._local_symbols):
1533
return False
1534
for i in range(t):
1535
if self._local_symbols[i] != other._local_symbols[i]:
1536
return False
1537
return True
1538
1539
1540
1541
def __ne__(self, other):
1542
"""
1543
Determines if two global genus symbols are unequal (not just inequivalent!).
1544
1545
INPUT:
1546
a GenusSymbol_global_ring object
1547
1548
OUTPUT:
1549
boolean
1550
1551
EXAMPLES:
1552
sage: from sage.quadratic_forms.genera.genus import GenusSymbol_global_ring
1553
1554
sage: A1 = DiagonalQuadraticForm(ZZ, [1,2,3,4]).Hessian_matrix()
1555
sage: GS1 = GenusSymbol_global_ring(A1)
1556
sage: A2 = DiagonalQuadraticForm(ZZ, [1,2,3,5]).Hessian_matrix()
1557
sage: GS2 = GenusSymbol_global_ring(A2)
1558
1559
sage: GS1 != GS2
1560
True
1561
1562
sage: GS2 != GS1
1563
True
1564
1565
sage: GS1 != GS1
1566
False
1567
1568
sage: GS2 != GS2
1569
False
1570
1571
"""
1572
return not self.__eq__(other)
1573
1574
1575
def signature_pair_of_matrix(self):
1576
"""
1577
Returns the signature pair (p, n) of the (non-degenerate)
1578
global genus symbol, where p is the number of positive
1579
eigenvalues and n is the number of negative eigenvalues.
1580
1581
INPUT:
1582
None
1583
1584
OUTPUT:
1585
a pair of integers (p, n) each >= 0
1586
1587
EXAMPLES:
1588
sage: from sage.quadratic_forms.genera.genus import GenusSymbol_global_ring
1589
1590
sage: A = DiagonalQuadraticForm(ZZ, [1,-2,3,4,8,-11]).Hessian_matrix()
1591
sage: GS = GenusSymbol_global_ring(A)
1592
sage: GS.signature_pair_of_matrix()
1593
(4, 2)
1594
1595
"""
1596
return self._signature
1597
1598
1599
def determinant(self):
1600
"""
1601
Returns the determinant of this genus, where the determinant
1602
is the Hessian determinant of the quadratic form whose Gram
1603
matrix is the Gram matrix giving rise to this global genus
1604
symbol.
1605
1606
INPUT:
1607
None
1608
1609
OUTPUT:
1610
an integer
1611
1612
EXAMPLES:
1613
sage: from sage.quadratic_forms.genera.genus import GenusSymbol_global_ring
1614
1615
sage: A = DiagonalQuadraticForm(ZZ, [1,-2,3,4]).Hessian_matrix()
1616
sage: GS = GenusSymbol_global_ring(A)
1617
sage: GS.determinant()
1618
-384
1619
1620
"""
1621
r, s = self.signature_pair_of_matrix()
1622
return (-1)**s*misc.prod([ G.determinant() for G in self._local_symbols ])
1623
1624