Book a Demo!
CoCalc Logo Icon
StoreFeaturesDocsShareSupportNewsAboutPoliciesSign UpSign In
sagemath
GitHub Repository: sagemath/sagesmc
Path: blob/master/src/sage/rings/finite_rings/element_base.pyx
8820 views
1
"""
2
Base class for finite field elements
3
4
AUTHORS::
5
6
- David Roe (2010-1-14) -- factored out of sage.structure.element
7
8
"""
9
include "sage/ext/stdsage.pxi"
10
11
from sage.structure.element cimport Element
12
from sage.structure.parent cimport Parent
13
from sage.rings.integer import Integer
14
15
def is_FiniteFieldElement(x):
16
"""
17
Returns if x is a finite field element.
18
19
EXAMPLE::
20
21
sage: from sage.rings.finite_rings.element_ext_pari import is_FiniteFieldElement
22
sage: is_FiniteFieldElement(1)
23
False
24
sage: is_FiniteFieldElement(IntegerRing())
25
False
26
sage: is_FiniteFieldElement(GF(5)(2))
27
True
28
"""
29
from sage.rings.finite_rings.finite_field_base import is_FiniteField
30
return isinstance(x, Element) and is_FiniteField(x.parent())
31
32
cdef class FiniteRingElement(CommutativeRingElement):
33
def _nth_root_common(self, n, all, algorithm, cunningham):
34
"""
35
This function exists to reduce code duplication between finite field
36
nth roots and integer_mod nth roots.
37
38
The inputs are described there.
39
40
TESTS::
41
42
sage: a = Zmod(17)(13)
43
sage: a._nth_root_common(4, True, "Johnston", False)
44
[3, 5, 14, 12]
45
sage: a._nth_root_common(4, True, "Johnston", cunningham = True) # optional - cunningham
46
[3, 5, 14, 12]
47
"""
48
K = self.parent()
49
q = K.order()
50
if self.is_one():
51
gcd = n.gcd(q-1)
52
if gcd == 1:
53
if all: return [self]
54
else: return self
55
else:
56
# the following may eventually be improved to not need a multiplicative generator.
57
g = K.multiplicative_generator()
58
q1overn = (q-1)//gcd
59
nthroot = g**q1overn
60
return [nthroot**a for a in range(gcd)] if all else nthroot
61
n = n % (q-1)
62
if n == 0:
63
if all: return []
64
else: raise ValueError, "no nth root"
65
gcd, alpha, beta = n.xgcd(q-1) # gcd = alpha*n + beta*(q-1), so 1/n = alpha/gcd (mod q-1)
66
if gcd == 1:
67
return [self**alpha] if all else self**alpha
68
n = gcd
69
q1overn = (q-1)//n
70
if self**q1overn != 1:
71
if all: return []
72
else: raise ValueError, "no nth root"
73
self = self**alpha
74
if cunningham:
75
from sage.rings.factorint import factor_cunningham
76
F = factor_cunningham(n)
77
else:
78
F = n.factor()
79
from sage.groups.generic import discrete_log
80
if algorithm is None or algorithm == 'Johnston':
81
g = K.multiplicative_generator()
82
for r, v in F:
83
k, h = (q-1).val_unit(r)
84
z = h * (-h).inverse_mod(r**v)
85
x = (1 + z) // r**v
86
if k == 1:
87
self = self**x
88
else:
89
t = discrete_log(self**h, g**(r**v*h), r**(k-v), operation='*')
90
self = self**x * g**(-z*t)
91
if all:
92
nthroot = g**q1overn
93
L = [self]
94
for i in range(1,n):
95
self *= nthroot
96
L.append(self)
97
return L
98
else:
99
return self
100
else:
101
raise ValueError, "unknown algorithm"
102
103
cdef class FinitePolyExtElement(FiniteRingElement):
104
"""
105
Elements represented as polynomials modulo a given ideal.
106
107
TESTS::
108
109
sage: k.<a> = GF(64)
110
sage: TestSuite(a).run()
111
"""
112
def _im_gens_(self, codomain, im_gens):
113
"""
114
Used for applying homomorphisms of finite fields.
115
116
EXAMPLES::
117
118
sage: k.<a> = FiniteField(73^2, 'a')
119
sage: K.<b> = FiniteField(73^4, 'b')
120
sage: phi = k.hom([ b^(73*73+1) ]) # indirect doctest
121
sage: phi(0)
122
0
123
sage: phi(a)
124
7*b^3 + 13*b^2 + 65*b + 71
125
126
sage: phi(a+3)
127
7*b^3 + 13*b^2 + 65*b + 1
128
"""
129
## NOTE: see the note in sage/rings/number_field_element.pyx,
130
## in the comments for _im_gens_ there -- something analogous
131
## applies here.
132
return codomain(self.polynomial()(im_gens[0]))
133
134
def minpoly(self,var='x'):
135
"""
136
Returns the minimal polynomial of this element
137
(over the corresponding prime subfield).
138
139
EXAMPLES::
140
141
sage: k.<a> = FiniteField(19^2)
142
sage: parent(a)
143
Finite Field in a of size 19^2
144
sage: b=a**20;p=b.charpoly("x");p
145
x^2 + 15*x + 4
146
sage: factor(p)
147
(x + 17)^2
148
sage: b.minpoly('x')
149
x + 17
150
"""
151
p=self.charpoly(var);
152
for q in p.factor():
153
if q[0](self)==0:
154
return q[0]
155
# This shouldn't be reached, but you never know!
156
raise ArithmeticError("Could not find the minimal polynomial")
157
158
## We have two names for the same method
159
## for compatibility with sage.matrix
160
def minimal_polynomial(self,var='x'):
161
"""
162
Returns the minimal polynomial of this element
163
(over the corresponding prime subfield).
164
165
EXAMPLES::
166
167
sage: k.<a> = FiniteField(3^4)
168
sage: parent(a)
169
Finite Field in a of size 3^4
170
sage: b=a**20;p=charpoly(b,"y");p
171
y^4 + 2*y^2 + 1
172
sage: factor(p)
173
(y^2 + 1)^2
174
sage: b.minimal_polynomial('y')
175
y^2 + 1
176
"""
177
return self.minpoly(var)
178
179
def vector(self, reverse=False):
180
r"""
181
See :meth:`_vector_`.
182
183
EXAMPLE::
184
185
sage: k.<a> = GF(2^16)
186
sage: e = a^2 + 1
187
sage: e.vector() # random-ish error message
188
doctest:1: DeprecationWarning:The function vector is replaced by _vector_.
189
(1, 0, 1, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0)
190
"""
191
from sage.misc.superseded import deprecation
192
deprecation(8218, "The function vector is replaced by _vector_.")
193
return self._vector_()
194
195
def _vector_(self, reverse=False):
196
"""
197
Return a vector in self.parent().vector_space() matching
198
self. The most significant bit is to the right.
199
200
INPUT::
201
202
- ``reverse`` -- reverse the order of the bits
203
from little endian to big endian.
204
205
EXAMPLES::
206
207
sage: k.<a> = GF(2^16)
208
sage: e = a^2 + 1
209
sage: v = vector(e)
210
sage: v
211
(1, 0, 1, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0)
212
sage: k(v)
213
a^2 + 1
214
215
sage: k.<a> = GF(3^16)
216
sage: e = 2*a^2 + 1
217
sage: v = vector(e)
218
sage: v
219
(1, 0, 2, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0)
220
sage: k(v)
221
2*a^2 + 1
222
223
You can also compute the vector in the other order::
224
225
sage: e._vector_(reverse=True)
226
(0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 2, 0, 1)
227
"""
228
#vector(foo) might pass in ZZ
229
if PY_TYPE_CHECK(reverse, Parent):
230
raise TypeError, "Base field is fixed to prime subfield."
231
232
k = self.parent()
233
234
v = self.polynomial().list()
235
236
ret = [v[i] for i in range(len(v))]
237
238
for i in range(k.degree() - len(ret)):
239
ret.append(0)
240
241
if reverse:
242
ret = list(reversed(ret))
243
return k.vector_space()(ret)
244
245
def matrix(self, reverse=False):
246
r"""
247
See :meth:`_matrix_`.
248
249
EXAMPLE::
250
251
sage: k.<a> = GF(2^16)
252
sage: e = a^2 + 1
253
sage: e.matrix() # random-ish error message
254
doctest:1: DeprecationWarning:The function matrix is replaced by _matrix_.
255
[1 0 0 0 0 0 0 0 0 0 0 0 0 0 1 0]
256
[0 1 0 0 0 0 0 0 0 0 0 0 0 0 0 1]
257
[1 0 1 0 0 0 0 0 0 0 0 0 0 0 1 0]
258
[0 1 0 1 0 0 0 0 0 0 0 0 0 0 1 1]
259
[0 0 1 0 1 0 0 0 0 0 0 0 0 0 0 1]
260
[0 0 0 1 0 1 0 0 0 0 0 0 0 0 1 0]
261
[0 0 0 0 1 0 1 0 0 0 0 0 0 0 0 1]
262
[0 0 0 0 0 1 0 1 0 0 0 0 0 0 0 0]
263
[0 0 0 0 0 0 1 0 1 0 0 0 0 0 0 0]
264
[0 0 0 0 0 0 0 1 0 1 0 0 0 0 0 0]
265
[0 0 0 0 0 0 0 0 1 0 1 0 0 0 0 0]
266
[0 0 0 0 0 0 0 0 0 1 0 1 0 0 0 0]
267
[0 0 0 0 0 0 0 0 0 0 1 0 1 0 0 0]
268
[0 0 0 0 0 0 0 0 0 0 0 1 0 1 0 0]
269
[0 0 0 0 0 0 0 0 0 0 0 0 1 0 1 0]
270
[0 0 0 0 0 0 0 0 0 0 0 0 0 1 0 1]
271
"""
272
from sage.misc.superseded import deprecation
273
deprecation(8218, "The function matrix is replaced by _matrix_.")
274
return self._matrix_()
275
276
277
def _matrix_(self, reverse=False):
278
"""
279
Return the matrix of right multiplication by the element on
280
the power basis `1, x, x^2, \ldots, x^{d-1}` for the field
281
extension. Thus the \emph{rows} of this matrix give the images
282
of each of the `x^i`.
283
284
INPUT:
285
286
- ``reverse`` - if True act on vectors in reversed order
287
288
EXAMPLE::
289
290
sage: k.<a> = GF(2^4)
291
sage: a._vector_(reverse=True), a._matrix_(reverse=True) * a._vector_(reverse=True)
292
((0, 0, 1, 0), (0, 1, 0, 0))
293
sage: vector(a), matrix(a) * vector(a)
294
((0, 1, 0, 0), (0, 0, 1, 0))
295
"""
296
import sage.matrix.matrix_space
297
298
K = self.parent()
299
a = K.gen()
300
x = K(1)
301
d = K.degree()
302
303
columns = []
304
305
if not reverse:
306
l = xrange(d)
307
else:
308
l = reversed(range(d))
309
310
for i in l:
311
columns.append( (self * x)._vector_() )
312
x *= a
313
314
k = K.base_ring()
315
M = sage.matrix.matrix_space.MatrixSpace(k, d)
316
317
if reverse:
318
return M(columns)
319
else:
320
return M(columns).transpose()
321
def _latex_(self):
322
r"""
323
Return the latex representation of self, which is just the
324
latex representation of the polynomial representation of self.
325
326
EXAMPLES::
327
328
sage: k.<b> = GF(5^2); k
329
Finite Field in b of size 5^2
330
sage: b._latex_()
331
'b'
332
sage: (b^2+1)._latex_()
333
'b + 4'
334
"""
335
if self.parent().degree()>1:
336
return self.polynomial()._latex_()
337
else:
338
return str(self)
339
340
def _pari_init_(self, var=None):
341
r"""
342
Return a string that defines this element when evaluated in PARI.
343
344
INPUT:
345
346
- ``var`` - default: ``None`` - a string for a new variable name to use.
347
348
EXAMPLES::
349
350
sage: S.<b> = GF(5^2); S
351
Finite Field in b of size 5^2
352
sage: b._pari_init_()
353
'Mod(Mod(1, 5)*b, Mod(1, 5)*b^2 + Mod(4, 5)*b + Mod(2, 5))'
354
sage: (2*b+3)._pari_init_()
355
'Mod(Mod(2, 5)*b + Mod(3, 5), Mod(1, 5)*b^2 + Mod(4, 5)*b + Mod(2, 5))'
356
357
TESTS:
358
359
The following tests against a bug fixed in trac ticket #11530. ::
360
361
sage: F.<d> = GF(3^4)
362
sage: F.modulus()
363
x^4 + 2*x^3 + 2
364
sage: d._pari_init_()
365
'Mod(Mod(1, 3)*d, Mod(1, 3)*d^4 + Mod(2, 3)*d^3 + Mod(2, 3))'
366
sage: (d^2+2*d+1)._pari_init_("p")
367
'Mod(Mod(1, 3)*p^2 + Mod(2, 3)*p + Mod(1, 3), Mod(1, 3)*p^4 + Mod(2, 3)*p^3 + Mod(2, 3))'
368
sage: d._pari_()
369
Mod(Mod(1, 3)*d, Mod(1, 3)*d^4 + Mod(2, 3)*d^3 + Mod(2, 3))
370
371
sage: K.<M> = GF(2^8)
372
sage: K.modulus()
373
x^8 + x^4 + x^3 + x^2 + 1
374
sage: (M^3+1)._pari_init_()
375
'Mod(Mod(1, 2)*M^3 + Mod(1, 2), Mod(1, 2)*M^8 + Mod(1, 2)*M^4 + Mod(1, 2)*M^3 + Mod(1, 2)*M^2 + Mod(1, 2))'
376
sage: M._pari_init_(var='foo')
377
'Mod(Mod(1, 2)*foo, Mod(1, 2)*foo^8 + Mod(1, 2)*foo^4 + Mod(1, 2)*foo^3 + Mod(1, 2)*foo^2 + Mod(1, 2))'
378
"""
379
if var is None:
380
var = self.parent().variable_name()
381
g = self.parent().modulus()._pari_with_name(var)
382
f = self.polynomial()._pari_with_name(var)
383
return 'Mod({0}, {1})'.format(f, g)
384
385
def charpoly(self, var='x', algorithm='matrix'):
386
"""
387
Return the characteristic polynomial of self as a polynomial with given variable.
388
389
INPUT:
390
391
- ``var`` - string (default: 'x')
392
393
- ``algorithm`` - string (default: 'matrix')
394
395
- 'matrix' - return the charpoly computed from the matrix of
396
left multiplication by self
397
398
- 'pari' -- use pari's charpoly routine on polymods, which
399
is not very good except in small cases
400
401
The result is not cached.
402
403
EXAMPLES::
404
405
sage: k.<a> = GF(19^2)
406
sage: parent(a)
407
Finite Field in a of size 19^2
408
sage: a.charpoly('X')
409
X^2 + 18*X + 2
410
sage: a^2 + 18*a + 2
411
0
412
sage: a.charpoly('X', algorithm='pari')
413
X^2 + 18*X + 2
414
"""
415
if algorithm == 'matrix':
416
return self._matrix_().charpoly(var)
417
elif algorithm == 'pari':
418
from sage.rings.polynomial.polynomial_ring_constructor import PolynomialRing
419
R = PolynomialRing(self.parent().prime_subfield(), var)
420
return R(self._pari_().charpoly('x').lift())
421
else:
422
raise ValueError, "unknown algorithm '%s'"%algorithm
423
424
def norm(self):
425
"""
426
Return the norm of self down to the prime subfield.
427
428
This is the product of the Galois conjugates of self.
429
430
EXAMPLES::
431
432
sage: S.<b> = GF(5^2); S
433
Finite Field in b of size 5^2
434
sage: b.norm()
435
2
436
sage: b.charpoly('t')
437
t^2 + 4*t + 2
438
439
Next we consider a cubic extension::
440
441
sage: S.<a> = GF(5^3); S
442
Finite Field in a of size 5^3
443
sage: a.norm()
444
2
445
sage: a.charpoly('t')
446
t^3 + 3*t + 3
447
sage: a * a^5 * (a^25)
448
2
449
"""
450
f = self.charpoly('x')
451
n = f[0]
452
if f.degree() % 2 != 0:
453
return -n
454
else:
455
return n
456
457
def trace(self):
458
"""
459
Return the trace of this element, which is the sum of the
460
Galois conjugates.
461
462
EXAMPLES::
463
464
sage: S.<a> = GF(5^3); S
465
Finite Field in a of size 5^3
466
sage: a.trace()
467
0
468
sage: a.charpoly('t')
469
t^3 + 3*t + 3
470
sage: a + a^5 + a^25
471
0
472
sage: z = a^2 + a + 1
473
sage: z.trace()
474
2
475
sage: z.charpoly('t')
476
t^3 + 3*t^2 + 2*t + 2
477
sage: z + z^5 + z^25
478
2
479
"""
480
return self.parent().prime_subfield()(self._pari_().trace().lift())
481
482
def multiplicative_order(self):
483
r"""
484
Return the multiplicative order of this field element.
485
486
EXAMPLE::
487
488
sage: S.<a> = GF(5^3); S
489
Finite Field in a of size 5^3
490
sage: a.multiplicative_order()
491
124
492
sage: (a^8).multiplicative_order()
493
31
494
sage: S(0).multiplicative_order()
495
Traceback (most recent call last):
496
...
497
ArithmeticError: Multiplicative order of 0 not defined.
498
"""
499
import sage.rings.arith
500
501
if self.is_zero():
502
raise ArithmeticError("Multiplicative order of 0 not defined.")
503
n = self._parent.order() - 1
504
F = self._parent.factored_unit_order()[0]
505
order = 1
506
for p, e in F:
507
# Determine the power of p that divides the order.
508
a = self**(n//(p**e))
509
while a != 1:
510
order = order * p
511
a = a**p
512
513
return order
514
515
def additive_order(self):
516
"""
517
Return the additive order of this finite field element.
518
519
EXAMPLES::
520
521
sage: k.<a> = FiniteField(2^12, 'a')
522
sage: b = a^3 + a + 1
523
sage: b.additive_order()
524
2
525
sage: k(0).additive_order()
526
1
527
"""
528
if self.is_zero():
529
from sage.rings.integer import Integer
530
return Integer(1)
531
return self.parent().characteristic()
532
533
def nth_root(self, n, extend = False, all = False, algorithm=None, cunningham=False):
534
r"""
535
Returns an `n`\th root of ``self``.
536
537
INPUT:
538
539
- ``n`` - integer `\geq 1`
540
541
- ``extend`` - bool (default: ``False``); if ``True``, return an `n`\th
542
root in an extension ring, if necessary. Otherwise, raise a
543
ValueError if the root is not in the base ring. Warning:
544
this option is not implemented!
545
546
- ``all`` - bool (default: ``False``); if ``True``, return all `n`\th
547
roots of ``self``, instead of just one.
548
549
- ``algorithm`` - string (default: ``None``); 'Johnston' is the only
550
currently supported option. For IntegerMod elements, the problem
551
is reduced to the prime modulus case using CRT and `p`-adic logs,
552
and then this algorithm used.
553
554
OUTPUT:
555
556
If self has an `n`\th root, returns one (if ``all`` is ``False``) or a
557
list of all of them (if ``all`` is ``True``).
558
Otherwise, raises a ``ValueError`` (if ``extend`` is ``False``)
559
or a ``NotImplementedError`` (if ``extend`` is ``True``).
560
561
.. warning::
562
563
The ``extend`` option is not implemented (yet).
564
565
EXAMPLES::
566
567
sage: K = GF(31)
568
sage: a = K(22)
569
sage: K(22).nth_root(7)
570
13
571
sage: K(25).nth_root(5)
572
5
573
sage: K(23).nth_root(3)
574
29
575
576
sage: K.<a> = GF(625)
577
sage: (3*a^2+a+1).nth_root(13)**13
578
3*a^2 + a + 1
579
580
sage: k.<a> = GF(29^2)
581
sage: b = a^2 + 5*a + 1
582
sage: b.nth_root(11)
583
3*a + 20
584
sage: b.nth_root(5)
585
Traceback (most recent call last):
586
...
587
ValueError: no nth root
588
sage: b.nth_root(5, all = True)
589
[]
590
sage: b.nth_root(3, all = True)
591
[14*a + 18, 10*a + 13, 5*a + 27]
592
593
sage: k.<a> = GF(29^5)
594
sage: b = a^2 + 5*a + 1
595
sage: b.nth_root(5)
596
19*a^4 + 2*a^3 + 2*a^2 + 16*a + 3
597
sage: b.nth_root(7)
598
Traceback (most recent call last):
599
...
600
ValueError: no nth root
601
sage: b.nth_root(4, all=True)
602
[]
603
604
TESTS::
605
606
sage: for p in [2,3,5,7,11]: # long time, random because of PARI warnings
607
....: for n in [2,5,10]:
608
....: q = p^n
609
....: K.<a> = GF(q)
610
....: for r in (q-1).divisors():
611
....: if r == 1: continue
612
....: x = K.random_element()
613
....: y = x^r
614
....: assert y.nth_root(r)^r == y
615
....: assert (y^41).nth_root(41*r)^(41*r) == y^41
616
....: assert (y^307).nth_root(307*r)^(307*r) == y^307
617
sage: k.<a> = GF(4)
618
sage: a.nth_root(0,all=True)
619
[]
620
sage: k(1).nth_root(0,all=True)
621
[a, a + 1, 1]
622
623
ALGORITHMS:
624
625
- The default is currently an algorithm described in the following paper:
626
627
Johnston, Anna M. A generalized qth root algorithm. Proceedings of the tenth annual ACM-SIAM symposium on Discrete algorithms. Baltimore, 1999: pp 929-930.
628
629
AUTHOR:
630
631
- David Roe (2010-02-13)
632
"""
633
if self.is_zero():
634
if n <= 0:
635
if all: return []
636
else: raise ValueError
637
if all: return [self]
638
else: return self
639
if n < 0:
640
self = ~self
641
n = -n
642
elif n == 0:
643
if self == 1:
644
if all: return [a for a in self.parent().list() if a != 0]
645
else: return self
646
else:
647
if all: return []
648
else: raise ValueError
649
if extend:
650
raise NotImplementedError
651
from sage.rings.integer import Integer
652
n = Integer(n)
653
return self._nth_root_common(n, all, algorithm, cunningham)
654
655
def pth_power(self, int k = 1):
656
"""
657
Return the `(p^k)^{th}` power of self, where `p` is the
658
characteristic of the field.
659
660
INPUT:
661
662
- ``k`` - integer (default: 1, must fit in C int type)
663
664
Note that if `k` is negative, then this computes the appropriate root.
665
666
EXAMPLES::
667
668
sage: F.<a> = GF(29^2)
669
sage: z = a^2 + 5*a + 1
670
sage: z.pth_power()
671
19*a + 20
672
sage: z.pth_power(10)
673
10*a + 28
674
sage: z.pth_power(-10) == z
675
True
676
sage: F.<b> = GF(2^12)
677
sage: y = b^3 + b + 1
678
sage: y == (y.pth_power(-3))^(2^3)
679
True
680
sage: y.pth_power(2)
681
b^7 + b^6 + b^5 + b^4 + b^3 + b
682
"""
683
p = self.additive_order()
684
n = self.parent().degree()
685
return self**(p**(k % n))
686
687
frobenius = pth_power
688
689
def pth_root(self, int k = 1):
690
"""
691
Return the `(p^k)^{th}` root of self, where `p` is the characteristic
692
of the field.
693
694
INPUT:
695
696
- ``k`` - integer (default: 1, must fit in C int type)
697
698
Note that if `k` is negative, then this computes the appropriate power.
699
700
EXAMPLES::
701
702
sage: F.<b> = GF(2^12)
703
sage: y = b^3 + b + 1
704
sage: y == (y.pth_root(3))^(2^3)
705
True
706
sage: y.pth_root(2)
707
b^11 + b^10 + b^9 + b^7 + b^5 + b^4 + b^2 + b
708
"""
709
return self.pth_power(-k)
710
711
712