Book a Demo!
CoCalc Logo Icon
StoreFeaturesDocsShareSupportNewsAboutPoliciesSign UpSign In
sagemath
GitHub Repository: sagemath/sagelib
Path: blob/master/sage/modular/modsym/p1list_nf.py
4069 views
1
r"""
2
Lists of Manin symbols (elements of `\mathbb{P}^1(R/N)`) over number fields
3
4
Lists of elements of `\mathbb{P}^1(R/N)` where `R` is the ring of integers of a number
5
field `K` and `N` is an integral ideal.
6
7
AUTHORS:
8
9
- Maite Aranes (2009): Initial version
10
11
EXAMPLES:
12
13
We define a P1NFList:
14
15
::
16
17
sage: k.<a> = NumberField(x^3 + 11)
18
sage: N = k.ideal(5, a^2 - a + 1)
19
sage: P = P1NFList(N); P
20
The projective line over the ring of integers modulo the Fractional ideal (5, a^2 - a + 1)
21
22
List operations with the P1NFList:
23
24
::
25
26
sage: len(P)
27
26
28
sage: [p for p in P]
29
[M-symbol (0: 1) of level Fractional ideal (5, a^2 - a + 1),
30
...
31
M-symbol (1: 2*a^2 + 2*a) of level Fractional ideal (5, a^2 - a + 1)]
32
33
The elements of the P1NFList are M-symbols:
34
35
::
36
37
sage: type(P[2])
38
<class 'sage.modular.modsym.p1list_nf.MSymbol'>
39
40
Definition of MSymbols:
41
42
::
43
44
sage: alpha = MSymbol(N, 3, a^2); alpha
45
M-symbol (3: a^2) of level Fractional ideal (5, a^2 - a + 1)
46
47
Find the index of the class of an M-Symbol `(c: d)` in the list:
48
49
::
50
51
sage: i = P.index(alpha)
52
sage: P[i].c*alpha.d - P[i].d*alpha.c in N
53
True
54
55
Lift an MSymbol to a matrix in `SL(2, R)`:
56
57
::
58
59
sage: alpha = MSymbol(N, a + 2, 3*a^2)
60
sage: alpha.lift_to_sl2_Ok()
61
[1, -4*a^2 + 9*a - 21, a + 2, a^2 - 3*a + 3]
62
sage: Ok = k.ring_of_integers()
63
sage: M = Matrix(Ok, 2, alpha.lift_to_sl2_Ok())
64
sage: det(M)
65
1
66
sage: M[1][1] - alpha.d in N
67
True
68
69
Lift an MSymbol from P1NFList to a matrix in `SL(2, R)`
70
71
::
72
73
sage: P[3]
74
M-symbol (1: -2*a) of level Fractional ideal (5, a^2 - a + 1)
75
sage: P.lift_to_sl2_Ok(3)
76
[0, -1, 1, -2*a]
77
"""
78
#*****************************************************************************
79
# Copyright (C) 2009, Maite Aranes <[email protected]>
80
#
81
# Distributed under the terms of the GNU General Public License (GPL)
82
# http://www.gnu.org/licenses/
83
#*****************************************************************************
84
85
86
from sage.structure.sage_object import SageObject
87
88
from sage.misc.search import search
89
90
_level_cache = {} # The info stored here is used in the normalization of MSymbols.
91
92
def P1NFList_clear_level_cache():
93
"""
94
Clear the global cache of data for the level ideals.
95
96
EXAMPLES::
97
98
sage: k.<a> = NumberField(x^3 + 11)
99
sage: N = k.ideal(a+1)
100
sage: alpha = MSymbol(N, 2*a^2, 5)
101
sage: alpha.normalize()
102
M-symbol (-4*a^2: 5*a^2) of level Fractional ideal (a + 1)
103
sage: sage.modular.modsym.p1list_nf._level_cache
104
{Fractional ideal (a + 1): (...)}
105
sage: sage.modular.modsym.p1list_nf.P1NFList_clear_level_cache()
106
sage: sage.modular.modsym.p1list_nf._level_cache
107
{}
108
"""
109
global _level_cache
110
_level_cache = {}
111
112
113
class MSymbol(SageObject):
114
"""
115
The constructor for an M-symbol over a number field.
116
117
INPUT:
118
119
- ``N`` -- integral ideal (the modulus or level).
120
121
- ``c`` -- integral element of the underlying number field or an MSymbol of
122
level N.
123
124
- ``d`` -- (optional) when present, it must be an integral element such
125
that <c> + <d> + N = R, where R is the corresponding ring of integers.
126
127
- ``check`` -- bool (default True). If ``check=False`` the constructor does
128
not check the condition <c> + <d> + N = R.
129
130
OUTPUT:
131
132
An M-symbol modulo the given ideal N, i.e. an element of the
133
projective line `\\mathbb{P}^1(R/N)`, where R is the ring of integers of
134
the underlying number field.
135
136
EXAMPLES::
137
138
sage: k.<a> = NumberField(x^3 + 11)
139
sage: N = k.ideal(a + 1, 2)
140
sage: MSymbol(N, 3, a^2 + 1)
141
M-symbol (3: a^2 + 1) of level Fractional ideal (2, a + 1)
142
143
We can give a tuple as input:
144
145
::
146
147
sage: MSymbol(N, (1, 0))
148
M-symbol (1: 0) of level Fractional ideal (2, a + 1)
149
150
We get an error if <c>, <d> and N are not coprime:
151
152
::
153
154
sage: MSymbol(N, 2*a, a - 1)
155
Traceback (most recent call last):
156
...
157
ValueError: (2*a, a - 1) is not an element of P1(R/N).
158
sage: MSymbol(N, (0, 0))
159
Traceback (most recent call last):
160
...
161
ValueError: (0, 0) is not an element of P1(R/N).
162
163
Saving and loading works:
164
165
::
166
167
sage: alpha = MSymbol(N, 3, a^2 + 1)
168
sage: loads(dumps(alpha))==alpha
169
True
170
"""
171
def __init__(self, N, c, d=None, check=True):
172
"""
173
See ``MSymbol`` for full documentation.
174
175
EXAMPLES::
176
177
sage: k.<a> = NumberField(x^4 + 13*x - 7)
178
sage: N = k.ideal(5)
179
sage: MSymbol(N, 0, 6*a)
180
M-symbol (0: 6*a) of level Fractional ideal (5)
181
sage: MSymbol(N, a^2 + 3, 7)
182
M-symbol (a^2 + 3: 7) of level Fractional ideal (5)
183
"""
184
k = N.number_field()
185
R = k.ring_of_integers()
186
self.__N = N
187
if d is None: # if we give a list (c, d) or an MSymbol as input
188
if isinstance(c, MSymbol):
189
if c.N() is N:
190
c1 = R(c[0])
191
d1 = R(c[1])
192
else:
193
raise ValueError, "Cannot change level of an MSymbol"
194
else:
195
try:
196
c1 = R(c[0])
197
d1 = R(c[1])
198
except (ValueError, TypeError):
199
raise TypeError, "Unable to create a Manin symbol from %s"%c
200
else:
201
try:
202
c1 = R(c)
203
d1 = R(d)
204
except (ValueError, TypeError):
205
raise TypeError, "Unable to create a Manin symbol from (%s, %s)"%(c, d)
206
if check:
207
if (c1.is_zero() and d1.is_zero()) or not N.is_coprime(k.ideal(c1, d1)):
208
raise ValueError, "(%s, %s) is not an element of P1(R/N)."%(c1, d1)
209
self.__c, self.__d = (c1, d1)
210
211
def __repr__(self):
212
"""
213
Returns the string representation of this MSymbol.
214
215
EXAMPLES::
216
217
sage: k.<a> = NumberField(x^2 + 23)
218
sage: N = k.ideal(3, a - 1)
219
sage: MSymbol(N, 3, a)
220
M-symbol (3: a) of level Fractional ideal (3, 1/2*a - 1/2)
221
"""
222
return "M-symbol (%s: %s) of level %s"%(self.__c, self.__d, self.__N)
223
224
def _latex_(self):
225
r"""
226
Return latex representation of self.
227
228
EXAMPLES::
229
230
sage: k.<a> = NumberField(x^4 + 13*x - 7)
231
sage: N = k.ideal(a^3 - 1)
232
sage: alpha = MSymbol(N, 3, 5*a^2 - 1)
233
sage: latex(alpha) # indirect doctest
234
\(3: 5 a^{2} - 1\)
235
"""
236
return "\\(%s: %s\)"%(self.c._latex_(), self.d._latex_())
237
238
def __cmp__(self, other):
239
"""
240
Comparison function for objects of the class MSymbol.
241
242
The order is the same as for the underlying lists of lists.
243
244
EXAMPLES::
245
246
sage: k.<a> = NumberField(x^2 + 23)
247
sage: N = k.ideal(3, a - 1)
248
sage: alpha = MSymbol(N, 3, a)
249
sage: beta = MSymbol(N, 1, 0)
250
sage: alpha < beta
251
False
252
sage: beta = MSymbol(N, 3, a + 1)
253
sage: alpha < beta
254
True
255
"""
256
if not isinstance(other, MSymbol):
257
raise ValueError, "You can only compare with another M-symbol"
258
return cmp([self.__c.list(), self.__d.list()],
259
[other.__c.list(), other.__d.list()])
260
261
def N(self):
262
"""
263
Returns the level or modulus of this MSymbol.
264
265
EXAMPLES::
266
267
sage: k.<a> = NumberField(x^2 + 23)
268
sage: N = k.ideal(3, a - 1)
269
sage: alpha = MSymbol(N, 3, a)
270
sage: alpha.N()
271
Fractional ideal (3, 1/2*a - 1/2)
272
"""
273
return self.__N
274
275
def tuple(self):
276
"""
277
Returns the MSymbol as a list (c, d).
278
279
EXAMPLES::
280
281
sage: k.<a> = NumberField(x^2 + 23)
282
sage: N = k.ideal(3, a - 1)
283
sage: alpha = MSymbol(N, 3, a); alpha
284
M-symbol (3: a) of level Fractional ideal (3, 1/2*a - 1/2)
285
sage: alpha.tuple()
286
(3, a)
287
"""
288
return self.__c, self.__d
289
290
def __getitem__(self, n):
291
"""
292
Indexing function for the list defined by an M-symbol.
293
294
INPUT:
295
296
- ``n`` -- integer (0 or 1, since the list defined by an M-symbol has
297
length 2)
298
299
EXAMPLES::
300
301
sage: k.<a> = NumberField(x^2 + 23)
302
sage: N = k.ideal(3, a - 1)
303
sage: alpha = MSymbol(N, 3, a); alpha
304
M-symbol (3: a) of level Fractional ideal (3, 1/2*a - 1/2)
305
sage: alpha[0]
306
3
307
sage: alpha[1]
308
a
309
"""
310
return self.tuple()[n]
311
312
def __get_c(self):
313
"""
314
Returns the first coefficient of the M-symbol.
315
316
EXAMPLES::
317
318
sage: k.<a> = NumberField(x^3 + 11)
319
sage: N = k.ideal(a + 1, 2)
320
sage: alpha = MSymbol(N, 3, a^2 + 1)
321
sage: alpha.c # indirect doctest
322
3
323
"""
324
return self.__c
325
c = property(__get_c)
326
327
def __get_d(self):
328
"""
329
Returns the second coefficient of the M-symbol.
330
331
EXAMPLES::
332
333
sage: k.<a> = NumberField(x^3 + 11)
334
sage: N = k.ideal(a + 1, 2)
335
sage: alpha = MSymbol(N, 3, a^2 + 1)
336
sage: alpha.d # indirect doctest
337
a^2 + 1
338
"""
339
return self.__d
340
d = property(__get_d)
341
342
def lift_to_sl2_Ok(self):
343
"""
344
Lift the MSymbol to an element of `SL(2, Ok)`, where `Ok` is the ring
345
of integers of the corresponding number field.
346
347
OUTPUT:
348
349
A list of integral elements `[a, b, c', d']` that are the entries of
350
a 2x2 matrix with determinant 1. The lower two entries are congruent
351
(modulo the level) to the coefficients `c, d` of the MSymbol self.
352
353
EXAMPLES::
354
355
sage: k.<a> = NumberField(x^2 + 23)
356
sage: N = k.ideal(3, a - 1)
357
sage: alpha = MSymbol(N, 3*a + 1, a)
358
sage: alpha.lift_to_sl2_Ok()
359
[0, -1, 1, a]
360
"""
361
return lift_to_sl2_Ok(self.__N, self.__c, self.__d)
362
363
def normalize(self, with_scalar=False):
364
"""
365
Returns a normalized MSymbol (a canonical representative of an element
366
of `\mathbb{P}^1(R/N)` ) equivalent to ``self``.
367
368
INPUT:
369
370
- ``with_scalar`` -- bool (default False)
371
372
OUTPUT:
373
374
- (only if ``with_scalar=True``) a transforming scalar `u`, such that
375
`(u*c', u*d')` is congruent to `(c: d)` (mod `N`), where `(c: d)`
376
are the coefficients of ``self`` and `N` is the level.
377
378
- a normalized MSymbol (c': d') equivalent to ``self``.
379
380
EXAMPLES::
381
382
sage: k.<a> = NumberField(x^2 + 23)
383
sage: N = k.ideal(3, a - 1)
384
sage: alpha1 = MSymbol(N, 3, a); alpha1
385
M-symbol (3: a) of level Fractional ideal (3, 1/2*a - 1/2)
386
sage: alpha1.normalize()
387
M-symbol (0: 1) of level Fractional ideal (3, 1/2*a - 1/2)
388
sage: alpha2 = MSymbol(N, 4, a + 1)
389
sage: alpha2.normalize()
390
M-symbol (1: -a) of level Fractional ideal (3, 1/2*a - 1/2)
391
392
We get the scaling factor by setting ``with_scalar=True``:
393
394
::
395
396
sage: alpha1.normalize(with_scalar=True)
397
(a, M-symbol (0: 1) of level Fractional ideal (3, 1/2*a - 1/2))
398
sage: r, beta1 = alpha1.normalize(with_scalar=True)
399
sage: r*beta1.c - alpha1.c in N
400
True
401
sage: r*beta1.d - alpha1.d in N
402
True
403
sage: r, beta2 = alpha2.normalize(with_scalar=True)
404
sage: r*beta2.c - alpha2.c in N
405
True
406
sage: r*beta2.d - alpha2.d in N
407
True
408
"""
409
N = self.__N
410
k = N.number_field()
411
R = k.ring_of_integers()
412
413
if self.__c in N:
414
if with_scalar:
415
return N.reduce(self.d), MSymbol(N, 0, 1)
416
else:
417
return MSymbol(N, 0, 1)
418
if self.d in N:
419
if with_scalar:
420
return N.reduce(self.c), MSymbol(N, 1, 0)
421
else:
422
return MSymbol(N, 1, 0)
423
if N.is_coprime(self.c):
424
cinv = R(self.c).inverse_mod(N)
425
if with_scalar:
426
return N.reduce(self.c), MSymbol(N, 1, N.reduce(self.d*cinv))
427
else:
428
return MSymbol(N, 1, N.reduce(self.d*cinv))
429
430
if _level_cache.has_key(N):
431
Lfacs, Lxs = _level_cache[N]
432
else:
433
Lfacs = [p**e for p, e in N.factor()]
434
Lxs = [(N/p).element_1_mod(p) for p in Lfacs]
435
# Lfacs, Lxs only depend of the ideal: same lists every time we
436
# call normalize for a given level, so we store the lists.
437
_level_cache[N] = (Lfacs, Lxs)
438
u = 0 # normalizer factor
439
p_i = 0
440
for p in Lfacs:
441
if p.is_coprime(self.c):
442
inv = self.c.inverse_mod(p)
443
else:
444
inv = self.d.inverse_mod(p)
445
u = u + inv*Lxs[p_i]
446
p_i = p_i + 1
447
c, d = (N.reduce(u*self.c), N.reduce(u*self.d))
448
if (c - 1) in N:
449
c = R(1)
450
if with_scalar:
451
return u.inverse_mod(N), MSymbol(N, c, d)
452
else:
453
return MSymbol(N, c, d)
454
455
456
#**************************************************************************
457
#* P1NFList class *
458
#**************************************************************************
459
class P1NFList(SageObject):
460
"""
461
The class for `\mathbb{P}^1(R/N)`, the projective line modulo `N`, where
462
`R` is the ring of integers of a number field `K` and `N` is an integral ideal.
463
464
INPUT:
465
466
- ``N`` - integral ideal (the modulus or level).
467
468
OUTPUT:
469
470
A P1NFList object representing `\mathbb{P}^1(R/N)`.
471
472
EXAMPLES::
473
474
sage: k.<a> = NumberField(x^3 + 11)
475
sage: N = k.ideal(5, a + 1)
476
sage: P = P1NFList(N); P
477
The projective line over the ring of integers modulo the Fractional ideal (5, a + 1)
478
479
Saving and loading works.
480
481
::
482
483
sage: loads(dumps(P)) == P
484
True
485
"""
486
def __init__(self, N):
487
"""
488
The constructor for the class P1NFList. See ``P1NFList`` for full
489
documentation.
490
491
EXAMPLES::
492
493
sage: k.<a> = NumberField(x^2 + 5)
494
sage: N = k.ideal(3, a - 1)
495
sage: P = P1NFList(N); P
496
The projective line over the ring of integers modulo the Fractional ideal (3, a + 2)
497
"""
498
self.__N = N
499
self.__list = p1NFlist(N)
500
self.__list.sort()
501
502
def __cmp__(self, other):
503
"""
504
Comparison function for objects of the class P1NFList.
505
506
The order is the same as for the underlying modulus.
507
508
EXAMPLES::
509
510
sage: k.<a> = NumberField(x^2 + 23)
511
sage: N1 = k.ideal(3, a + 1)
512
sage: P1 = P1NFList(N1)
513
sage: N2 = k.ideal(a + 2)
514
sage: P2 = P1NFList(N2)
515
sage: P1 < P2
516
False
517
sage: P1 > P2
518
True
519
sage: P1 == P1NFList(N1)
520
True
521
"""
522
if not isinstance(other, P1NFList):
523
raise ValueError, "You can only compare with another P1NFList"
524
return cmp(self.__N, other.__N)
525
526
def __getitem__(self, n):
527
"""
528
Standard indexing function for the class P1NFList.
529
530
EXAMPLES::
531
532
sage: k.<a> = NumberField(x^3 + 11)
533
sage: N = k.ideal(a)
534
sage: P = P1NFList(N)
535
sage: list(P) == P._P1NFList__list
536
True
537
sage: j = randint(0,len(P)-1)
538
sage: P[j] == P._P1NFList__list[j]
539
True
540
"""
541
return self.__list[n]
542
543
def __len__(self):
544
"""
545
Returns the length of this P1NFList.
546
547
EXAMPLES::
548
549
sage: k.<a> = NumberField(x^3 + 11)
550
sage: N = k.ideal(5, a^2 - a + 1)
551
sage: P = P1NFList(N)
552
sage: len(P)
553
26
554
"""
555
return len(self.__list)
556
557
def __repr__(self):
558
"""
559
Returns the string representation of this P1NFList.
560
561
EXAMPLES::
562
563
sage: k.<a> = NumberField(x^3 + 11)
564
sage: N = k.ideal(5, a+1)
565
sage: P = P1NFList(N); P
566
The projective line over the ring of integers modulo the Fractional ideal (5, a + 1)
567
568
"""
569
return "The projective line over the ring of integers modulo the %s"%self.__N
570
571
def list(self):
572
"""
573
Returns the underlying list of this P1NFList object.
574
575
EXAMPLES::
576
577
sage: k.<a> = NumberField(x^3 + 11)
578
sage: N = k.ideal(5, a+1)
579
sage: P = P1NFList(N)
580
sage: type(P)
581
<class 'sage.modular.modsym.p1list_nf.P1NFList'>
582
sage: type(P.list())
583
<type 'list'>
584
"""
585
return self.__list
586
587
def normalize(self, c, d=None, with_scalar=False):
588
"""
589
Returns a normalised element of `\mathbb{P}^1(R/N)`.
590
591
INPUT:
592
593
- ``c`` -- integral element of the underlying number field, or an
594
MSymbol.
595
596
- ``d`` -- (optional) when present, it must be an integral element of
597
the number field such that `(c, d)` defines an M-symbol of level `N`.
598
599
- ``with_scalar`` -- bool (default False)
600
601
OUTPUT:
602
603
- (only if ``with_scalar=True``) a transforming scalar `u`, such that
604
`(u*c', u*d')` is congruent to `(c: d)` (mod `N`).
605
606
- a normalized MSymbol (c': d') equivalent to `(c: d)`.
607
608
EXAMPLES::
609
610
sage: k.<a> = NumberField(x^2 + 31)
611
sage: N = k.ideal(5, a + 3)
612
sage: P = P1NFList(N)
613
sage: P.normalize(3, a)
614
M-symbol (1: 2*a) of level Fractional ideal (5, 1/2*a + 3/2)
615
616
We can use an MSymbol as input:
617
618
::
619
620
sage: alpha = MSymbol(N, 3, a)
621
sage: P.normalize(alpha)
622
M-symbol (1: 2*a) of level Fractional ideal (5, 1/2*a + 3/2)
623
624
If we are interested in the normalizing scalar:
625
626
::
627
628
sage: P.normalize(alpha, with_scalar=True)
629
(-a, M-symbol (1: 2*a) of level Fractional ideal (5, 1/2*a + 3/2))
630
sage: r, beta = P.normalize(alpha, with_scalar=True)
631
sage: (r*beta.c - alpha.c in N) and (r*beta.d - alpha.d in N)
632
True
633
"""
634
if d is None:
635
try:
636
c = MSymbol(self.__N, c) # check that c is an MSymbol
637
except ValueError: # catch special case of wrong level
638
raise ValueError, "The MSymbol is of a different level"
639
return c.normalize(with_scalar)
640
return MSymbol(self.N(), c, d).normalize(with_scalar)
641
642
def N(self):
643
"""
644
Returns the level or modulus of this P1NFList.
645
646
EXAMPLES::
647
648
sage: k.<a> = NumberField(x^2 + 31)
649
sage: N = k.ideal(5, a + 3)
650
sage: P = P1NFList(N)
651
sage: P.N()
652
Fractional ideal (5, 1/2*a + 3/2)
653
"""
654
return self.__N
655
656
def index(self, c, d=None, with_scalar=False):
657
"""
658
Returns the index of the class of the pair `(c, d)` in the fixed list
659
of representatives of `\mathbb{P}^1(R/N)`.
660
661
INPUT:
662
663
- ``c`` -- integral element of the corresponding number field, or an
664
MSymbol.
665
666
- ``d`` -- (optional) when present, it must be an integral element of
667
the number field such that `(c, d)` defines an M-symbol of level `N`.
668
669
- ``with_scalar`` -- bool (default False)
670
671
OUTPUT:
672
673
- ``u`` - the normalizing scalar (only if ``with_scalar=True``)
674
675
- ``i`` - the index of `(c, d)` in the list.
676
677
EXAMPLES::
678
679
sage: k.<a> = NumberField(x^2 + 31)
680
sage: N = k.ideal(5, a + 3)
681
sage: P = P1NFList(N)
682
sage: P.index(3,a)
683
5
684
sage: P[5]==MSymbol(N, 3, a).normalize()
685
True
686
687
We can give an MSymbol as input:
688
689
::
690
691
sage: alpha = MSymbol(N, 3, a)
692
sage: P.index(alpha)
693
5
694
695
We cannot look for the class of an MSymbol of a different level:
696
697
::
698
699
sage: M = k.ideal(a + 1)
700
sage: beta = MSymbol(M, 0, 1)
701
sage: P.index(beta)
702
Traceback (most recent call last):
703
...
704
ValueError: The MSymbol is of a different level
705
706
If we are interested in the transforming scalar:
707
708
::
709
710
sage: alpha = MSymbol(N, 3, a)
711
sage: P.index(alpha, with_scalar=True)
712
(-a, 5)
713
sage: u, i = P.index(alpha, with_scalar=True)
714
sage: (u*P[i].c - alpha.c in N) and (u*P[i].d - alpha.d in N)
715
True
716
"""
717
if d is None:
718
try:
719
c = MSymbol(self.__N, c) # check that c is an MSymbol
720
except ValueError: # catch special case of wrong level
721
raise ValueError, "The MSymbol is of a different level"
722
if with_scalar:
723
u, norm_c = c.normalize(with_scalar=True)
724
else:
725
norm_c = c.normalize()
726
else:
727
if with_scalar:
728
u, norm_c = MSymbol(self.__N, c, d).normalize(with_scalar=True)
729
else:
730
norm_c = MSymbol(self.__N, c, d).normalize()
731
t, i = search(self.__list, norm_c)
732
if t:
733
if with_scalar:
734
return u, i
735
else:
736
return i
737
return False
738
739
def index_of_normalized_pair(self, c, d=None):
740
"""
741
Returns the index of the class `(c, d)` in the fixed list of
742
representatives of `\mathbb(P)^1(R/N)`.
743
744
INPUT:
745
746
- ``c`` -- integral element of the corresponding number field, or a
747
normalized MSymbol.
748
749
- ``d`` -- (optional) when present, it must be an integral element of
750
the number field such that `(c, d)` defines a normalized M-symbol of
751
level `N`.
752
753
OUTPUT:
754
755
- ``i`` - the index of `(c, d)` in the list.
756
757
EXAMPLES::
758
759
sage: k.<a> = NumberField(x^2 + 31)
760
sage: N = k.ideal(5, a + 3)
761
sage: P = P1NFList(N)
762
sage: P.index_of_normalized_pair(1, 0)
763
3
764
sage: j = randint(0,len(P)-1)
765
sage: P.index_of_normalized_pair(P[j])==j
766
True
767
"""
768
if d is None:
769
try:
770
c = MSymbol(self.__N, c) # check that c is an MSymbol
771
except ValueError: # catch special case of wrong level
772
raise ValueError, "The MSymbol is of a different level"
773
t, i = search(self.__list, c)
774
else:
775
t, i = search(self.__list, MSymbol(self.__N, c, d))
776
if t: return i
777
return False
778
779
def lift_to_sl2_Ok(self, i):
780
"""
781
Lift the `i`-th element of this P1NFList to an element of `SL(2, R)`,
782
where `R` is the ring of integers of the corresponding number field.
783
784
INPUT:
785
786
- ``i`` - integer (index of the element to lift)
787
788
OUTPUT:
789
790
If the `i`-th element is `(c : d)`, the function returns a list of
791
integral elements `[a, b, c', d']` that defines a 2x2 matrix with
792
determinant 1 and such that `c=c'` (mod `N`) and `d=d'` (mod `N`).
793
794
EXAMPLES::
795
796
sage: k.<a> = NumberField(x^2 + 23)
797
sage: N = k.ideal(3)
798
sage: P = P1NFList(N)
799
sage: len(P)
800
16
801
sage: P[5]
802
M-symbol (1/2*a + 1/2: -a) of level Fractional ideal (3)
803
sage: P.lift_to_sl2_Ok(5)
804
[1, -2, 1/2*a + 1/2, -a]
805
806
::
807
808
sage: Ok = k.ring_of_integers()
809
sage: L = [Matrix(Ok, 2, P.lift_to_sl2_Ok(i)) for i in range(len(P))]
810
sage: all([det(L[i]) == 1 for i in range(len(L))])
811
True
812
"""
813
return self[i].lift_to_sl2_Ok()
814
815
def apply_S(self, i):
816
"""
817
Applies the matrix S = [0, -1, 1, 0] to the i-th M-Symbol of the list.
818
819
INPUT:
820
821
- ``i`` -- integer
822
823
OUTPUT:
824
825
integer -- the index of the M-Symbol obtained by the right action of
826
the matrix S = [0, -1, 1, 0] on the i-th M-Symbol.
827
828
EXAMPLES::
829
830
sage: k.<a> = NumberField(x^3 + 11)
831
sage: N = k.ideal(5, a + 1)
832
sage: P = P1NFList(N)
833
sage: j = P.apply_S(P.index_of_normalized_pair(1, 0))
834
sage: P[j]
835
M-symbol (0: 1) of level Fractional ideal (5, a + 1)
836
837
We test that S has order 2:
838
839
::
840
841
sage: j = randint(0,len(P)-1)
842
sage: P.apply_S(P.apply_S(j))==j
843
True
844
"""
845
c, d = self.__list[i].tuple()
846
t, j = search(self.__list, self.normalize(d, -c))
847
return j
848
849
def apply_TS(self, i):
850
"""
851
Applies the matrix TS = [1, -1, 0, 1] to the i-th M-Symbol of the list.
852
853
INPUT:
854
855
- ``i`` -- integer
856
857
OUTPUT:
858
859
integer -- the index of the M-Symbol obtained by the right action of
860
the matrix TS = [1, -1, 0, 1] on the i-th M-Symbol.
861
862
EXAMPLES::
863
864
sage: k.<a> = NumberField(x^3 + 11)
865
sage: N = k.ideal(5, a + 1)
866
sage: P = P1NFList(N)
867
sage: P.apply_TS(3)
868
2
869
870
We test that TS has order 3:
871
872
::
873
874
sage: j = randint(0,len(P)-1)
875
sage: P.apply_TS(P.apply_TS(P.apply_TS(j)))==j
876
True
877
"""
878
c, d = self.__list[i].tuple()
879
t, j = search(self.__list, self.normalize(c + d, -c))
880
return j
881
882
def apply_T_alpha(self, i, alpha=1):
883
"""
884
Applies the matrix T_alpha = [1, alpha, 0, 1] to the i-th M-Symbol of
885
the list.
886
887
INPUT:
888
889
- ``i`` -- integer
890
891
- ``alpha`` -- element of the corresponding ring of integers(default 1)
892
893
OUTPUT:
894
895
integer -- the index of the M-Symbol obtained by the right action of
896
the matrix T_alpha = [1, alpha, 0, 1] on the i-th M-Symbol.
897
898
EXAMPLES::
899
900
sage: k.<a> = NumberField(x^3 + 11)
901
sage: N = k.ideal(5, a + 1)
902
sage: P = P1NFList(N)
903
sage: P.apply_T_alpha(4, a^ 2 - 2)
904
3
905
906
We test that T_a*T_b = T_(a+b):
907
908
::
909
910
sage: P.apply_T_alpha(3, a^2 - 2)==P.apply_T_alpha(P.apply_T_alpha(3,a^2),-2)
911
True
912
"""
913
c, d = self.__list[i].tuple()
914
t, j = search(self.__list, self.normalize(c, alpha*c + d))
915
return j
916
917
def apply_J_epsilon(self, i, e1, e2=1):
918
"""
919
Applies the matrix `J_{\epsilon}` = [e1, 0, 0, e2] to the i-th
920
M-Symbol of the list.
921
922
e1, e2 are units of the underlying number field.
923
924
INPUT:
925
926
- ``i`` -- integer
927
928
- ``e1`` -- unit
929
930
- ``e2`` -- unit (default 1)
931
932
OUTPUT:
933
934
integer -- the index of the M-Symbol obtained by the right action of
935
the matrix `J_{\epsilon}` = [e1, 0, 0, e2] on the i-th M-Symbol.
936
937
EXAMPLES::
938
939
sage: k.<a> = NumberField(x^3 + 11)
940
sage: N = k.ideal(5, a + 1)
941
sage: P = P1NFList(N)
942
sage: u = k.unit_group().gens(); u
943
[-1, 2*a^2 + 4*a - 1]
944
sage: P.apply_J_epsilon(4, -1)
945
2
946
sage: P.apply_J_epsilon(4, u[0], u[1])
947
1
948
949
::
950
951
sage: k.<a> = NumberField(x^4 + 13*x - 7)
952
sage: N = k.ideal(a + 1)
953
sage: P = P1NFList(N)
954
sage: u = k.unit_group().gens(); u
955
[-1, a^3 + a^2 + a + 12, a^3 + 3*a^2 - 1]
956
sage: P.apply_J_epsilon(3, u[2]^2)==P.apply_J_epsilon(P.apply_J_epsilon(3, u[2]),u[2])
957
True
958
"""
959
c, d = self.__list[i].tuple()
960
t, j = search(self.__list, self.normalize(c*e1, d*e2))
961
return j
962
963
964
#**************************************************************************
965
# Global functions:
966
# - p1NFList --compute list of M-symbols
967
# - lift_to_sl2_Ok
968
# - make_coprime -- need it for ``lift_to_sl2_Ok``
969
# - psi -- useful to check cardinality of the M-symbols list
970
#**************************************************************************
971
972
def p1NFlist(N):
973
"""
974
Returns a list of the normalized elements of `\\mathbb{P}^1(R/N)`, where
975
`N` is an integral ideal.
976
977
INPUT:
978
979
- ``N`` - integral ideal (the level or modulus).
980
981
EXAMPLES::
982
983
sage: k.<a> = NumberField(x^2 + 23)
984
sage: N = k.ideal(3)
985
sage: from sage.modular.modsym.p1list_nf import p1NFlist, psi
986
sage: len(p1NFlist(N))==psi(N)
987
True
988
"""
989
k = N.number_field()
990
991
L = [MSymbol(N, k(0),k(1), check=False)]
992
#N.residues() = iterator through the residues mod N
993
L = L+[MSymbol(N, k(1), r, check=False) for r in N.residues()]
994
995
from sage.rings.arith import divisors
996
for D in divisors(N):
997
if not D.is_trivial() and D!=N:
998
#we find Dp ideal coprime to N, in inverse class to D
999
if D.is_principal():
1000
Dp = k.ideal(1)
1001
c = D.gens_reduced()[0]
1002
else:
1003
it = k.primes_of_degree_one_iter()
1004
Dp = it.next()
1005
while not Dp.is_coprime(N) or not (Dp*D).is_principal():
1006
Dp = it.next()
1007
c = (D*Dp).gens_reduced()[0]
1008
#now we find all the (c,d)'s which have associated divisor D
1009
I = D + N/D
1010
for d in (N/D).residues():
1011
if I.is_coprime(d):
1012
M = D.prime_to_idealM_part(N/D)
1013
u = (Dp*M).element_1_mod(N/D)
1014
d1 = u*d + (1-u)
1015
L.append(MSymbol(N, c, d1, check=False).normalize())
1016
return L
1017
1018
def lift_to_sl2_Ok(N, c, d):
1019
"""
1020
Lift a pair (c, d) to an element of `SL(2, O_k)`, where `O_k` is the ring
1021
of integers of the corresponding number field.
1022
1023
INPUT:
1024
1025
- ``N`` -- number field ideal
1026
1027
- ``c`` -- integral element of the number field
1028
1029
- ``d`` -- integral element of the number field
1030
1031
OUTPUT:
1032
1033
A list [a, b, c', d'] of integral elements that are the entries of
1034
a 2x2 matrix with determinant 1. The lower two entries are congruent to
1035
c, d modulo the ideal `N`.
1036
1037
1038
EXAMPLES::
1039
1040
sage: from sage.modular.modsym.p1list_nf import lift_to_sl2_Ok
1041
sage: k.<a> = NumberField(x^2 + 23)
1042
sage: Ok = k.ring_of_integers(k)
1043
sage: N = k.ideal(3)
1044
sage: M = Matrix(Ok, 2, lift_to_sl2_Ok(N, 1, a))
1045
sage: det(M)
1046
1
1047
sage: M = Matrix(Ok, 2, lift_to_sl2_Ok(N, 0, a))
1048
sage: det(M)
1049
1
1050
sage: (M[1][0] in N) and (M[1][1] - a in N)
1051
True
1052
sage: M = Matrix(Ok, 2, lift_to_sl2_Ok(N, 0, 0))
1053
Traceback (most recent call last):
1054
...
1055
ValueError: Cannot lift (0, 0) to an element of Sl2(Ok).
1056
1057
::
1058
1059
sage: k.<a> = NumberField(x^3 + 11)
1060
sage: Ok = k.ring_of_integers(k)
1061
sage: N = k.ideal(3, a - 1)
1062
sage: M = Matrix(Ok, 2, lift_to_sl2_Ok(N, 2*a, 0))
1063
sage: det(M)
1064
1
1065
sage: (M[1][0] - 2*a in N) and (M[1][1] in N)
1066
True
1067
sage: M = Matrix(Ok, 2, lift_to_sl2_Ok(N, 4*a^2, a + 1))
1068
sage: det(M)
1069
1
1070
sage: (M[1][0] - 4*a^2 in N) and (M[1][1] - (a+1) in N)
1071
True
1072
1073
::
1074
1075
sage: k.<a> = NumberField(x^4 - x^3 -21*x^2 + 17*x + 133)
1076
sage: Ok = k.ring_of_integers(k)
1077
sage: N = k.ideal(7, a)
1078
sage: M = Matrix(Ok, 2, lift_to_sl2_Ok(N, 0, a^2 - 1))
1079
sage: det(M)
1080
1
1081
sage: (M[1][0] in N) and (M[1][1] - (a^2-1) in N)
1082
True
1083
sage: M = Matrix(Ok, 2, lift_to_sl2_Ok(N, 0, 7))
1084
Traceback (most recent call last):
1085
...
1086
ValueError: <0> + <7> and the Fractional ideal (7, a) are not coprime.
1087
"""
1088
k = N.number_field()
1089
#check the input
1090
if c.is_zero() and d.is_zero():
1091
raise ValueError, "Cannot lift (%s, %s) to an element of Sl2(Ok)."%(c, d)
1092
if not N.is_coprime(k.ideal(c, d)):
1093
raise ValueError, "<%s> + <%s> and the %s are not coprime."%(c, d, N)
1094
#a few special cases
1095
if c - 1 in N:
1096
return [k(0), k(-1), 1, d]
1097
if d - 1 in N:
1098
return [k(1), k(0), c, 1]
1099
if c.is_zero(): # and d!=1, so won't happen for normalized M-symbols (c: d)
1100
it = k.primes_of_degree_one_iter()
1101
q = k.ideal(1)
1102
while not (q.is_coprime(d) and (q*N).is_principal()):
1103
q = it.next()
1104
m = (q*N).gens_reduced()[0]
1105
B = k.ideal(m).element_1_mod(k.ideal(d))
1106
return [(1-B)/d, -B/m, m, d]
1107
if d.is_zero(): # and c!=1, so won't happen for normalized M-symbols (c: d)
1108
it = k.primes_of_degree_one_iter()
1109
q = k.ideal(1)
1110
while not (q.is_coprime(c) and (q*N).is_principal()):
1111
q = it.next()
1112
m = (q*N).gens_reduced()[0]
1113
B = k.ideal(c).element_1_mod(k.ideal(m))
1114
return [(1-B)/m, -B/c, c, m]
1115
1116
c, d = make_coprime(N, c, d)
1117
1118
B = k.ideal(c).element_1_mod(k.ideal(d))
1119
b = -B/c
1120
a = (1-B)/d
1121
return [a, b, c, d]
1122
1123
def make_coprime(N, c, d):
1124
"""
1125
Returns (c, d') so d' is congruent to d modulo N, and such that c and d' are
1126
coprime (<c> + <d'> = R).
1127
1128
INPUT:
1129
1130
- ``N`` -- number field ideal
1131
1132
- ``c`` -- integral element of the number field
1133
1134
- ``d`` -- integral element of the number field
1135
1136
OUTPUT:
1137
1138
A pair (c, d') where c, d' are integral elements of the corresponding
1139
number field, with d' congruent to d mod N, and such that <c> + <d'> = R
1140
(R being the corresponding ring of integers).
1141
1142
EXAMPLES::
1143
1144
sage: from sage.modular.modsym.p1list_nf import make_coprime
1145
sage: k.<a> = NumberField(x^2 + 23)
1146
sage: N = k.ideal(3, a - 1)
1147
sage: c = 2*a; d = a + 1
1148
sage: N.is_coprime(k.ideal(c, d))
1149
True
1150
sage: k.ideal(c).is_coprime(d)
1151
False
1152
sage: c, dp = make_coprime(N, c, d)
1153
sage: k.ideal(c).is_coprime(dp)
1154
True
1155
"""
1156
k = N.number_field()
1157
if k.ideal(c).is_coprime(d):
1158
return c, d
1159
else:
1160
q = k.ideal(c).prime_to_idealM_part(d)
1161
it = k.primes_of_degree_one_iter()
1162
r = k.ideal(1)
1163
qN = q*N
1164
while not (r.is_coprime(c) and (r*qN).is_principal()):
1165
r = it.next()
1166
m = (r*qN).gens_reduced()[0]
1167
d1 = d + m
1168
return c, d1
1169
1170
def psi(N):
1171
"""
1172
The index `[\Gamma : \Gamma_0(N)]`, where `\Gamma = GL(2, R)` for `R` the
1173
corresponding ring of integers, and `\Gamma_0(N)` standard congruence
1174
subgroup.
1175
1176
EXAMPLES::
1177
1178
sage: from sage.modular.modsym.p1list_nf import psi
1179
sage: k.<a> = NumberField(x^2 + 23)
1180
sage: N = k.ideal(3, a - 1)
1181
sage: psi(N)
1182
4
1183
1184
::
1185
1186
sage: k.<a> = NumberField(x^2 + 23)
1187
sage: N = k.ideal(5)
1188
sage: psi(N)
1189
26
1190
"""
1191
if not N.is_integral():
1192
raise ValueError, "psi only defined for integral ideals"
1193
1194
from sage.misc.misc import prod
1195
return prod([(np+1)*np**(e-1) \
1196
for np,e in [(p.absolute_norm(),e) \
1197
for p,e in N.factor()]])
1198
1199