-
Notifications
You must be signed in to change notification settings - Fork 0
/
sqrt5.py
1490 lines (1264 loc) · 47.1 KB
/
sqrt5.py
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
296
297
298
299
300
301
302
303
304
305
306
307
308
309
310
311
312
313
314
315
316
317
318
319
320
321
322
323
324
325
326
327
328
329
330
331
332
333
334
335
336
337
338
339
340
341
342
343
344
345
346
347
348
349
350
351
352
353
354
355
356
357
358
359
360
361
362
363
364
365
366
367
368
369
370
371
372
373
374
375
376
377
378
379
380
381
382
383
384
385
386
387
388
389
390
391
392
393
394
395
396
397
398
399
400
401
402
403
404
405
406
407
408
409
410
411
412
413
414
415
416
417
418
419
420
421
422
423
424
425
426
427
428
429
430
431
432
433
434
435
436
437
438
439
440
441
442
443
444
445
446
447
448
449
450
451
452
453
454
455
456
457
458
459
460
461
462
463
464
465
466
467
468
469
470
471
472
473
474
475
476
477
478
479
480
481
482
483
484
485
486
487
488
489
490
491
492
493
494
495
496
497
498
499
500
501
502
503
504
505
506
507
508
509
510
511
512
513
514
515
516
517
518
519
520
521
522
523
524
525
526
527
528
529
530
531
532
533
534
535
536
537
538
539
540
541
542
543
544
545
546
547
548
549
550
551
552
553
554
555
556
557
558
559
560
561
562
563
564
565
566
567
568
569
570
571
572
573
574
575
576
577
578
579
580
581
582
583
584
585
586
587
588
589
590
591
592
593
594
595
596
597
598
599
600
601
602
603
604
605
606
607
608
609
610
611
612
613
614
615
616
617
618
619
620
621
622
623
624
625
626
627
628
629
630
631
632
633
634
635
636
637
638
639
640
641
642
643
644
645
646
647
648
649
650
651
652
653
654
655
656
657
658
659
660
661
662
663
664
665
666
667
668
669
670
671
672
673
674
675
676
677
678
679
680
681
682
683
684
685
686
687
688
689
690
691
692
693
694
695
696
697
698
699
700
701
702
703
704
705
706
707
708
709
710
711
712
713
714
715
716
717
718
719
720
721
722
723
724
725
726
727
728
729
730
731
732
733
734
735
736
737
738
739
740
741
742
743
744
745
746
747
748
749
750
751
752
753
754
755
756
757
758
759
760
761
762
763
764
765
766
767
768
769
770
771
772
773
774
775
776
777
778
779
780
781
782
783
784
785
786
787
788
789
790
791
792
793
794
795
796
797
798
799
800
801
802
803
804
805
806
807
808
809
810
811
812
813
814
815
816
817
818
819
820
821
822
823
824
825
826
827
828
829
830
831
832
833
834
835
836
837
838
839
840
841
842
843
844
845
846
847
848
849
850
851
852
853
854
855
856
857
858
859
860
861
862
863
864
865
866
867
868
869
870
871
872
873
874
875
876
877
878
879
880
881
882
883
884
885
886
887
888
889
890
891
892
893
894
895
896
897
898
899
900
901
902
903
904
905
906
907
908
909
910
911
912
913
914
915
916
917
918
919
920
921
922
923
924
925
926
927
928
929
930
931
932
933
934
935
936
937
938
939
940
941
942
943
944
945
946
947
948
949
950
951
952
953
954
955
956
957
958
959
960
961
962
963
964
965
966
967
968
969
970
971
972
973
974
975
976
977
978
979
980
981
982
983
984
985
986
987
988
989
990
991
992
993
994
995
996
997
998
999
1000
#########################################################################
# Copyright (C) 2010-2012 William Stein <wstein@gmail.com>
#
# Distributed under the terms of the GNU General Public License (GPL)
#
# http://www.gnu.org/licenses/
#########################################################################
"""
Hilbert Modular Forms over Q(sqrt(5)) of Weight (2,2)
AUTHORS:
- William Stein (2010, 2011, 2012)
This file contains an implementation of Dembele's quaternion algebra
algorithm (from "Explicit Computations of Hilbert Modular Forms on
Q(sqrt(5))") for computing Hilbert modular forms of weight (2,2) and
ramified or split prime level. The modules sqrt5_fast contains a much
faster version.
EXAMPLES:
LEVEL 31::
sage: from sage.modular.hilbert.sqrt5 import THETA, hecke_ops, F
sage: B.<i,j,k> = QuaternionAlgebra(F,-1,-1)
sage: c = F.factor(31)[1][0]
sage: P = F.primes_above(5)[0]
sage: TH = THETA(20) # about 1 minute
pi = [...]
Sorting through 22440 elements
sage: T = hecke_ops(c, TH); T # random output do to choice of basis
[(5, a + 2, [1 5]
[3 3]), (9, 3, [5 5]
[3 7]), (11, a + 3, [ 2 10]
[ 6 6]), (11, 2*a + 3, [7 5]
[3 9]), (19, a + 4, [10 10]
[ 6 14]), (19, 3*a + 4, [ 5 15]
[ 9 11])]
sage: for nm,p,t in T:
... print nm, p, t.charpoly().factor()
5 a + 2 (x - 6) * (x + 2)
9 3 (x - 10) * (x - 2)
11 a + 3 (x - 12) * (x + 4)
11 2*a + 3 (x - 12) * (x - 4)
19 a + 4 (x - 20) * (x - 4)
19 3*a + 4 (x - 20) * (x + 4)
LEVEL 41::
sage: from sage.modular.hilbert.sqrt5 import THETA, hecke_ops, F
sage: B.<i,j,k> = QuaternionAlgebra(F,-1,-1)
sage: F.primes_above(41)
[Fractional ideal (a - 7), Fractional ideal (a + 6)]
sage: c = F.primes_above(41)[0]
sage: TH = THETA(11) # about 30 seconds
pi = [...]
Sorting through 6660 elements
sage: T = hecke_ops(c, TH); T # random output do to choice of basis
[(5, a + 2, [4 2]
[5 1]), (9, 3, [ 6 4]
[10 0]), (11, a + 3, [10 2]
[ 5 7]), (11, 2*a + 3, [ 8 4]
[10 2])]
sage: for nm,p,t in T:
... print nm, p, t.charpoly().factor()
5 a + 2 (x - 6) * (x + 1)
9 3 (x - 10) * (x + 4)
11 a + 3 (x - 12) * (x - 5)
11 2*a + 3 (x - 12) * (x + 2)
LEVEL 389!:
This relies on having TH from above (say from the level 31 block above)::
sage: F.primes_above(389)
[Fractional ideal (18*a - 5), Fractional ideal (-18*a + 13)]
sage: c = F.primes_above(389)[0]
sage: T = hecke_ops(c, TH)
sage: for nm,p,t in T:
... print nm, p, t.charpoly().factor()
5 a + 2 (x - 6) * (x^2 + 4*x - 1) * (x^2 - x - 4)^2
9 3 (x - 10) * (x^2 + 3*x - 9) * (x^4 - 5*x^3 + 3*x^2 + 6*x - 4)
11 a + 3 (x - 12) * (x + 3)^2 * (x^4 - 17*x^2 + 68)
11 2*a + 3 (x - 12) * (x^2 + 5*x + 5) * (x^4 - x^3 - 23*x^2 + 18*x + 52)
"""
from sage.all import (NumberField, polygen, QQ, ZZ, QuaternionAlgebra,
cached_function, disk_cached_function)
x = polygen(QQ,'x')
F = NumberField(x**2 - x -1, 'a')
O_F = F.ring_of_integers()
B = QuaternionAlgebra(F, -1, -1, 'i,j,k')
def modp_splitting(p):
"""
INPUT:
- p -- ideal of the number field K = B.base() with ring O of integers.
OUTPUT:
- matrices I, J in M_2(O/p) such that i |--> I and j |--> J defines
an algebra morphism, i.e., I^2=a, J^2=b, I*J=-J*I.
EXAMPLES::
sage: from sage.modular.hilbert.sqrt5 import F, B, modp_splitting
sage: c = F.factor(31)[0][0]
sage: modp_splitting(c)
(
[ 0 30] [18 4]
[ 1 0], [ 4 13]
)
sage: c = F.factor(37)[0][0]; c
Fractional ideal (37)
sage: I, J = modp_splitting(c); I, J
(
[ 0 36] [23*abar + 21 36*abar + 8]
[ 1 0], [ 36*abar + 8 14*abar + 16]
)
sage: I^2
[36 0]
[ 0 36]
sage: J^2
[36 0]
[ 0 36]
sage: I*J == -J*I
True
AUTHOR: William Stein
"""
global B, F
# Inspired by the code in the function
# modp_splitting_data in algebras/quatalg/quaternion_algebra.py
if p.number_field() != B.base():
raise ValueError, "p must be a prime ideal in the base field of the quaternion algebra"
if not p.is_prime():
raise ValueError, "p must be prime"
if F is not p.number_field():
raise ValueError, "p must be a prime of %s"%F
k = p.residue_field()
from sage.all import MatrixSpace
M = MatrixSpace(k, 2)
i2, j2 = B.invariants()
i2 = k(i2); j2 = k(j2)
if k.characteristic() == 2:
if i2 == 0 or j2 == 0:
raise NotImplementedError
return M([0,1,1,0]), M([1,1,0,1])
# Find I -- just write it down
I = M([0,i2,1,0])
# Find J -- I figured this out by just writing out the general case
# and seeing what the parameters have to satisfy
i2inv = 1/i2
a = None
for b in list(k):
if not b: continue
c = j2 + i2inv * b*b
if c.is_square():
a = -c.sqrt()
break
if a is None:
# do a fallback search; needed in char 3 sometimes.
for J in M:
K = I*J
if J*J == j2 and K == -J*I:
return I, J, K
J = M([a,b,(j2-a*a)/b, -a])
K = I*J
assert K == -J*I, "bug in that I,J don't skew commute"
return I, J
def modp_splitting_map(p):
"""
Return a map from subset of B to 2x2 matrix space isomorphic
to R tensor OF/p.
INPUT:
- `B` -- quaternion algebra over F=Q(sqrt(5))
with invariants -1, -1.
- `p` -- prime ideal of F=Q(sqrt(5))
EXAMPLES::
sage: from sage.modular.hilbert.sqrt5 import F, B, modp_splitting_map
sage: i,j,k = B.gens()
sage: theta = modp_splitting_map(F.primes_above(5)[0])
sage: theta(i + j - k)
[2 1]
[3 3]
sage: s = 2 + 3*i - 2*j - 2*k
sage: theta(s)
[1 3]
[4 3]
sage: s.reduced_characteristic_polynomial()
x^2 - 4*x + 21
sage: theta(s).charpoly()
x^2 + x + 1
sage: s.reduced_characteristic_polynomial().change_ring(GF(5))
x^2 + x + 1
sage: theta = modp_splitting_map(F.primes_above(3)[0])
sage: smod = theta(s); smod
[2*abar + 1 abar + 1]
[ abar + 1 abar]
sage: smod^2 - 4*smod + 21
[0 0]
[0 0]
"""
I, J = modp_splitting(p)
F = p.residue_field()
def f(x):
return F(x[0]) + I*F(x[1]) + J*F(x[2]) + I*J*F(x[3])
return f
def icosian_gens():
"""
Return generators of the icosian group, as elements of the
Hamilton quaternion algebra B over Q(sqrt(5)).
AUTHOR: William Stein
EXAMPLES::
sage: from sage.modular.hilbert.sqrt5 import icosian_gens
sage: icosian_gens()
[i, j, k, -1/2 + 1/2*i + 1/2*j + 1/2*k, 1/2*i + 1/2*a*j + (-1/2*a + 1/2)*k]
sage: [a.reduced_norm() for a in icosian_gens()]
[1, 1, 1, 1, 1]
"""
global B, F
omega = F.gen() # (1+sqrt(5))/2
omega_bar = 1 - F.gen() # (1-sqrt(5))/2 = 1 - (1+sqrt(5))/2
return [B(v)/2 for v in [(0,2,0,0), (0,0,2,0),
(0,0,0,2), (-1,1,1,1), (0,1,omega,omega_bar)]]
def compute_all_icosians():
"""
Return a list of the elements of the Icosian group of order 120,
which we compute by generating enough products of icosian
generators.
EXAMPLES::
sage: from sage.modular.hilbert.sqrt5 import compute_all_icosians, all_icosians
sage: v = compute_all_icosians()
sage: len(v)
120
sage: v
[1/2 + 1/2*a*i + (-1/2*a + 1/2)*k, 1/2 + (-1/2*a + 1/2)*i + 1/2*a*j,..., -k, i, j, -i]
sage: assert set(v) == set(all_icosians()) # double check
"""
from sage.all import permutations, cartesian_product_iterator
Icos = []
ig = icosian_gens()
per = permutations(range(5))
exp = cartesian_product_iterator([range(1,i) for i in [5,5,5,4,5]])
for f in exp:
for p in per:
e0 = ig[p[0]]**f[0]
e1 = ig[p[1]]**f[1]
e2 = ig[p[2]]**f[2]
e3 = ig[p[3]]**f[3]
e4 = ig[p[4]]**f[4]
elt = e0*e1*e2*e3*e4
if elt not in Icos:
Icos.append(elt)
if len(Icos) == 120:
return Icos
@cached_function
def all_icosians():
"""
Return a list of all 120 icosians, from a precomputed table.
EXAMPLES::
sage: from sage.modular.hilbert.sqrt5 import all_icosians
sage: v = all_icosians()
sage: len(v)
120
"""
s = '[1+a*i+(-a+1)*k,1+(-a+1)*i+a*j,-a+i+(a-1)*j,1+(-a)*i+(-a+1)*k,-a-j+(-a+1)*k,1+(a-1)*i+(-a)*j,-1+(-a)*i+(a-1)*k,-1+(a-1)*i+(-a)*j,-a+1+(-a)*j-k,-1+a*i+(-a+1)*k,-a+1+i+(-a)*k,-1+(-a+1)*i+(-a)*j,(-a+1)*i+j+(-a)*k,a-1+(-a)*j+k,(a-1)*i-j+a*k,a+i+(-a+1)*j,1+a*i+(a-1)*k,-1+(-a)*i+(-a+1)*k,a*i+(-a+1)*j-k,a-1+i+a*k,(-a)*i+(a-1)*j+k,a+j+(-a+1)*k,1+(-a+1)*i+(-a)*j,-1+(a-1)*i+a*j,a-i+(-a+1)*j,-1+a*i+(a-1)*k,a+j+(a-1)*k,-1+(-a+1)*i+a*j,(-a+1)*i+j+a*k,a-1+a*j+k,-a+i+(-a+1)*j,1+(-a)*i+(a-1)*k,a*i+(-a+1)*j+k,a-1-i+a*k,-a+j+(a-1)*k,1+(a-1)*i+a*j,(a-1)*i+j+a*k,-a+1+a*j+k,(-a)*i+(-a+1)*j+k,-a+1+i+a*k,-a-i+(-a+1)*j,-a+1+(-a)*j+k,(-a+1)*i-j+a*k,(a-1)*i+j+(-a)*k,a+i+(a-1)*j,a-1+a*j-k,-a+j+(-a+1)*k,-a+1-i+a*k,a*i+(a-1)*j+k,(-a)*i+(-a+1)*j-k,a-j+(a-1)*k,a-1+i+(-a)*k,-1+i+j+k,-1-i-j+k,1-i-j-k,1+i-j+k,a-j+(-a+1)*k,-1+i-j-k,1-i+j+k,(a-1)*i-j+(-a)*k,-a-j+(a-1)*k,1+i+j-k,a*i+(a-1)*j-k,-1-i+j-k,a-1+(-a)*j-k,-a+1+a*j-k,(-a)*i+(a-1)*j-k,-a+1-i+(-a)*k,-a-i+(a-1)*j,a-1-i+(-a)*k,(-a+1)*i-j+(-a)*k,a-i+(a-1)*j,-i+(-a)*j+(a-1)*k,i+a*j+(a-1)*k,i+a*j+(-a+1)*k,-i+a*j+(a-1)*k,-i+a*j+(-a+1)*k,i+(-a)*j+(a-1)*k,-i+(-a)*j+(-a+1)*k,i+(-a)*j+(-a+1)*k,-1-i-j-k,1+i+j+k,2,-a+1+a*i-j,-2,-a+(-a+1)*i-k,a-1+a*i+j,a+(-a+1)*i+k,a-1+(-a)*i+j,1+(-a+1)*j+(-a)*k,-a+1+a*i+j,-1+(-a+1)*j+a*k,a+(a-1)*i+k,-1+(a-1)*j+a*k,-a+(-a+1)*i+k,1+(-a+1)*j+a*k,-a+1+(-a)*i+j,-a+(a-1)*i+k,a-1+a*i-j,1+(a-1)*j+a*k,a+(-a+1)*i-k,-1+(-a+1)*j+(-a)*k,-a+1+(-a)*i-j,-a+(a-1)*i-k,a-1+(-a)*i-j,1+(a-1)*j+(-a)*k,a+(a-1)*i-k,-1+(a-1)*j+(-a)*k,-1+i+j-k,1-i+j-k,1-i-j+k,-1-i+j+k,-1+i-j+k,1+i-j-k,2*k,(-2)*j,(-2)*k,2*i,2*j,(-2)*i]'
v = eval(s, {'a':F.gen(), 'i':B.gen(0), 'j':B.gen(1), 'k':B.gen(0)*B.gen(1)})
return [B(x)/B(2) for x in v]
def icosian_ring_gens():
"""
Return ring generators for the icosian ring (a maximal order) in the
quaternion algebra ramified only at infinity over F=Q(sqrt(5)).
These are generators over the ring of integers of F.
OUTPUT:
- list
EXAMPLES::
sage: from sage.modular.hilbert.sqrt5 import icosian_ring_gens
sage: icosian_ring_gens()
[1/2 + (1/2*a - 1/2)*i + 1/2*a*j, (1/2*a - 1/2)*i + 1/2*j + 1/2*a*k, 1/2*a*i + (1/2*a - 1/2)*j + 1/2*k, 1/2*i + 1/2*a*j + (1/2*a - 1/2)*k]
"""
global B, F
# See page 6 of Dembele.
# DO NOT CHANGE THESE! You'll break, e.g., the optimized code for
# writing elements in terms of these...
omega = F.gen()
omega_bar = 1 - F.gen()
return [B(v)/2 for v in [(1,-omega_bar,omega,0),
(0,-omega_bar,1,omega),
(0,omega,-omega_bar,1),
(0,1,omega,-omega_bar)]]
def icosian_ring_gens_over_ZZ():
"""
Return basis over ZZ for the icosian ring, which has ZZ-rank 8.
EXAMPLES::
sage: from sage.modular.hilbert.sqrt5 import icosian_ring_gens_over_ZZ
sage: v = icosian_ring_gens_over_ZZ(); v
[1/2 + (1/2*a - 1/2)*i + 1/2*a*j, (1/2*a - 1/2)*i + 1/2*j + 1/2*a*k, 1/2*a*i + (1/2*a - 1/2)*j + 1/2*k, 1/2*i + 1/2*a*j + (1/2*a - 1/2)*k, 1/2*a + 1/2*i + (1/2*a + 1/2)*j, 1/2*i + 1/2*a*j + (1/2*a + 1/2)*k, (1/2*a + 1/2)*i + 1/2*j + 1/2*a*k, 1/2*a*i + (1/2*a + 1/2)*j + 1/2*k]
sage: len(v)
8
"""
global B, F
I = icosian_ring_gens()
omega = F.gen()
return I + [omega*x for x in I]
def tensor_over_QQ_with_RR(prec=53):
"""
Return map from the quaternion algebra B to the tensor product of
B over QQ with RealField(prec), viewed as an 8-dimensional real
vector space.
INPUT:
- prec -- (integer: default 53); bits of real precision
OUTPUT:
- a Python function
EXAMPLES::
sage: from sage.modular.hilbert.sqrt5 import tensor_over_QQ_with_RR, B
sage: f = tensor_over_QQ_with_RR()
sage: B.gens()
[i, j, k]
sage: f(B.0)
(0.000000000000000, 1.00000000000000, 0.000000000000000, 0.000000000000000, 0.000000000000000, 1.00000000000000, 0.000000000000000, 0.000000000000000)
sage: f(B.1)
(0.000000000000000, 0.000000000000000, 1.00000000000000, 0.000000000000000, 0.000000000000000, 0.000000000000000, 1.00000000000000, 0.000000000000000)
sage: f = tensor_over_QQ_with_RR(20)
sage: f(B.0 - (1/9)*B.1)
(0.00000, 1.0000, -0.11111, 0.00000, 0.00000, 1.0000, -0.11111, 0.00000)
"""
global B, F
from sage.all import RealField
RR = RealField(prec=prec)
V = RR**8
S = F.embeddings(RR)
def f(x):
return V(sum([[sigma(a) for a in x] for sigma in S],[]))
return f
def modp_icosians(p):
"""
Return matrices of images of all 120 icosians modulo p.
INPUT:
- p -- *split* or ramified prime ideal of real quadratic field F=Q(sqrt(5))
OUTPUT:
- list of matrices modulo p.
EXAMPLES::
sage: from sage.modular.hilbert.sqrt5 import F, modp_icosians
sage: len(modp_icosians(F.primes_above(5)[0]))
120
sage: v = modp_icosians(F.primes_above(11)[0])
sage: len(v)
120
sage: v[0]
[10 8]
[ 8 1]
sage: v[-1]
[0 3]
[7 0]
"""
I, J = modp_splitting(p); K = I*J
k = p.residue_field()
G = [k(g[0]) + k(g[1])*I + k(g[2])*J + k(g[3])*K for g in icosian_gens()]
from sage.all import MatrixGroup
return [g.matrix() for g in MatrixGroup(G)]
class P1ModList(object):
"""
Object the represents the elements of the projective line modulo
a nonzero *prime* ideal of the ring of integers of Q(sqrt(5)).
Elements of the projective line are represented by elements of a 2-dimension
vector space over the residue field.
EXAMPLES::
We construct the projective line modulo the ideal (2), and illustrate
all the standard operations with it::
sage: from sage.modular.hilbert.sqrt5 import F, P1ModList
sage: P1 = P1ModList(F.primes_above(2)[0]); P1
Projective line over Residue field in abar of Fractional ideal (2)
sage: len(P1)
5
sage: list(P1)
[(0, 1), (1, 0), (1, abar), (1, abar + 1), (1, 1)]
sage: P1.random_element() # random output
(1, abar + 1)
sage: z = P1.random_element(); z # random output
(1, 0)
sage: z[0].parent()
Residue field in abar of Fractional ideal (2)
sage: g = z[0].parent().gen()
sage: P1.normalize((g,g))
(1, 1)
sage: P1((g,g))
(1, 1)
"""
def __init__(self, c):
"""
INPUT:
- c -- a nonzero prime of the ring of integers of Q(sqrt(5))
EXAMPLES::
sage: from sage.modular.hilbert.sqrt5 import F, P1ModList
sage: P1ModList(F.primes_above(3)[0])
Projective line over Residue field in abar of Fractional ideal (3)
sage: P1ModList(F.primes_above(11)[1])
Projective line over Residue field of Fractional ideal (-3*a + 1)
sage: list(P1ModList(F.primes_above(5)[0]))
[(0, 1), (1, 0), (1, 1), (1, 2), (1, 3), (1, 4)]
"""
self._c = c
F = c.residue_field()
V = F**2
self._V = V
self._F = F
self._list = [ V([0,1]) ] + [V([1,a]) for a in F]
for a in self._list:
a.set_immutable()
def random_element(self):
"""
Return a random element of this projective line.
sage: from sage.modular.hilbert.sqrt5 import F, P1ModList
sage: P1 = P1ModList(F.primes_above(13)[0]); P1
Projective line over Residue field in abar of Fractional ideal (13)
sage: P1.random_element() # random output
(1, 10*abar + 5)
"""
import random
return random.choice(self._list)
def normalize(self, uv):
"""
Normalize a representative element so it is either of the form
(1,*) if the first entry is nonzero, or of the form (0,1)
otherwise.
EXAMPLES::
sage: from sage.modular.hilbert.sqrt5 import F, P1ModList
sage: p = F.primes_above(13)[0]
sage: P1 = P1ModList(p)
sage: k = p.residue_field()
sage: g = k.gen()
sage: P1.normalize([3,4])
(1, 10)
sage: P1.normalize([g,g])
(1, 1)
sage: P1.normalize([0,g])
(0, 1)
"""
w = self._V(uv)
if w[0]:
w = (~w[0]) * w
w.set_immutable()
#assert w in self._list
return w
else:
return self._list[0]
def __len__(self):
"""
Return number of elements of this P1.
EXAMPLES::
sage: from sage.modular.hilbert.sqrt5 import F, P1ModList
sage: len(P1ModList(F.primes_above(3)[0]))
10
sage: len(P1ModList(F.primes_above(5)[0]))
6
sage: len(P1ModList(F.primes_above(19)[0]))
20
sage: len(P1ModList(F.primes_above(19)[1]))
20
"""
return len(self._list)
def __getitem__(self, i):
"""
Return i-th element.
EXAMPLES::
sage: from sage.modular.hilbert.sqrt5 import F, P1ModList
sage: P = P1ModList(F.primes_above(3)[0]); list(P)
[(0, 1), (1, 0), (1, 2*abar), (1, abar + 1), (1, abar + 2), (1, 2), (1, abar), (1, 2*abar + 2), (1, 2*abar + 1), (1, 1)]
sage: P[2]
(1, 2*abar)
"""
return self._list[i]
def __call__(self, x):
"""
Coerce x into this P1 list. Here x is anything that coerces to
the 2-dimensional vector space over the residue field. The
result is normalized (in fact this function is just an alias
for the normalize function).
EXAMPLES::
sage: from sage.modular.hilbert.sqrt5 import F, P1ModList
sage: p = F.primes_above(13)[0]
sage: k = p.residue_field(); g = k.gen()
sage: P1 = P1ModList(p)
sage: P1([3,4])
(1, 10)
sage: P1([g,g])
(1, 1)
sage: P1(P1([g,g]))
(1, 1)
"""
return self.normalize(x)
def __repr__(self):
"""
EXAMPLES::
sage: from sage.modular.hilbert.sqrt5 import F, P1ModList
sage: P1ModList(F.primes_above(19)[1]).__repr__()
'Projective line over Residue field of Fractional ideal (-4*a + 3)'
"""
return 'Projective line over %s'%self._F
def P1_orbits(p):
"""
INPUT:
- p -- a split or ramified prime of the integers of Q(sqrt(5)).
OUTPUT:
- ``orbits`` -- dictionary mapping elements of P1 to a choice of orbit rep
- ``reps`` -- list of representatives for the orbits
- ``P1`` -- the P1ModList object
AUTHOR: William Stein
EXAMPLES::
sage: from sage.modular.hilbert.sqrt5 import P1_orbits, F
sage: orbits, reps, P1 = P1_orbits(F.primes_above(5)[0])
sage: orbits # random output
{(1, 2): (1, 0), (0, 1): (1, 0), (1, 3): (1, 0), (1, 4): (1, 0), (1, 0): (1, 0), (1, 1): (1, 0)}
sage: reps # random output
[(1, 1)]
sage: len(reps)
1
sage: P1
Projective line over Residue field of Fractional ideal (2*a - 1)
sage: orbits, reps, P1 = P1_orbits(F.primes_above(41)[0])
sage: reps # random output
[(1, 40), (1, 5)]
sage: len(reps)
2
"""
global B
P1 = P1ModList(p)
ICO = modp_icosians(p)
def act(u, t):
return P1(u*t)
cur = P1.random_element()
reps = [cur]
orbits = {cur:cur}
while len(orbits) < len(P1):
for u in ICO:
s = act(u, cur)
if not orbits.has_key(s):
orbits[s] = cur
if len(orbits) < len(P1):
# choose element of P1 not a key
while True:
c = P1.random_element()
if c not in orbits:
cur = c
reps.append(cur)
orbits[cur] = cur
break
# done
return orbits, reps, P1
def P1_orbits2(p):
"""
INPUT:
- p -- a split or ramified prime of the integers of Q(sqrt(5)).
OUTPUT:
- list of disjoint sets of elements of P1 that are orbits
AUTHOR: William Stein
EXAMPLES::
sage: from sage.modular.hilbert.sqrt5 import P1_orbits2, F
sage: P1_orbits2(F.primes_above(5)[0]) # random output
[set([(1, 2), (0, 1), (1, 3), (1, 4), (1, 0), (1, 1)])]
sage: P1_orbits2(F.primes_above(11)[0]) # random output
[set([(0, 1), (1, 2), (1, 3), (1, 4), (1, 5), (1, 8), (1, 6), (1, 9), (1, 7), (1, 10), (1, 0), (1, 1)])]
sage: len(P1_orbits2(F.primes_above(41)[0]))
2
"""
global B
P1 = P1ModList(p)
ICO = modp_icosians(p)
orbits = []
while sum(len(x) for x in orbits) < len(P1):
v = P1.random_element()
skip = False
for O in orbits:
if v in O:
skip = True
break
if skip: continue
O = set([P1(g*v) for g in ICO])
orbits.append(O)
# Now make a dictionary
return orbits
def totally_pos_gen(p):
"""
Given a prime ideal p of a narrow class number 1 real quadratic
field, return a totally positive generator for p.
INPUT:
- p -- prime ideal of narrow class number 1 real
quadratic field
OUTPUT:
- generator of p that is totally positive
AUTHOR: William Stein
EXAMPLES::
sage: from sage.modular.hilbert.sqrt5 import F, totally_pos_gen
sage: g = totally_pos_gen(F.factor(19)[0][0]); g
3*a + 4
sage: g.norm()
19
sage: g.complex_embeddings()
[2.14589803375032, 8.85410196624968]
sage: for p in primes(14):
... for P, e in F.factor(p):
... g = totally_pos_gen(P)
... print P, g, g.complex_embeddings()
Fractional ideal (2) 2 [2.00000000000000, 2.00000000000000]
Fractional ideal (3) 3 [3.00000000000000, 3.00000000000000]
Fractional ideal (-2*a + 1) a + 2 [1.38196601125011, 3.61803398874989]
Fractional ideal (7) 7 [7.00000000000000, 7.00000000000000]
Fractional ideal (-3*a + 2) a + 3 [2.38196601125011, 4.61803398874989]
Fractional ideal (-3*a + 1) 2*a + 3 [1.76393202250021, 6.23606797749979]
Fractional ideal (13) 13 [13.0000000000000, 13.0000000000000]
"""
F = p.number_field()
assert F.degree() == 2 and F.discriminant() > 0
G = p.gens_reduced()
if len(G) != 1:
raise ValueError, "ideal not principal"
g = G[0]
from sage.all import RR
sigma = F.embeddings(RR)
e = [s(g) for s in sigma]
u = F.unit_group().gen(1)
ue = [s(u) for s in sigma]
if ue[0] > 0 and ue[1] < 0:
u *= -1
if e[0] < 0 and e[1] < 0:
return -g
elif e[0] < 0 and e[1] > 0:
if ue[0] < 0 and ue[1] > 0:
return u*g
else:
raise ValueError, "no totally positive generator"
elif e[0] > 0 and e[1] > 0:
return g
elif e[0] > 0 and e[1] < 0:
if ue[0] < 0 and ue[1] > 0:
return -u*g
else:
raise ValueError, "no totally positive generator"
assert False, "bug"
def gram_matrix_of_maximal_order(R):
"""
Return 8x8 Gram matrix of maximal order defined by R, which is assumed to be
a basis for maximal order over ZZ.
EXAMPLES::
sage: from sage.modular.hilbert.sqrt5 import icosian_ring_gens_over_ZZ, gram_matrix_of_maximal_order
sage: R = icosian_ring_gens_over_ZZ()
sage: gram_matrix_of_maximal_order(R)
[4 2 2 1 2 1 1 3]
[2 4 1 1 1 2 3 3]
[2 1 4 1 1 3 2 3]
[1 1 1 4 3 3 3 2]
[2 1 1 3 6 3 3 4]
[1 2 3 3 3 6 4 4]
[1 3 2 3 3 4 6 4]
[3 3 3 2 4 4 4 6]
"""
G = [[(R[i]*R[j].conjugate()).reduced_trace().trace()
for i in range(8)] for j in range(8)]
from sage.all import matrix
return matrix(ZZ, G)
def qfminim(qf, N):
"""
Call the PARI qfminim method on qf and 2*N, with smaller and
and smaller search range, until a MemoryError is *not* raised.
On a large-memory machine this will succeed the first time.
EXAMPLES::
sage: from sage.modular.hilbert.sqrt5 import icosian_ring_gens_over_ZZ, gram_matrix_of_maximal_order
sage: R = icosian_ring_gens_over_ZZ()
sage: G = gram_matrix_of_maximal_order(R)
sage: qf = pari(G)
sage: from sage.modular.hilbert.sqrt5 import icosian_ring_gens_over_ZZ, gram_matrix_of_maximal_order, qfminim
sage: n, m, v = qfminim(qf, 2)
sage: n
120
sage: m
4
sage: v[0]
[0, 0, 0, -1, 1, 0, 0, 0]~
"""
i = 32
while i>10:
try:
return qf.qfminim(2*N, 2**i)
except MemoryError:
i -= 1
def bounded_elements(N):
"""
Return elements in maximal order of B that have reduced norm
whose trace is at most N.
EXAMPLES::
sage: from sage.modular.hilbert.sqrt5 import bounded_elements
sage: X = bounded_elements(3)
sage: len(X)
180
sage: rnX = [a.reduced_norm() for a in X]
sage: set([a.trace() for a in rnX])
set([2, 3])
sage: set([a.norm() for a in rnX])
set([1])
sage: X = bounded_elements(5)
sage: len(X)
1200
sage: rnX = [a.reduced_norm() for a in X]
sage: set([a.trace() for a in rnX])
set([2, 3, 4, 5])
sage: set([a.norm() for a in rnX])
set([1, 4, 5])
"""
# Get our maximal order
R = icosian_ring_gens_over_ZZ()
# Compute Gram matrix of R
G = gram_matrix_of_maximal_order(R)
# Make PARI quadratic form
from sage.all import pari
qf = pari(G)
# Get the vectors of norm up to N.
# The 2 is because we had to scale by 2 to get
# rid of denominator in Gram matrix.
Z = qfminim(qf, N)
Z2 = Z[2].sage().transpose()
# For each vector, make the corresponding element of B.
# TODO: This step massively dominates the runtime, and can be
# easily made trivial with careful thought.
V = []
for i in range(Z2.nrows()):
w = Z2[i]
V.append(sum(w[j]*R[j] for j in range(8)))
return V
from sqrt5_tables import primes_of_bounded_norm
def THETA(N):
r"""
Return representative elements of the maximal order of `R` of
reduced norm `\pi_p` up to `N` modulo the left action of the units
of `R` (the icosians). Here `\pi_p` runs through totally positive
generators of the odd prime ideals with norm up to `N`.
INPUT:
- `N` -- a positive integer >= 4.
OUTPUT:
- dictionary with keys the totally positive generators of the
odd prime ideals with norm up to and including `N`, and values
a dictionary with values the actual elements of reduced norm
`\pi_p`.
AUTHOR: William Stein
EXAMPLES::
sage: from sage.modular.hilbert.sqrt5 import THETA
sage: d = THETA(9)
pi = [2, a + 2, 3]
Sorting through 2400 elements
sage: d.keys()
[a + 2, 3]
sage: d[3]
{(0, 1): a + 1/2*i + (-1/2*a + 1)*j + (-1/2*a + 1/2)*k, (1, 2): a - 1 + a*j, (1, abar): a - 1/2 + (1/2*a + 1/2)*i + (-1/2*a + 1)*j, (1, abar + 1): a - 1/2 + 1/2*i + 1/2*j + (-a + 1/2)*k, (1, abar + 2): a - 1 + (1/2*a + 1/2)*i + 1/2*j + (-1/2*a)*k, (1, 2*abar + 2): a - 1 + 1/2*a*i + (1/2*a + 1/2)*j + 1/2*k, (1, 2*abar): a - 1/2 + i + (-1/2*a + 1/2)*j + (-1/2*a)*k, (1, 2*abar + 1): a - 1/2 + (-1/2*a)*i + j + (-1/2*a + 1/2)*k, (1, 0): a + (-a + 1)*k, (1, 1): a - 1/2 + 1/2*a*i + j + (-1/2*a + 1/2)*k}
sage: k = d[3].keys(); k
[(0, 1), (1, 2), (1, abar), (1, abar + 1), (1, abar + 2), (1, 2*abar + 2), (1, 2*abar), (1, 2*abar + 1), (1, 0), (1, 1)]
sage: d[3][k[0]]
a + 1/2*i + (-1/2*a + 1)*j + (-1/2*a + 1/2)*k
sage: d[3][k[0]].reduced_norm()
3
"""
# ** NOTE: This algorithm will not scale up well, because there
# are so many vectors of bounded norm.
####################################################################
# Algorithm:
# * Enumerate set S of primes of norm <= N.
# * For each prime p in S:
# - Find a totally positive generator pi_p for p.
# - Compute mod-p local splitting theta_p.
# * Compute set X of elements in maximal order of B of norm <= N
# using the Gram matrix of the icosian ring (maximal order).
# * For each element z of X, compute the reduced norm of z.
# If it equals pi_p for one of the pi_p, compute the
# top row v of the reduced row echelon form of theta_p(z).
# Store v:z if there isn't already something for v.
# Also, of course, store this stuff separately for each p.
####################################################################
global B, F
# Get primes of norm up to N.
S = primes_of_bounded_norm(N)
# Find totally positive generators pi_p
pi = [totally_pos_gen(p) for p in S]
print "pi =",pi
# Compute traces of the generators, since that's what
# the bounded_elements command computes up to.
tr = [abs(x.trace()) for x in pi]
N = max(tr)
# A list that at least includes all elements (up to -1) whose
# reduced norm has trace at most N.
X = bounded_elements(N)
# Compute mod-p local splitting maps
theta_map = {}
for i, p in enumerate(S):
theta_map[pi[i]] = modp_splitting_map(p)
# Sort through the elements of X.
pi_set = set(pi)
# TODO: We skip the prime 2, since the mod-p splitting map is
# broken there.
pi_set.remove(2)
Theta = {}
for pi_p in pi_set:
Theta[pi_p] = {}
# The dictionary Theta will have keys the pi_p and
# the dictionary for pi_p has keys reduced vectors
# in (OF/p)^2. Here "reduced" just means "reduced
# row echelon form", so scaled so first entry is 1.
print "Sorting through %s elements"%len(X)
for a in X:
nrm = a.reduced_norm()
if nrm in pi_set:
# this is: mod right action of R^* acting on the right,
# so column echelon form
v = theta_map[nrm](a).transpose().echelon_form()[0]
## for reference, this is: mod left action of R^*,
## which is wrong, I guess:
# v = theta_map[nrm](a).echelon_form()[0]
z = Theta[nrm]
if z.has_key(v):
pass
else:
z[v] = a
return Theta
def hecke_ops(c, X):
orbits, reps, P1 = P1_orbits(c)
theta_c = modp_splitting_map(c)
def Tp(pi):
z = X[pi]
mat = []
for x in reps:
row = [0]*len(reps)
for _, w in z.iteritems():
w_c = theta_c(w)
y = w_c**(-1) * x
y_red = orbits[P1.normalize(y)]
row[reps.index(y_red)] += 1
mat.append(row)
from sage.all import matrix
return matrix(ZZ, mat)
ans = [(pi.norm(), pi, Tp(pi)) for pi in X.keys()]
ans.sort()
return ans
def hecke_ops2(c, X):
reduce, reps, P1 = P1_orbits(c)
theta_c = modp_splitting_map(c)
def Tp(pi):
z = X[pi]
mat = []
for x in reps:
print "x = %s, card = %s"%(x, len([M for M in reduce.keys() if reduce[M]==x]))
row = [0]*len(reps)
for _, w in z.iteritems():
w_c = theta_c(w)
y = w_c**(-1) * x
print "y =", y
y_red = reduce[P1(y)]
row[reps.index(y_red)] += 1
print "y_red =", y_red
mat.append(row)
from sage.all import matrix
return matrix(ZZ, mat)
ans = [(pi.norm(), pi, Tp(pi)) for pi in X.keys()]
ans.sort()
return ans
class AlphaZ:
def __init__(self, P):
"""
Computing elements with norm pi, where P=(pi).