Finite Fields of characteristic 2.

This implementation uses NTL’s GF2E class to perform the arithmetic and is the standard implementation for GF(2^n) for n >= 16.

AUTHORS:

class sage.rings.finite_rings.element_ntl_gf2e.Cache_ntl_gf2e

Bases: sage.rings.finite_rings.element_base.Cache_base

This class stores information for an NTL finite field in a Cython class so that elements can access it quickly.

It’s modeled on NativeIntStruct, but includes many functions that were previously included in the parent (see trac ticket #12062).

degree()

If the field has cardinality \(2^n\) this method returns \(n\).

EXAMPLES:

sage: k.<a> = GF(2^64)
sage: k._cache.degree()
64
fetch_int(number)

Given an integer less than \(p^n\) with base \(2\) representation \(a_0 + a_1 \cdot 2 + \cdots + a_k 2^k\), this returns \(a_0 + a_1 x + \cdots + a_k x^k\), where \(x\) is the generator of this finite field.

INPUT:

  • number – an integer, of size less than the cardinality

EXAMPLES:

sage: k.<a> = GF(2^48)
sage: k._cache.fetch_int(2^33 + 2 + 1)
a^33 + a + 1
import_data(e)

EXAMPLES:

sage: k.<a> = GF(2^17)
sage: V = k.vector_space(map=False)
sage: v = [1,0,0,0,0,1,0,0,1,0,0,0,0,1,0,0,0]
sage: k._cache.import_data(v)
a^13 + a^8 + a^5 + 1
sage: k._cache.import_data(V(v))
a^13 + a^8 + a^5 + 1
order()

Return the cardinality of the field.

EXAMPLES:

sage: k.<a> = GF(2^64)
sage: k._cache.order()
18446744073709551616
polynomial()

Returns the list of 0’s and 1’s giving the defining polynomial of the field.

EXAMPLES:

sage: k.<a> = GF(2^20,modulus="minimal_weight")
sage: k._cache.polynomial()
[1, 0, 0, 1, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 1]
class sage.rings.finite_rings.element_ntl_gf2e.FiniteField_ntl_gf2eElement

Bases: sage.rings.finite_rings.element_base.FinitePolyExtElement

An element of an NTL:GF2E finite field.

charpoly(var='x')

Return the characteristic polynomial of self as a polynomial in var over the prime subfield.

INPUT:

  • var – string (default: 'x')

OUTPUT:

polynomial

EXAMPLES:

sage: k.<a> = GF(2^8, impl="ntl")
sage: b = a^3 + a
sage: b.minpoly()
x^4 + x^3 + x^2 + x + 1
sage: b.charpoly()
x^8 + x^6 + x^4 + x^2 + 1
sage: b.charpoly().factor()
(x^4 + x^3 + x^2 + x + 1)^2
sage: b.charpoly('Z')
Z^8 + Z^6 + Z^4 + Z^2 + 1
integer_representation()

Return the int representation of self. When self is in the prime subfield, the integer returned is equal to self and not to log_repr.

Elements of this field are represented as ints in as follows: for \(e \in \GF{p}[x]\) with \(e = a_0 + a_1 x + a_2 x^2 + \cdots\), \(e\) is represented as: \(n = a_0 + a_1 p + a_2 p^2 + \cdots\).

EXAMPLES:

sage: k.<a> = GF(2^20)
sage: a.integer_representation()
2
sage: (a^2 + 1).integer_representation()
5
sage: k.<a> = GF(2^70)
sage: (a^65 + a^64 + 1).integer_representation()
55340232221128654849L
is_one()

Return True if self == k(1).

Equivalent to self != k(0).

EXAMPLES:

sage: k.<a> = GF(2^20)
sage: a.is_one() # indirect doctest
False
sage: k(1).is_one()
True
is_square()

Return True as every element in \(\GF{2^n}\) is a square.

EXAMPLES:

sage: k.<a> = GF(2^18)
sage: e = k.random_element()
sage: e
a^15 + a^14 + a^13 + a^11 + a^10 + a^9 + a^6 + a^5 + a^4 + 1
sage: e.is_square()
True
sage: e.sqrt()
a^16 + a^15 + a^14 + a^11 + a^9 + a^8 + a^7 + a^6 + a^4 + a^3 + 1
sage: e.sqrt()^2 == e
True
is_unit()

Return True if self is nonzero, so it is a unit as an element of the finite field.

EXAMPLES:

sage: k.<a> = GF(2^17)
sage: a.is_unit()
True
sage: k(0).is_unit()
False
log(base)

Return \(x\) such that \(b^x = a\), where \(x\) is \(a\) and \(b\) is the base.

INPUT:

  • base – finite field element that generates the multiplicative group.

OUTPUT:

Integer \(x\) such that \(a^x = b\), if it exists. Raises a ValueError exception if no such \(x\) exists.

EXAMPLES:

sage: F = GF(17)
sage: F(3^11).log(F(3))
11
sage: F = GF(113)
sage: F(3^19).log(F(3))
19
sage: F = GF(next_prime(10000))
sage: F(23^997).log(F(23))
997

sage: F = FiniteField(2^10, 'a')
sage: g = F.gen()
sage: b = g; a = g^37
sage: a.log(b)
37
sage: b^37; a
a^8 + a^7 + a^4 + a + 1
a^8 + a^7 + a^4 + a + 1

AUTHOR: David Joyner and William Stein (2005-11)

minpoly(var='x')

Return the minimal polynomial of self, which is the smallest degree polynomial \(f \in \GF{2}[x]\) such that f(self) == 0.

INPUT:

  • var – string (default: 'x')

OUTPUT:

polynomial

EXAMPLES:

sage: K.<a> = GF(2^100)
sage: f = a.minpoly(); f
x^100 + x^57 + x^56 + x^55 + x^52 + x^48 + x^47 + x^46 + x^45 + x^44 + x^43 + x^41 + x^37 + x^36 + x^35 + x^34 + x^31 + x^30 + x^27 + x^25 + x^24 + x^22 + x^20 + x^19 + x^16 + x^15 + x^11 + x^9 + x^8 + x^6 + x^5 + x^3 + 1
sage: f(a)
0
sage: g = K.random_element()
sage: g.minpoly()(g)
0
polynomial(name=None)

Return self viewed as a polynomial over self.parent().prime_subfield().

INPUT:

  • name – (optional) variable name

EXAMPLES:

sage: k.<a> = GF(2^17)
sage: e = a^15 + a^13 + a^11 + a^10 + a^9 + a^8 + a^7 + a^6 + a^4 + a + 1
sage: e.polynomial()
a^15 + a^13 + a^11 + a^10 + a^9 + a^8 + a^7 + a^6 + a^4 + a + 1

sage: from sage.rings.polynomial.polynomial_element import is_Polynomial
sage: is_Polynomial(e.polynomial())
True

sage: e.polynomial('x')
x^15 + x^13 + x^11 + x^10 + x^9 + x^8 + x^7 + x^6 + x^4 + x + 1
sqrt(all=False, extend=False)

Return a square root of this finite field element in its parent.

EXAMPLES:

sage: k.<a> = GF(2^20)
sage: a.is_square()
True
sage: a.sqrt()
a^19 + a^15 + a^14 + a^12 + a^9 + a^7 + a^4 + a^3 + a + 1
sage: a.sqrt()^2 == a
True

This failed before trac ticket #4899:

sage: GF(2^16,'a')(1).sqrt()
1
trace()

Return the trace of self.

EXAMPLES:

sage: K.<a> = GF(2^25)
sage: a.trace()
0
sage: a.charpoly()
x^25 + x^8 + x^6 + x^2 + 1
sage: parent(a.trace())
Finite Field of size 2

sage: b = a+1
sage: b.trace()
1
sage: b.charpoly()[1]
1
weight()

Returns the number of non-zero coefficients in the polynomial representation of self.

EXAMPLES:

sage: K.<a> = GF(2^21)
sage: a.weight()
1
sage: (a^5+a^2+1).weight()
3
sage: b = 1/(a+1); b
a^20 + a^19 + a^18 + a^17 + a^16 + a^15 + a^14 + a^13 + a^12 + a^11 + a^10 + a^9 + a^8 + a^7 + a^6 + a^4 + a^3 + a^2
sage: b.weight()
18
sage.rings.finite_rings.element_ntl_gf2e.unpickleFiniteField_ntl_gf2eElement(parent, elem)

EXAMPLES:

sage: k.<a> = GF(2^20)
sage: e = k.random_element()
sage: f = loads(dumps(e)) # indirect doctest
sage: e == f
True