Hot-keys on this page
r m x p toggle line displays
j k next/prev highlighted chunk
0 (zero) top of page
1 (one) first highlighted chunk
""" Univariate Polynomials over GF(2) via NTL's GF2X.
AUTHOR: - Martin Albrecht (2008-10) initial implementation """
# We need to define this stuff before including the templating stuff # to make sure the function get_cparent is found since it is used in # 'polynomial_template.pxi'.
cdef inline cparent get_cparent(parent):
# first we include the definitions include "sage/libs/ntl/decl.pxi" include "sage/libs/ntl/ntl_GF2X_linkage.pxi"
# and then the interface include "polynomial_template.pxi"
from sage.libs.m4ri cimport mzd_write_bit, mzd_read_bit from sage.matrix.matrix_mod2_dense cimport Matrix_mod2_dense
cdef class Polynomial_GF2X(Polynomial_template): """ Univariate Polynomials over GF(2) via NTL's GF2X.
EXAMPLES::
sage: P.<x> = GF(2)[] sage: x^3 + x^2 + 1 x^3 + x^2 + 1 """ def __init__(self, parent, x=None, check=True, is_gen=False, construct=False): """ Create a new univariate polynomials over GF(2).
EXAMPLES::
sage: P.<x> = GF(2)[] sage: x^3 + x^2 + 1 x^3 + x^2 + 1
We check that the bug noted at :trac:`12724` is fixed::
sage: R.<x> = Zmod(2)[] sage: R([2^80]) 0 """ pass
cdef get_unsafe(self, Py_ssize_t i): """ Return the `i`-th coefficient of ``self``.
EXAMPLES::
sage: P.<x> = GF(2)[] sage: f = x^3 + x^2 + 1; f x^3 + x^2 + 1 sage: f[0] 1 sage: f[1] 0 sage: f[:50] == f True sage: f[:3] x^2 + 1 """
def __pari__(self, variable=None): """ EXAMPLES::
sage: P.<x> = GF(2)[] sage: f = x^3 + x^2 + 1 sage: pari(f) Mod(1, 2)*x^3 + Mod(1, 2)*x^2 + Mod(1, 2) """ #TODO: put this in a superclass
def modular_composition(Polynomial_GF2X self, Polynomial_GF2X g, Polynomial_GF2X h, algorithm=None): """ Compute `f(g) \pmod h`.
Both implementations use Brent-Kung's Algorithm 2.1 (*Fast Algorithms for Manipulation of Formal Power Series*, JACM 1978).
INPUT:
- ``g`` -- a polynomial - ``h`` -- a polynomial - ``algorithm`` -- either 'native' or 'ntl' (default: 'native')
EXAMPLES::
sage: P.<x> = GF(2)[] sage: r = 279 sage: f = x^r + x +1 sage: g = x^r sage: g.modular_composition(g, f) == g(g) % f True
sage: P.<x> = GF(2)[] sage: f = x^29 + x^24 + x^22 + x^21 + x^20 + x^16 + x^15 + x^14 + x^10 + x^9 + x^8 + x^7 + x^6 + x^5 + x^2 sage: g = x^31 + x^30 + x^28 + x^26 + x^24 + x^21 + x^19 + x^18 + x^11 + x^10 + x^9 + x^8 + x^5 + x^2 + 1 sage: h = x^30 + x^28 + x^26 + x^25 + x^24 + x^22 + x^21 + x^18 + x^17 + x^15 + x^13 + x^12 + x^11 + x^10 + x^9 + x^4 sage: f.modular_composition(g,h) == f(g) % h True
AUTHORS:
- Paul Zimmermann (2008-10) initial implementation - Martin Albrecht (2008-10) performance improvements """ raise TypeError("Parents of the first three parameters must match.")
cdef Polynomial_GF2X res cdef GF2XModulus_c modulus
t = cputime() sig_on() GF2X_CompMod(res.x, self.x, g.x, modulus) sig_off() verbose("NTL %5.3f s"%cputime(t),level=1) return res
cdef Py_ssize_t i, j, k, l, n, maxlength cdef Matrix_mod2_dense F, G, H
cdef GF2X_c gpow, g2, tt
# we store all matrices transposed for performance reasons
# first compute g^j mod h, 2 <= j < k # first deal with j=0 # precompute g^2 # we now process 2j, 4j, 8j, ... by squaring each time # we need that gpow = g^k at the end else: # k is even, last j is k-1
# split f in chunks of degree < k else:
# H is a n x l matrix now H[i,j] = sum(G[i,m]*F[m,j], # m=0..k-1) = sum(g^m[i] * f[j*k+m], m=0..k-1) where g^m[i] is # the coefficient of degree i in g^m and f[j*k+m] is the # coefficient of degree j*k+m in f thus f[j*k+m]*g^m[i] should # be multiplied by g^(j*k) gpow = (g^k) % h
#res = (res * gpow) % h
# res = res + parent([H[j,i] for i in range(0,n)])
r""" Return whether this polynomial is irreducible over `\GF{2}`.`
EXAMPLES::
sage: R.<x> = GF(2)[] sage: (x^2 + 1).is_irreducible() False sage: (x^3 + x + 1).is_irreducible() True
Test that caching works::
sage: R.<x> = GF(2)[] sage: f = x^2 + 1 sage: f.is_irreducible() False sage: f.is_irreducible.cache False
"""
# The three functions below are used in polynomial_ring.py, but are in # this Cython file since they call C++ functions. They return # polynomials as lists so that no variable has to be specified. # AUTHOR: Peter Bruin (June 2013)
""" Return the list of coefficients of the lexicographically smallest irreducible polynomial of degree `n` over the field of 2 elements.
EXAMPLES::
sage: from sage.rings.polynomial.polynomial_gf2x import GF2X_BuildIrred_list sage: GF2X_BuildIrred_list(2) [1, 1, 1] sage: GF2X_BuildIrred_list(3) [1, 1, 0, 1] sage: GF2X_BuildIrred_list(4) [1, 1, 0, 0, 1] sage: GF(2)['x'](GF2X_BuildIrred_list(33)) x^33 + x^6 + x^3 + x + 1 """ cdef GF2X_c f
""" Return the list of coefficients of an irreducible polynomial of degree `n` of minimal weight over the field of 2 elements.
EXAMPLES::
sage: from sage.rings.polynomial.polynomial_gf2x import GF2X_BuildIrred_list, GF2X_BuildSparseIrred_list sage: all([GF2X_BuildSparseIrred_list(n) == GF2X_BuildIrred_list(n) ....: for n in range(1,33)]) True sage: GF(2)['x'](GF2X_BuildSparseIrred_list(33)) x^33 + x^10 + 1 """ cdef GF2X_c f
""" Return the list of coefficients of an irreducible polynomial of degree `n` of minimal weight over the field of 2 elements.
EXAMPLES::
sage: from sage.rings.polynomial.polynomial_gf2x import GF2X_BuildRandomIrred_list sage: GF2X_BuildRandomIrred_list(2) [1, 1, 1] sage: GF2X_BuildRandomIrred_list(3) in [[1, 1, 0, 1], [1, 0, 1, 1]] True """ cdef GF2X_c tmp, f |