Path: blob/master/src/sage/rings/finite_rings/finite_field_ext_pari.py
8820 views
"""1Finite Extension Fields implemented via PARI.23AUTHORS:45- William Stein: initial version67- Jeroen Demeyer (2010-12-16): fix formatting of docstrings (:trac:`10487`)8"""9#*****************************************************************************10# Copyright (C) 2005,2007 William Stein <[email protected]>11# Copyright (C) 2010 Jeroen Demeyer <[email protected]>12#13# Distributed under the terms of the GNU General Public License (GPL)14# as published by the Free Software Foundation; either version 2 of15# the License, or (at your option) any later version.16# http://www.gnu.org/licenses/17#*****************************************************************************181920import sage.rings.polynomial.polynomial_element as polynomial_element21import sage.rings.polynomial.multi_polynomial_element as multi_polynomial_element2223import sage.rings.integer as integer24import sage.rings.rational as rational2526import sage.libs.pari.all as pari2728import element_ext_pari2930from sage.rings.finite_rings.finite_field_base import FiniteField as FiniteField_generic313233import sage.interfaces.gap3435class FiniteField_ext_pari(FiniteField_generic):36r"""37Finite Field of order `q`, where `q` is a prime power (not a prime),38implemented using PARI ``POLMOD``. This implementation is the default39implementation for `q \geq 2^{16}`.4041INPUT:4243- ``q`` -- integer, size of the finite field, not prime4445- ``name`` -- variable used for printing element of the finite46field. Also, two finite fields are considered equal if they47have the same variable name, and not otherwise.4849- ``modulus`` -- you may provide a polynomial to use for reduction or50a string:5152- ``'conway'`` -- force the use of a Conway polynomial, will53raise a ``RuntimeError`` if none is found in the database54- ``'random'`` -- use a random irreducible polynomial55- ``'default'`` -- a Conway polynomial is used if found. Otherwise56a random polynomial is used5758OUTPUT:5960A finite field of order `q` with the given variable name6162EXAMPLES::6364sage: from sage.rings.finite_rings.finite_field_ext_pari import FiniteField_ext_pari65sage: k = FiniteField_ext_pari(9, 'a')66sage: k67Finite Field in a of size 3^268sage: k.is_field()69True70sage: k.characteristic()71372sage: a = k.gen()73sage: a74a75sage: a.parent()76Finite Field in a of size 3^277sage: a.charpoly('x')78x^2 + 2*x + 279sage: [a^i for i in range(8)]80[1, a, a + 1, 2*a + 1, 2, 2*a, 2*a + 2, a + 2]8182Fields can be coerced into sets or list and iterated over::8384sage: list(k)85[0, 1, 2, a, a + 1, a + 2, 2*a, 2*a + 1, 2*a + 2]8687The following is a native Python set::8889sage: set(k)90set([0, 1, 2, a, a + 1, a + 2, 2*a, 2*a + 1, 2*a + 2])9192And the following is a Sage set::9394sage: Set(k)95{0, 1, 2, a, a + 1, a + 2, 2*a, 2*a + 1, 2*a + 2}9697We can also make a list via comprehension:98sage: [x for x in k]99[0, 1, 2, a, a + 1, a + 2, 2*a, 2*a + 1, 2*a + 2]100101Next we compute with the finite field of order 16, where102the name is named ``b``::103104sage: from sage.rings.finite_rings.finite_field_ext_pari import FiniteField_ext_pari105sage: k16 = FiniteField_ext_pari(16, "b")106sage: z = k16.gen()107sage: z108b109sage: z.charpoly('x')110x^4 + x + 1111sage: k16.is_field()112True113sage: k16.characteristic()1142115sage: z.multiplicative_order()11615117118Of course one can also make prime finite fields::119120sage: k = FiniteField(7)121122Note that the generator is 1::123124sage: k.gen()1251126sage: k.gen().multiplicative_order()1271128129Prime finite fields are implemented elsewhere, they cannot be130constructed using :class:`FiniteField_ext_pari`::131132sage: k = FiniteField_ext_pari(7, 'a')133Traceback (most recent call last):134...135ValueError: The size of the finite field must not be prime.136137Illustration of dumping and loading::138139sage: K = FiniteField(7)140sage: loads(K.dumps()) == K141True142sage: K = FiniteField(7^10, 'b', impl='pari_mod')143sage: loads(K.dumps()) == K144True145sage: K = FiniteField(7^10, 'a', impl='pari_mod')146sage: loads(K.dumps()) == K147True148149In this example `K` is large enough that Conway polynomials are not150used. Note that when the field is dumped the defining polynomial `f`151is also dumped. Since `f` is determined by a random algorithm, it's152important that `f` is dumped as part of `K`. If you quit Sage and153restart and remake a finite field of the same order (and the order is154large enough so that there is no Conway polynomial), then defining155polynomial is probably different. However, if you load a previously156saved field, that will have the same defining polynomial. ::157158sage: K = GF(10007^10, 'a', impl='pari_mod')159sage: loads(K.dumps()) == K160True161162.. NOTE::163164We do NOT yet define natural consistent inclusion maps165between different finite fields.166"""167def __init__(self, q, name, modulus=None):168"""169Create finite field of order `q` with variable printed as name.170171EXAMPLES::172173sage: from sage.rings.finite_rings.finite_field_ext_pari import FiniteField_ext_pari174sage: k = FiniteField_ext_pari(9, 'a'); k175Finite Field in a of size 3^2176"""177if element_ext_pari.dynamic_FiniteField_ext_pariElement is None: element_ext_pari._late_import()178from constructor import FiniteField as GF179q = integer.Integer(q)180if q < 2:181raise ArithmeticError, "q must be a prime power"182from sage.structure.proof.all import arithmetic183proof = arithmetic()184if proof:185F = q.factor()186else:187from sage.rings.arith import is_pseudoprime_small_power188F = is_pseudoprime_small_power(q, get_data=True)189if len(F) != 1:190raise ArithmeticError, "q must be a prime power"191192if F[0][1] > 1:193base_ring = GF(F[0][0])194else:195raise ValueError, "The size of the finite field must not be prime."196#base_ring = self197198FiniteField_generic.__init__(self, base_ring, name, normalize=True)199200self._kwargs = {}201self.__char = F[0][0]202self.__pari_one = pari.pari(1).Mod(self.__char)203self.__degree = integer.Integer(F[0][1])204self.__order = q205self.__is_field = True206207if modulus is None or modulus == "default":208from conway_polynomials import exists_conway_polynomial209if exists_conway_polynomial(self.__char, self.__degree):210modulus = "conway"211else:212modulus = "random"213214if isinstance(modulus,str):215if modulus == "conway":216from conway_polynomials import conway_polynomial217modulus = conway_polynomial(self.__char, self.__degree)218elif modulus == "random":219# The following is fast/deterministic, but has serious problems since220# it crashes on 64-bit machines, and I can't figure out why:221# self.__pari_modulus = pari.pari.finitefield_init(self.__char, self.__degree, self.variable_name())222# So instead we iterate through random polys until we find an irreducible one.223224R = GF(self.__char)['x']225while True:226modulus = R.random_element(self.__degree)227modulus = modulus.monic()228if modulus.degree() == self.__degree and modulus.is_irreducible():229break230else:231raise ValueError("Modulus parameter not understood")232233elif isinstance(modulus, (list, tuple)):234modulus = GF(self.__char)['x'](modulus)235elif sage.rings.polynomial.polynomial_element.is_Polynomial(modulus):236if modulus.parent() is not base_ring:237modulus = modulus.change_ring(base_ring)238else:239raise ValueError("Modulus parameter not understood")240241self.__modulus = modulus242f = pari.pari(str(modulus))243self.__pari_modulus = f.subst(modulus.parent().variable_name(), 'a') * self.__pari_one244self.__gen = element_ext_pari.FiniteField_ext_pariElement(self, pari.pari('a'))245246self._zero_element = self._element_constructor_(0)247self._one_element = self._element_constructor_(1)248249def __reduce__(self):250"""251For pickling.252253EXAMPLES::254255sage: k.<b> = GF(5^20, impl='pari_mod'); type(k)256<class 'sage.rings.finite_rings.finite_field_ext_pari.FiniteField_ext_pari_with_category'>257sage: k is loads(dumps(k))258True259"""260return self._factory_data[0].reduce_data(self)261262def __cmp__(self, other):263"""264Compare ``self`` to ``other``.265266EXAMPLES::267268sage: k = GF(7^20, 'a', impl='pari_mod')269sage: k == loads(dumps(k))270True271"""272if not isinstance(other, FiniteField_ext_pari):273return cmp(type(self), type(other))274return cmp((self.__order, self.variable_name(), self.__modulus), (other.__order, other.variable_name(), other.__modulus))275276def __richcmp__(left, right, op):277r"""278Compare ``left`` with ``right``.279280EXAMPLE::281282sage: k.<a> = GF(2^17, impl='pari_mod')283sage: j.<b> = GF(2^18, impl='pari_mod')284sage: k == j285False286287sage: GF(2^17, 'a', impl='pari_mod') == copy(GF(2^17, 'a', impl='pari_mod'))288True289"""290return left._richcmp_helper(right, op)291292def _pari_one(self):293r"""294The PARI object ``Mod(1,p)``. This is implementation specific295and should be ignored by users.296297EXAMPLES::298299sage: k = GF(7^20, 'a', impl='pari_mod')300sage: k._pari_one()301Mod(1, 7)302"""303return self.__pari_one304305def _pari_modulus(self):306"""307The polynomial mod `p` that defines the finite field, as a PARI308object. This is implementation specific, and some finite fields309might not be implemented using PARI, so you should avoid using310this function.311312OUTPUT:313314- ``gen`` -- a PARI polynomial gen315316EXAMPLES::317318sage: from sage.rings.finite_rings.finite_field_ext_pari import FiniteField_ext_pari319sage: FiniteField_ext_pari(19^2, 'a')._pari_modulus()320Mod(1, 19)*a^2 + Mod(18, 19)*a + Mod(2, 19)321322sage: FiniteField_ext_pari(13^3, 'a')._pari_modulus()323Mod(1, 13)*a^3 + Mod(2, 13)*a + Mod(11, 13)324325Note that the PARI modulus is always in terms of a, even if326the field variable isn't. This is because the specific choice327of variable name has meaning in PARI, i.e., it can't be328arbitrary. ::329330sage: from sage.rings.finite_rings.finite_field_ext_pari import FiniteField_ext_pari331sage: FiniteField_ext_pari(2^4, "b")._pari_modulus()332Mod(1, 2)*a^4 + Mod(1, 2)*a + Mod(1, 2)333"""334return self.__pari_modulus335336def gen(self, n=0):337"""338Return a generator of the finite field.339340This generator is a root of the defining polynomial of the finite341field, and might differ between different runs or different342architectures.343344.. WARNING::345346This generator is not guaranteed to be a generator347for the multiplicative group. To obtain the latter, use348:meth:`~sage.rings.finite_rings.finite_field_base.FiniteFields.multiplicative_generator()`.349350INPUT:351352- ``n`` -- ignored353354OUTPUT:355356Field generator of finite field357358EXAMPLES::359360sage: from sage.rings.finite_rings.finite_field_ext_pari import FiniteField_ext_pari361sage: FiniteField_ext_pari(2^4, "b").gen()362b363sage: k = FiniteField_ext_pari(3^4, "alpha")364sage: a = k.gen()365sage: a366alpha367sage: a^4368alpha^3 + 1369370"""371return self.__gen372373def characteristic(self):374"""375Returns the characteristic of the finite field, which is a376prime number.377378EXAMPLES::379380sage: from sage.rings.finite_rings.finite_field_ext_pari import FiniteField_ext_pari381sage: k = FiniteField_ext_pari(3^4, 'a')382sage: k.characteristic()3833384"""385return self.__char386387def modulus(self):388r"""389Return the minimal polynomial of the generator of ``self`` in390``self.polynomial_ring('x')``.391392EXAMPLES::393394sage: F.<a> = GF(7^20, 'a', impl='pari_mod')395sage: f = F.modulus(); f396x^20 + x^12 + 6*x^11 + 2*x^10 + 5*x^9 + 2*x^8 + 3*x^7 + x^6 + 3*x^5 + 3*x^3 + x + 3397398sage: f(a)3990400"""401return self.__modulus402403def degree(self):404"""405Returns the degree of the finite field, which is a positive406integer.407408EXAMPLES::409410sage: from sage.rings.finite_rings.finite_field_ext_pari import FiniteField_ext_pari411sage: FiniteField_ext_pari(3^20, 'a').degree()41220413"""414return self.__degree415416def _element_constructor_(self, x):417r"""418Coerce ``x`` into the finite field.419420INPUT:421422- ``x`` -- object423424OUTPUT:425426If possible, makes a finite field element from ``x``.427428EXAMPLES::429430sage: from sage.rings.finite_rings.finite_field_ext_pari import FiniteField_ext_pari431sage: k = FiniteField_ext_pari(3^4, 'a')432sage: b = k(5) # indirect doctest433sage: b.parent()434Finite Field in a of size 3^4435sage: a = k.gen()436sage: k(a + 2)437a + 2438439Univariate polynomials coerce into finite fields by evaluating440the polynomial at the field's generator::441442sage: from sage.rings.finite_rings.finite_field_ext_pari import FiniteField_ext_pari443sage: R.<x> = QQ[]444sage: k, a = FiniteField_ext_pari(5^2, 'a').objgen()445sage: k(R(2/3))4464447sage: k(x^2)448a + 3449sage: R.<x> = GF(5)[]450sage: k(x^3-2*x+1)4512*a + 4452453sage: x = polygen(QQ)454sage: k(x^25)455a456457sage: Q, q = FiniteField_ext_pari(5^7, 'q').objgen()458sage: L = GF(5)459sage: LL.<xx> = L[]460sage: Q(xx^2 + 2*xx + 4)461q^2 + 2*q + 4462463Multivariate polynomials only coerce if constant::464465sage: R = k['x,y,z']; R466Multivariate Polynomial Ring in x, y, z over Finite Field in a of size 5^2467sage: k(R(2))4682469sage: R = QQ['x,y,z']470sage: k(R(1/5))471Traceback (most recent call last):472...473TypeError: unable to coerce474475Gap elements can also be coerced into finite fields::476477sage: from sage.rings.finite_rings.finite_field_ext_pari import FiniteField_ext_pari478sage: F = FiniteField_ext_pari(8, 'a')479sage: a = F.multiplicative_generator(); a480a481sage: b = gap(a^3); b482Z(2^3)^3483sage: F(b)484a + 1485sage: a^3486a + 1487488sage: a = GF(13)(gap('0*Z(13)')); a4890490sage: a.parent()491Finite Field of size 13492493sage: F = GF(16, 'a', impl='pari_mod')494sage: F(gap('Z(16)^3'))495a^3496sage: F(gap('Z(16)^2'))497a^2498499You can also call a finite extension field with a string500to produce an element of that field, like this::501502sage: k = GF(2^8, 'a', impl='pari_mod')503sage: k('a^200')504a^4 + a^3 + a^2505506This is especially useful for fast conversions from Singular etc.507to ``FiniteField_ext_pariElements``.508509AUTHORS:510511- David Joyner (2005-11)512- Martin Albrecht (2006-01-23)513- Martin Albrecht (2006-03-06): added coercion from string514"""515if isinstance(x, element_ext_pari.FiniteField_ext_pariElement):516if x.parent() is self:517return x518elif x.parent() == self:519# canonically isomorphic finite fields520return element_ext_pari.FiniteField_ext_pariElement(self, x)521else:522# This is where we *would* do coercion from one finite field to another...523raise TypeError, "no coercion defined"524525elif sage.interfaces.gap.is_GapElement(x):526from sage.interfaces.gap import gfq_gap_to_sage527try:528return gfq_gap_to_sage(x, self)529except (ValueError, IndexError, TypeError):530raise TypeError, "no coercion defined"531532if isinstance(x, (int, long, integer.Integer, rational.Rational,533pari.pari_gen, list)):534535return element_ext_pari.FiniteField_ext_pariElement(self, x)536537elif isinstance(x, multi_polynomial_element.MPolynomial):538if x.is_constant():539return self(x.constant_coefficient())540else:541raise TypeError, "no coercion defined"542543elif isinstance(x, polynomial_element.Polynomial):544if x.is_constant():545return self(x.constant_coefficient())546else:547return x.change_ring(self)(self.gen())548549elif isinstance(x, str):550x = x.replace(self.variable_name(),'a')551x = pari.pari(x)552t = x.type()553if t == 't_POL':554if (x.variable() == 'a' \555and x.polcoeff(0).type()[2] == 'I'): #t_INT and t_INTMOD556return self(x)557if t[2] == 'I': #t_INT and t_INTMOD558return self(x)559raise TypeError, "string element does not match this finite field"560561try:562if x.parent() == self.vector_space():563x = pari.pari('+'.join(['%s*a^%s'%(x[i], i) for i in range(self.degree())]))564return element_ext_pari.FiniteField_ext_pariElement(self, x)565except AttributeError:566pass567try:568return element_ext_pari.FiniteField_ext_pariElement(self, integer.Integer(x))569except TypeError, msg:570raise TypeError, "%s\nno coercion defined"%msg571572def __len__(self):573"""574The number of elements of the finite field.575576EXAMPLES::577578sage: from sage.rings.finite_rings.finite_field_ext_pari import FiniteField_ext_pari579sage: k = FiniteField_ext_pari(2^10, 'a')580sage: k581Finite Field in a of size 2^10582sage: len(k)5831024584"""585return self.__order586587def order(self):588"""589The number of elements of the finite field.590591EXAMPLES::592593sage: from sage.rings.finite_rings.finite_field_ext_pari import FiniteField_ext_pari594sage: k = FiniteField_ext_pari(2^10,'a')595sage: k596Finite Field in a of size 2^10597sage: k.order()5981024599"""600return self.__order601602def polynomial(self, name=None):603"""604Return the irreducible characteristic polynomial of the605generator of this finite field, i.e., the polynomial `f(x)` so606elements of the finite field as elements modulo `f`.607608EXAMPLES::609610sage: k = FiniteField(17)611sage: k.polynomial('x')612x613sage: from sage.rings.finite_rings.finite_field_ext_pari import FiniteField_ext_pari614sage: k = FiniteField_ext_pari(9,'a')615sage: k.polynomial('x')616x^2 + 2*x + 2617"""618if name is None:619name = self.variable_name()620try:621return self.__polynomial[name]622except (AttributeError, KeyError):623from constructor import FiniteField as GF624R = GF(self.characteristic())[name]625f = R(self._pari_modulus())626try:627self.__polynomial[name] = f628except (KeyError, AttributeError):629self.__polynomial = {}630self.__polynomial[name] = f631return f632633def __hash__(self):634"""635Return the hash of this field.636637EXAMPLES::638639sage: {GF(9, 'a', impl='pari_mod'): 1} # indirect doctest640{Finite Field in a of size 3^2: 1}641sage: {GF(9, 'b', impl='pari_mod'): 1} # indirect doctest642{Finite Field in b of size 3^2: 1}643"""644try:645return self.__hash646except AttributeError:647self.__hash = hash((self.__order, self.variable_name(), self.__modulus))648return self.__hash649650651