Dense univariate polynomials over \(\ZZ/n\ZZ\), implemented using FLINT¶
This module gives a fast implementation of \((\ZZ/n\ZZ)[x]\) whenever \(n\) is at
most sys.maxsize
. We use it by default in preference to NTL when the modulus
is small, falling back to NTL if the modulus is too large, as in the example
below.
EXAMPLES:
sage: R.<a> = PolynomialRing(Integers(100))
sage: type(a)
<class 'sage.rings.polynomial.polynomial_zmod_flint.Polynomial_zmod_flint'>
sage: R.<a> = PolynomialRing(Integers(5*2^64))
sage: type(a)
<class 'sage.rings.polynomial.polynomial_modn_dense_ntl.Polynomial_dense_modn_ntl_ZZ'>
sage: R.<a> = PolynomialRing(Integers(5*2^64), implementation="FLINT")
Traceback (most recent call last):
...
ValueError: FLINT does not support modulus 92233720368547758080
>>> from sage.all import *
>>> R = PolynomialRing(Integers(Integer(100)), names=('a',)); (a,) = R._first_ngens(1)
>>> type(a)
<class 'sage.rings.polynomial.polynomial_zmod_flint.Polynomial_zmod_flint'>
>>> R = PolynomialRing(Integers(Integer(5)*Integer(2)**Integer(64)), names=('a',)); (a,) = R._first_ngens(1)
>>> type(a)
<class 'sage.rings.polynomial.polynomial_modn_dense_ntl.Polynomial_dense_modn_ntl_ZZ'>
>>> R = PolynomialRing(Integers(Integer(5)*Integer(2)**Integer(64)), implementation="FLINT", names=('a',)); (a,) = R._first_ngens(1)
Traceback (most recent call last):
...
ValueError: FLINT does not support modulus 92233720368547758080
AUTHORS:
Burcin Erocal (2008-11) initial implementation
Martin Albrecht (2009-01) another initial implementation
- class sage.rings.polynomial.polynomial_zmod_flint.Polynomial_template[source]¶
Bases:
Polynomial
Template for interfacing to external C / C++ libraries for implementations of polynomials.
AUTHORS:
Robert Bradshaw (2008-10): original idea for templating
Martin Albrecht (2008-10): initial implementation
This file implements a simple templating engine for linking univariate polynomials to their C/C++ library implementations. It requires a ‘linkage’ file which implements the
celement_
functions (seesage.libs.ntl.ntl_GF2X_linkage
for an example). Both parts are then plugged together by inclusion of the linkage file when inheriting from this class. Seesage.rings.polynomial.polynomial_gf2x
for an example.We illustrate the generic glueing using univariate polynomials over \(\mathop{\mathrm{GF}}(2)\).
Note
Implementations using this template MUST implement coercion from base ring elements and
get_unsafe()
. SeePolynomial_GF2X
for an example.- degree()[source]¶
EXAMPLES:
sage: P.<x> = GF(2)[] sage: x.degree() 1 sage: P(1).degree() 0 sage: P(0).degree() -1
>>> from sage.all import * >>> P = GF(Integer(2))['x']; (x,) = P._first_ngens(1) >>> x.degree() 1 >>> P(Integer(1)).degree() 0 >>> P(Integer(0)).degree() -1
- gcd(other)[source]¶
Return the greatest common divisor of
self
andother
.EXAMPLES:
sage: P.<x> = GF(2)[] sage: f = x*(x+1) sage: f.gcd(x+1) x + 1 sage: f.gcd(x^2) x
>>> from sage.all import * >>> P = GF(Integer(2))['x']; (x,) = P._first_ngens(1) >>> f = x*(x+Integer(1)) >>> f.gcd(x+Integer(1)) x + 1 >>> f.gcd(x**Integer(2)) x
- is_gen()[source]¶
EXAMPLES:
sage: P.<x> = GF(2)[] sage: x.is_gen() True sage: (x+1).is_gen() False
>>> from sage.all import * >>> P = GF(Integer(2))['x']; (x,) = P._first_ngens(1) >>> x.is_gen() True >>> (x+Integer(1)).is_gen() False
- is_one()[source]¶
EXAMPLES:
sage: P.<x> = GF(2)[] sage: P(1).is_one() True
>>> from sage.all import * >>> P = GF(Integer(2))['x']; (x,) = P._first_ngens(1) >>> P(Integer(1)).is_one() True
- is_zero()[source]¶
EXAMPLES:
sage: P.<x> = GF(2)[] sage: x.is_zero() False
>>> from sage.all import * >>> P = GF(Integer(2))['x']; (x,) = P._first_ngens(1) >>> x.is_zero() False
- list(copy=True)[source]¶
EXAMPLES:
sage: P.<x> = GF(2)[] sage: x.list() [0, 1] sage: list(x) [0, 1]
>>> from sage.all import * >>> P = GF(Integer(2))['x']; (x,) = P._first_ngens(1) >>> x.list() [0, 1] >>> list(x) [0, 1]
- quo_rem(right)[source]¶
EXAMPLES:
sage: P.<x> = GF(2)[] sage: f = x^2 + x + 1 sage: f.quo_rem(x + 1) (x, 1)
>>> from sage.all import * >>> P = GF(Integer(2))['x']; (x,) = P._first_ngens(1) >>> f = x**Integer(2) + x + Integer(1) >>> f.quo_rem(x + Integer(1)) (x, 1)
- shift(n)[source]¶
EXAMPLES:
sage: P.<x> = GF(2)[] sage: f = x^3 + x^2 + 1 sage: f.shift(1) x^4 + x^3 + x sage: f.shift(-1) x^2 + x
>>> from sage.all import * >>> P = GF(Integer(2))['x']; (x,) = P._first_ngens(1) >>> f = x**Integer(3) + x**Integer(2) + Integer(1) >>> f.shift(Integer(1)) x^4 + x^3 + x >>> f.shift(-Integer(1)) x^2 + x
- truncate(n)[source]¶
Return this polynomial mod \(x^n\).
EXAMPLES:
sage: R.<x> =GF(2)[] sage: f = sum(x^n for n in range(10)); f x^9 + x^8 + x^7 + x^6 + x^5 + x^4 + x^3 + x^2 + x + 1 sage: f.truncate(6) x^5 + x^4 + x^3 + x^2 + x + 1
>>> from sage.all import * >>> R = GF(Integer(2))['x']; (x,) = R._first_ngens(1) >>> f = sum(x**n for n in range(Integer(10))); f x^9 + x^8 + x^7 + x^6 + x^5 + x^4 + x^3 + x^2 + x + 1 >>> f.truncate(Integer(6)) x^5 + x^4 + x^3 + x^2 + x + 1
If the precision is higher than the degree of the polynomial then the polynomial itself is returned:
sage: f.truncate(10) is f True
>>> from sage.all import * >>> f.truncate(Integer(10)) is f True
If the precision is negative, the zero polynomial is returned:
sage: f.truncate(-1) 0
>>> from sage.all import * >>> f.truncate(-Integer(1)) 0
- xgcd(other)[source]¶
Compute extended gcd of
self
andother
.EXAMPLES:
sage: P.<x> = GF(7)[] sage: f = x*(x+1) sage: f.xgcd(x+1) (x + 1, 0, 1) sage: f.xgcd(x^2) (x, 1, 6)
>>> from sage.all import * >>> P = GF(Integer(7))['x']; (x,) = P._first_ngens(1) >>> f = x*(x+Integer(1)) >>> f.xgcd(x+Integer(1)) (x + 1, 0, 1) >>> f.xgcd(x**Integer(2)) (x, 1, 6)
- class sage.rings.polynomial.polynomial_zmod_flint.Polynomial_zmod_flint[source]¶
Bases:
Polynomial_template
Polynomial on \(\ZZ/n\ZZ\) implemented via FLINT.
- _add_(right)[source]¶
EXAMPLES:
sage: P.<x> = GF(2)[] sage: x + 1 x + 1
>>> from sage.all import * >>> P = GF(Integer(2))['x']; (x,) = P._first_ngens(1) >>> x + Integer(1) x + 1
- _sub_(right)[source]¶
EXAMPLES:
sage: P.<x> = GF(2)[] sage: x - 1 x + 1
>>> from sage.all import * >>> P = GF(Integer(2))['x']; (x,) = P._first_ngens(1) >>> x - Integer(1) x + 1
- _lmul_(left)[source]¶
EXAMPLES:
sage: P.<x> = GF(2)[] sage: t = x^2 + x + 1 sage: 0*t 0 sage: 1*t x^2 + x + 1 sage: R.<y> = GF(5)[] sage: u = y^2 + y + 1 sage: 3*u 3*y^2 + 3*y + 3 sage: 5*u 0 sage: (2^81)*u 2*y^2 + 2*y + 2 sage: (-2^81)*u 3*y^2 + 3*y + 3
>>> from sage.all import * >>> P = GF(Integer(2))['x']; (x,) = P._first_ngens(1) >>> t = x**Integer(2) + x + Integer(1) >>> Integer(0)*t 0 >>> Integer(1)*t x^2 + x + 1 >>> R = GF(Integer(5))['y']; (y,) = R._first_ngens(1) >>> u = y**Integer(2) + y + Integer(1) >>> Integer(3)*u 3*y^2 + 3*y + 3 >>> Integer(5)*u 0 >>> (Integer(2)**Integer(81))*u 2*y^2 + 2*y + 2 >>> (-Integer(2)**Integer(81))*u 3*y^2 + 3*y + 3
sage: P.<x> = GF(2)[] sage: t = x^2 + x + 1 sage: t*0 0 sage: t*1 x^2 + x + 1 sage: R.<y> = GF(5)[] sage: u = y^2 + y + 1 sage: u*3 3*y^2 + 3*y + 3 sage: u*5 0
>>> from sage.all import * >>> P = GF(Integer(2))['x']; (x,) = P._first_ngens(1) >>> t = x**Integer(2) + x + Integer(1) >>> t*Integer(0) 0 >>> t*Integer(1) x^2 + x + 1 >>> R = GF(Integer(5))['y']; (y,) = R._first_ngens(1) >>> u = y**Integer(2) + y + Integer(1) >>> u*Integer(3) 3*y^2 + 3*y + 3 >>> u*Integer(5) 0
- _rmul_(right)[source]¶
Multiply
self
on the right by a scalar.EXAMPLES:
sage: R.<x> = ZZ[] sage: f = (x^3 + x + 5) sage: f._rmul_(7) 7*x^3 + 7*x + 35 sage: f*7 7*x^3 + 7*x + 35
>>> from sage.all import * >>> R = ZZ['x']; (x,) = R._first_ngens(1) >>> f = (x**Integer(3) + x + Integer(5)) >>> f._rmul_(Integer(7)) 7*x^3 + 7*x + 35 >>> f*Integer(7) 7*x^3 + 7*x + 35
- _mul_(right)[source]¶
EXAMPLES:
sage: P.<x> = GF(2)[] sage: x*(x+1) x^2 + x
>>> from sage.all import * >>> P = GF(Integer(2))['x']; (x,) = P._first_ngens(1) >>> x*(x+Integer(1)) x^2 + x
- _mul_trunc_(right, n)[source]¶
Return the product of this polynomial and other truncated to the given length \(n\).
This function is usually more efficient than simply doing the multiplication and then truncating. The function is tuned for length \(n\) about half the length of a full product.
EXAMPLES:
sage: P.<a> = GF(7)[] sage: a = P(range(10)); b = P(range(5, 15)) sage: a._mul_trunc_(b, 5) 4*a^4 + 6*a^3 + 2*a^2 + 5*a
>>> from sage.all import * >>> P = GF(Integer(7))['a']; (a,) = P._first_ngens(1) >>> a = P(range(Integer(10))); b = P(range(Integer(5), Integer(15))) >>> a._mul_trunc_(b, Integer(5)) 4*a^4 + 6*a^3 + 2*a^2 + 5*a
- compose_mod(other, modulus)[source]¶
Compute \(f(g) \bmod h\).
To be precise about the order fo compostion, given
self
,other
andmodulus
as \(f(x)\), \(g(x)\) and \(h(x)\) compute \(f(g(x)) \bmod h(x)\).INPUT:
other
– a polynomial \(g(x)\)modulus
– a polynomial \(h(x)\)
EXAMPLES:
sage: R.<x> = GF(163)[] sage: f = R.random_element() sage: g = R.random_element() sage: g.compose_mod(g, f) == g(g) % f True sage: f = R([i for i in range(100)]) sage: g = R([i**2 for i in range(100)]) sage: h = 1 + x + x**5 sage: f.compose_mod(g, h) 82*x^4 + 56*x^3 + 45*x^2 + 60*x + 127 sage: f.compose_mod(g, h) == f(g) % h True
>>> from sage.all import * >>> R = GF(Integer(163))['x']; (x,) = R._first_ngens(1) >>> f = R.random_element() >>> g = R.random_element() >>> g.compose_mod(g, f) == g(g) % f True >>> f = R([i for i in range(Integer(100))]) >>> g = R([i**Integer(2) for i in range(Integer(100))]) >>> h = Integer(1) + x + x**Integer(5) >>> f.compose_mod(g, h) 82*x^4 + 56*x^3 + 45*x^2 + 60*x + 127 >>> f.compose_mod(g, h) == f(g) % h True
AUTHORS:
Giacomo Pope (2024-08) initial implementation
- factor()[source]¶
Return the factorization of the polynomial.
EXAMPLES:
sage: R.<x> = GF(5)[] sage: (x^2 + 1).factor() (x + 2) * (x + 3)
>>> from sage.all import * >>> R = GF(Integer(5))['x']; (x,) = R._first_ngens(1) >>> (x**Integer(2) + Integer(1)).factor() (x + 2) * (x + 3)
It also works for prime-power moduli:
sage: R.<x> = Zmod(23^5)[] sage: (x^3 + 1).factor() (x + 1) * (x^2 + 6436342*x + 1)
>>> from sage.all import * >>> R = Zmod(Integer(23)**Integer(5))['x']; (x,) = R._first_ngens(1) >>> (x**Integer(3) + Integer(1)).factor() (x + 1) * (x^2 + 6436342*x + 1)
- is_irreducible()[source]¶
Return whether this polynomial is irreducible.
EXAMPLES:
sage: R.<x> = GF(5)[] sage: (x^2 + 1).is_irreducible() False sage: (x^3 + x + 1).is_irreducible() True
>>> from sage.all import * >>> R = GF(Integer(5))['x']; (x,) = R._first_ngens(1) >>> (x**Integer(2) + Integer(1)).is_irreducible() False >>> (x**Integer(3) + x + Integer(1)).is_irreducible() True
Not implemented when the base ring is not a field:
sage: S.<s> = Zmod(10)[] sage: (s^2).is_irreducible() Traceback (most recent call last): ... NotImplementedError: checking irreducibility of polynomials over rings with composite characteristic is not implemented
>>> from sage.all import * >>> S = Zmod(Integer(10))['s']; (s,) = S._first_ngens(1) >>> (s**Integer(2)).is_irreducible() Traceback (most recent call last): ... NotImplementedError: checking irreducibility of polynomials over rings with composite characteristic is not implemented
- minpoly_mod(other)[source]¶
Thin wrapper for
sage.rings.polynomial.polynomial_modn_dense_ntl.Polynomial_dense_mod_n.minpoly_mod()
.EXAMPLES:
sage: R.<x> = GF(127)[] sage: type(x) <class 'sage.rings.polynomial.polynomial_zmod_flint.Polynomial_zmod_flint'> sage: (x^5 - 3).minpoly_mod(x^3 + 5*x - 1) x^3 + 34*x^2 + 125*x + 95
>>> from sage.all import * >>> R = GF(Integer(127))['x']; (x,) = R._first_ngens(1) >>> type(x) <class 'sage.rings.polynomial.polynomial_zmod_flint.Polynomial_zmod_flint'> >>> (x**Integer(5) - Integer(3)).minpoly_mod(x**Integer(3) + Integer(5)*x - Integer(1)) x^3 + 34*x^2 + 125*x + 95
- modular_composition(other, modulus)[source]¶
Compute \(f(g) \bmod h\).
To be precise about the order fo compostion, given
self
,other
andmodulus
as \(f(x)\), \(g(x)\) and \(h(x)\) compute \(f(g(x)) \bmod h(x)\).INPUT:
other
– a polynomial \(g(x)\)modulus
– a polynomial \(h(x)\)
EXAMPLES:
sage: R.<x> = GF(163)[] sage: f = R.random_element() sage: g = R.random_element() sage: g.compose_mod(g, f) == g(g) % f True sage: f = R([i for i in range(100)]) sage: g = R([i**2 for i in range(100)]) sage: h = 1 + x + x**5 sage: f.compose_mod(g, h) 82*x^4 + 56*x^3 + 45*x^2 + 60*x + 127 sage: f.compose_mod(g, h) == f(g) % h True
>>> from sage.all import * >>> R = GF(Integer(163))['x']; (x,) = R._first_ngens(1) >>> f = R.random_element() >>> g = R.random_element() >>> g.compose_mod(g, f) == g(g) % f True >>> f = R([i for i in range(Integer(100))]) >>> g = R([i**Integer(2) for i in range(Integer(100))]) >>> h = Integer(1) + x + x**Integer(5) >>> f.compose_mod(g, h) 82*x^4 + 56*x^3 + 45*x^2 + 60*x + 127 >>> f.compose_mod(g, h) == f(g) % h True
AUTHORS:
Giacomo Pope (2024-08) initial implementation
- monic()[source]¶
Return this polynomial divided by its leading coefficient.
Raises
ValueError
if the leading coefficient is not invertible in the base ring.EXAMPLES:
sage: R.<x> = GF(5)[] sage: (2*x^2 + 1).monic() x^2 + 3
>>> from sage.all import * >>> R = GF(Integer(5))['x']; (x,) = R._first_ngens(1) >>> (Integer(2)*x**Integer(2) + Integer(1)).monic() x^2 + 3
- rational_reconstruct(*args, **kwds)[source]¶
Deprecated: Use
rational_reconstruction()
instead. See Issue #12696 for details.
- rational_reconstruction(m, n_deg=0, d_deg=0)[source]¶
Construct a rational function \(n/d\) such that \(p*d\) is equivalent to \(n\) modulo \(m\) where \(p\) is this polynomial.
EXAMPLES:
sage: P.<x> = GF(5)[] sage: p = 4*x^5 + 3*x^4 + 2*x^3 + 2*x^2 + 4*x + 2 sage: n, d = p.rational_reconstruction(x^9, 4, 4); n, d (3*x^4 + 2*x^3 + x^2 + 2*x, x^4 + 3*x^3 + x^2 + x) sage: (p*d % x^9) == n True
>>> from sage.all import * >>> P = GF(Integer(5))['x']; (x,) = P._first_ngens(1) >>> p = Integer(4)*x**Integer(5) + Integer(3)*x**Integer(4) + Integer(2)*x**Integer(3) + Integer(2)*x**Integer(2) + Integer(4)*x + Integer(2) >>> n, d = p.rational_reconstruction(x**Integer(9), Integer(4), Integer(4)); n, d (3*x^4 + 2*x^3 + x^2 + 2*x, x^4 + 3*x^3 + x^2 + x) >>> (p*d % x**Integer(9)) == n True
Check that Issue #37169 is fixed - it does not throw an error:
sage: R.<x> = Zmod(4)[] sage: _.<z> = R.quotient_ring(x^2 - 1) sage: c = 2 * z + 1 sage: c * Zmod(2).zero() Traceback (most recent call last): ... RuntimeError: Aborted
>>> from sage.all import * >>> R = Zmod(Integer(4))['x']; (x,) = R._first_ngens(1) >>> _ = R.quotient_ring(x**Integer(2) - Integer(1), names=('z',)); (z,) = _._first_ngens(1) >>> c = Integer(2) * z + Integer(1) >>> c * Zmod(Integer(2)).zero() Traceback (most recent call last): ... RuntimeError: Aborted
- resultant(other)[source]¶
Return the resultant of
self
andother
, which must lie in the same polynomial ring.INPUT:
other
– a polynomial
OUTPUT: an element of the base ring of the polynomial ring
EXAMPLES:
sage: R.<x> = GF(19)['x'] sage: f = x^3 + x + 1; g = x^3 - x - 1 sage: r = f.resultant(g); r 11 sage: r.parent() is GF(19) True
>>> from sage.all import * >>> R = GF(Integer(19))['x']; (x,) = R._first_ngens(1) >>> f = x**Integer(3) + x + Integer(1); g = x**Integer(3) - x - Integer(1) >>> r = f.resultant(g); r 11 >>> r.parent() is GF(Integer(19)) True
The following example shows that Issue #11782 has been fixed:
sage: R.<x> = ZZ.quo(9)['x'] sage: f = 2*x^3 + x^2 + x; g = 6*x^2 + 2*x + 1 sage: f.resultant(g) 5
>>> from sage.all import * >>> R = ZZ.quo(Integer(9))['x']; (x,) = R._first_ngens(1) >>> f = Integer(2)*x**Integer(3) + x**Integer(2) + x; g = Integer(6)*x**Integer(2) + Integer(2)*x + Integer(1) >>> f.resultant(g) 5
- reverse(degree=None)[source]¶
Return a polynomial with the coefficients of this polynomial reversed.
If the optional argument
degree
is given, the coefficient list will be truncated or zero padded as necessary before computing the reverse.EXAMPLES:
sage: R.<x> = GF(5)[] sage: p = R([1,2,3,4]); p 4*x^3 + 3*x^2 + 2*x + 1 sage: p.reverse() x^3 + 2*x^2 + 3*x + 4 sage: p.reverse(degree=6) x^6 + 2*x^5 + 3*x^4 + 4*x^3 sage: p.reverse(degree=2) x^2 + 2*x + 3 sage: R.<x> = GF(101)[] sage: f = x^3 - x + 2; f x^3 + 100*x + 2 sage: f.reverse() 2*x^3 + 100*x^2 + 1 sage: f.reverse() == f(1/x) * x^f.degree() True
>>> from sage.all import * >>> R = GF(Integer(5))['x']; (x,) = R._first_ngens(1) >>> p = R([Integer(1),Integer(2),Integer(3),Integer(4)]); p 4*x^3 + 3*x^2 + 2*x + 1 >>> p.reverse() x^3 + 2*x^2 + 3*x + 4 >>> p.reverse(degree=Integer(6)) x^6 + 2*x^5 + 3*x^4 + 4*x^3 >>> p.reverse(degree=Integer(2)) x^2 + 2*x + 3 >>> R = GF(Integer(101))['x']; (x,) = R._first_ngens(1) >>> f = x**Integer(3) - x + Integer(2); f x^3 + 100*x + 2 >>> f.reverse() 2*x^3 + 100*x^2 + 1 >>> f.reverse() == f(Integer(1)/x) * x**f.degree() True
Note that if \(f\) has zero constant coefficient, its reverse will have lower degree.
sage: f = x^3 + 2*x sage: f.reverse() 2*x^2 + 1
>>> from sage.all import * >>> f = x**Integer(3) + Integer(2)*x >>> f.reverse() 2*x^2 + 1
In this case, reverse is not an involution unless we explicitly specify a degree.
sage: f x^3 + 2*x sage: f.reverse().reverse() x^2 + 2 sage: f.reverse(5).reverse(5) x^3 + 2*x
>>> from sage.all import * >>> f x^3 + 2*x >>> f.reverse().reverse() x^2 + 2 >>> f.reverse(Integer(5)).reverse(Integer(5)) x^3 + 2*x
- revert_series(n)[source]¶
Return a polynomial \(f\) such that
f(self(x)) = self(f(x)) = x
(mod \(x^n\)).EXAMPLES:
sage: R.<t> = GF(5)[] sage: f = t + 2*t^2 - t^3 - 3*t^4 sage: f.revert_series(5) 3*t^4 + 4*t^3 + 3*t^2 + t sage: f.revert_series(-1) Traceback (most recent call last): ... ValueError: argument n must be a nonnegative integer, got -1 sage: g = - t^3 + t^5 sage: g.revert_series(6) Traceback (most recent call last): ... ValueError: self must have constant coefficient 0 and a unit for coefficient t^1 sage: g = t + 2*t^2 - t^3 -3*t^4 + t^5 sage: g.revert_series(6) Traceback (most recent call last): ... ValueError: the integers 1 up to n=5 are required to be invertible over the base field
>>> from sage.all import * >>> R = GF(Integer(5))['t']; (t,) = R._first_ngens(1) >>> f = t + Integer(2)*t**Integer(2) - t**Integer(3) - Integer(3)*t**Integer(4) >>> f.revert_series(Integer(5)) 3*t^4 + 4*t^3 + 3*t^2 + t >>> f.revert_series(-Integer(1)) Traceback (most recent call last): ... ValueError: argument n must be a nonnegative integer, got -1 >>> g = - t**Integer(3) + t**Integer(5) >>> g.revert_series(Integer(6)) Traceback (most recent call last): ... ValueError: self must have constant coefficient 0 and a unit for coefficient t^1 >>> g = t + Integer(2)*t**Integer(2) - t**Integer(3) -Integer(3)*t**Integer(4) + t**Integer(5) >>> g.revert_series(Integer(6)) Traceback (most recent call last): ... ValueError: the integers 1 up to n=5 are required to be invertible over the base field
- small_roots(*args, **kwds)[source]¶
See
sage.rings.polynomial.polynomial_modn_dense_ntl.small_roots()
for the documentation of this function.EXAMPLES:
sage: N = 10001 sage: K = Zmod(10001) sage: P.<x> = PolynomialRing(K) sage: f = x^3 + 10*x^2 + 5000*x - 222 sage: f.small_roots() [4]
>>> from sage.all import * >>> N = Integer(10001) >>> K = Zmod(Integer(10001)) >>> P = PolynomialRing(K, names=('x',)); (x,) = P._first_ngens(1) >>> f = x**Integer(3) + Integer(10)*x**Integer(2) + Integer(5000)*x - Integer(222) >>> f.small_roots() [4]
- squarefree_decomposition()[source]¶
Return the squarefree decomposition of this polynomial.
EXAMPLES:
sage: R.<x> = GF(5)[] sage: ((x+1)*(x^2+1)^2*x^3).squarefree_decomposition() (x + 1) * (x^2 + 1)^2 * x^3
>>> from sage.all import * >>> R = GF(Integer(5))['x']; (x,) = R._first_ngens(1) >>> ((x+Integer(1))*(x**Integer(2)+Integer(1))**Integer(2)*x**Integer(3)).squarefree_decomposition() (x + 1) * (x^2 + 1)^2 * x^3