Multivariate Polynomials via libSINGULAR#
This module implements specialized and optimized implementations for multivariate polynomials over many coefficient rings, via a shared library interface to SINGULAR. In particular, the following coefficient rings are supported by this implementation:
the rational numbers \(\QQ\),
the ring of integers \(\ZZ\),
\(\ZZ/n\ZZ\) for any integer \(n\),
finite fields \(\GF{p^n}\) for \(p\) prime and \(n > 0\),
and absolute number fields \(\QQ(a)\).
EXAMPLES:
We show how to construct various multivariate polynomial rings:
sage: P.<x,y,z> = QQ[]
sage: P
Multivariate Polynomial Ring in x, y, z over Rational Field
sage: f = 27/113 * x^2 + y*z + 1/2; f
27/113*x^2 + y*z + 1/2
sage: P.term_order()
Degree reverse lexicographic term order
sage: P = PolynomialRing(GF(127),3,names='abc', order='lex')
sage: P
Multivariate Polynomial Ring in a, b, c over Finite Field of size 127
sage: a,b,c = P.gens()
sage: f = 57 * a^2*b + 43 * c + 1; f
57*a^2*b + 43*c + 1
sage: P.term_order()
Lexicographic term order
sage: z = QQ['z'].0
sage: K.<s> = NumberField(z^2 - 2)
sage: P.<x,y> = PolynomialRing(K, 2)
sage: 1/2*s*x^2 + 3/4*s
(1/2*s)*x^2 + (3/4*s)
sage: P.<x,y,z> = ZZ[]; P
Multivariate Polynomial Ring in x, y, z over Integer Ring
sage: P.<x,y,z> = Zmod(2^10)[]; P
Multivariate Polynomial Ring in x, y, z over Ring of integers modulo 1024
sage: P.<x,y,z> = Zmod(3^10)[]; P
Multivariate Polynomial Ring in x, y, z over Ring of integers modulo 59049
sage: P.<x,y,z> = Zmod(2^100)[]; P
Multivariate Polynomial Ring in x, y, z over Ring of integers modulo 1267650600228229401496703205376
sage: P.<x,y,z> = Zmod(2521352)[]; P
Multivariate Polynomial Ring in x, y, z over Ring of integers modulo 2521352
sage: type(P)
<class 'sage.rings.polynomial.multi_polynomial_libsingular.MPolynomialRing_libsingular'>
sage: P.<x,y,z> = Zmod(25213521351515232)[]; P
Multivariate Polynomial Ring in x, y, z over Ring of integers modulo 25213521351515232
sage: type(P)
<class 'sage.rings.polynomial.multi_polynomial_ring.MPolynomialRing_polydict_with_category'>
We construct the Frobenius morphism on \(\GF{5}[x,y,z]\) over \(\GF{5}\):
sage: R.<x,y,z> = PolynomialRing(GF(5), 3)
sage: frob = R.hom([x^5, y^5, z^5])
sage: frob(x^2 + 2*y - z^4)
-z^20 + x^10 + 2*y^5
sage: frob((x + 2*y)^3)
x^15 + x^10*y^5 + 2*x^5*y^10 - 2*y^15
sage: (x^5 + 2*y^5)^3
x^15 + x^10*y^5 + 2*x^5*y^10 - 2*y^15
We make a polynomial ring in one variable over a polynomial ring in two variables:
sage: R.<x, y> = PolynomialRing(QQ, 2)
sage: S.<t> = PowerSeriesRing(R)
sage: t*(x+y)
(x + y)*t
Todo
Implement Real, Complex coefficient rings via libSINGULAR
AUTHORS:
Martin Albrecht (2007-01): initial implementation
Joel Mohler (2008-01): misc improvements, polishing
Martin Albrecht (2008-08): added \(\QQ(a)\) and \(\ZZ\) support
Simon King (2009-04): improved coercion
Martin Albrecht (2009-05): added \(\ZZ/n\ZZ\) support, refactoring
Martin Albrecht (2009-06): refactored the code to allow better re-use
Simon King (2011-03): use a faster way of conversion from the base ring.
Volker Braun (2011-06): major cleanup, refcount singular rings, bugfixes.
- class sage.rings.polynomial.multi_polynomial_libsingular.MPolynomialRing_libsingular#
Bases:
MPolynomialRing_base
Construct a multivariate polynomial ring subject to the following conditions:
INPUT:
base_ring
- base ring (must be either GF(q), ZZ, ZZ/nZZ,QQ or absolute number field)
n
- number of variables (must be at least 1)names
- names of ring variables, may be string of list/tupleorder
- term order (default:degrevlex
)
EXAMPLES:
sage: P.<x,y,z> = QQ[] sage: P Multivariate Polynomial Ring in x, y, z over Rational Field sage: f = 27/113 * x^2 + y*z + 1/2; f 27/113*x^2 + y*z + 1/2 sage: P.term_order() Degree reverse lexicographic term order sage: P = PolynomialRing(GF(127),3,names='abc', order='lex') sage: P Multivariate Polynomial Ring in a, b, c over Finite Field of size 127 sage: a,b,c = P.gens() sage: f = 57 * a^2*b + 43 * c + 1; f 57*a^2*b + 43*c + 1 sage: P.term_order() Lexicographic term order sage: z = QQ['z'].0 sage: K.<s> = NumberField(z^2 - 2) sage: P.<x,y> = PolynomialRing(K, 2) sage: 1/2*s*x^2 + 3/4*s (1/2*s)*x^2 + (3/4*s) sage: P.<x,y,z> = ZZ[]; P Multivariate Polynomial Ring in x, y, z over Integer Ring sage: P.<x,y,z> = Zmod(2^10)[]; P Multivariate Polynomial Ring in x, y, z over Ring of integers modulo 1024 sage: P.<x,y,z> = Zmod(3^10)[]; P Multivariate Polynomial Ring in x, y, z over Ring of integers modulo 59049 sage: P.<x,y,z> = Zmod(2^100)[]; P Multivariate Polynomial Ring in x, y, z over Ring of integers modulo 1267650600228229401496703205376 sage: P.<x,y,z> = Zmod(2521352)[]; P Multivariate Polynomial Ring in x, y, z over Ring of integers modulo 2521352 sage: type(P) <class 'sage.rings.polynomial.multi_polynomial_libsingular.MPolynomialRing_libsingular'> sage: P.<x,y,z> = Zmod(25213521351515232)[]; P Multivariate Polynomial Ring in x, y, z over Ring of integers modulo 25213521351515232 sage: type(P) <class 'sage.rings.polynomial.multi_polynomial_ring.MPolynomialRing_polydict_with_category'> sage: P.<x,y,z> = PolynomialRing(Integers(2^32),order='lex') sage: P(2^32-1) 4294967295
- Element#
alias of
MPolynomial_libsingular
- gen(n=0)#
Returns the
n
-th generator of this multivariate polynomial ring.INPUT:
n
– an integer>= 0
EXAMPLES:
sage: P.<x,y,z> = QQ[] sage: P.gen(),P.gen(1) (x, y) sage: P = PolynomialRing(GF(127),1000,'x') sage: P.gen(500) x500 sage: P.<SAGE,SINGULAR> = QQ[] # weird names sage: P.gen(1) SINGULAR
- ideal(*gens, **kwds)#
Create an ideal in this polynomial ring.
INPUT:
*gens
- list or tuple of generators (or several input arguments)coerce
- bool (default:True
); this must be a keyword argument. Only set it toFalse
if you are certain that each generator is already in the ring.
EXAMPLES:
sage: P.<x,y,z> = QQ[] sage: sage.rings.ideal.Katsura(P) Ideal (x + 2*y + 2*z - 1, x^2 + 2*y^2 + 2*z^2 - x, 2*x*y + 2*y*z - y) of Multivariate Polynomial Ring in x, y, z over Rational Field sage: P.ideal([x + 2*y + 2*z-1, 2*x*y + 2*y*z-y, x^2 + 2*y^2 + 2*z^2-x]) Ideal (x + 2*y + 2*z - 1, 2*x*y + 2*y*z - y, x^2 + 2*y^2 + 2*z^2 - x) of Multivariate Polynomial Ring in x, y, z over Rational Field
- monomial_all_divisors(t)#
Return a list of all monomials that divide
t
.Coefficients are ignored.
INPUT:
t
- a monomial
- OUTPUT:
a list of monomials
EXAMPLES:
sage: P.<x,y,z> = QQ[] sage: P.monomial_all_divisors(x^2*z^3) [x, x^2, z, x*z, x^2*z, z^2, x*z^2, x^2*z^2, z^3, x*z^3, x^2*z^3]
ALGORITHM: addwithcarry idea by Toon Segers
- monomial_divides(a, b)#
Return
False
if a does not divide b andTrue
otherwise.Coefficients are ignored.
INPUT:
a
– monomialb
– monomial
EXAMPLES:
sage: P.<x,y,z> = QQ[] sage: P.monomial_divides(x*y*z, x^3*y^2*z^4) True sage: P.monomial_divides(x^3*y^2*z^4, x*y*z) False
- monomial_lcm(f, g)#
LCM for monomials. Coefficients are ignored.
INPUT:
f
- monomialg
- monomial
EXAMPLES:
sage: P.<x,y,z> = QQ[] sage: P.monomial_lcm(3/2*x*y,x) x*y
- monomial_pairwise_prime(g, h)#
Return
True
ifh
andg
are pairwise prime. Both are treated as monomials.Coefficients are ignored.
INPUT:
h
- monomialg
- monomial
EXAMPLES:
sage: P.<x,y,z> = QQ[] sage: P.monomial_pairwise_prime(x^2*z^3, y^4) True sage: P.monomial_pairwise_prime(1/2*x^3*y^2, 3/4*y^3) False
- monomial_quotient(f, g, coeff=False)#
Return
f/g
, where bothf
and``g
are treated as monomials.Coefficients are ignored by default.
INPUT:
f
- monomialg
- monomialcoeff
- divide coefficients as well (default:False
)
EXAMPLES:
sage: P.<x,y,z> = QQ[] sage: P.monomial_quotient(3/2*x*y,x) y sage: P.monomial_quotient(3/2*x*y,x,coeff=True) 3/2*y
Note, that \(\ZZ\) behaves different if
coeff=True
:sage: P.monomial_quotient(2*x,3*x) 1 sage: P.<x,y> = PolynomialRing(ZZ) sage: P.monomial_quotient(2*x,3*x,coeff=True) Traceback (most recent call last): ... ArithmeticError: Cannot divide these coefficients.
Warning
Assumes that the head term of f is a multiple of the head term of g and return the multiplicant m. If this rule is violated, funny things may happen.
- monomial_reduce(f, G)#
Try to find a
g
inG
whereg.lm()
dividesf
. If found(flt,g)
is returned,(0,0)
otherwise, whereflt
isf/g.lm()
.It is assumed that
G
is iterable and contains only elements in this polynomial ring.Coefficients are ignored.
INPUT:
f
- monomialG
- list/set of mpolynomials
EXAMPLES:
sage: P.<x,y,z> = QQ[] sage: f = x*y^2 sage: G = [ 3/2*x^3 + y^2 + 1/2, 1/4*x*y + 2/7, 1/2 ] sage: P.monomial_reduce(f,G) (y, 1/4*x*y + 2/7)
- ngens()#
Returns the number of variables in this multivariate polynomial ring.
EXAMPLES:
sage: P.<x,y> = QQ[] sage: P.ngens() 2 sage: k.<a> = GF(2^16) sage: P = PolynomialRing(k,1000,'x') sage: P.ngens() 1000
- class sage.rings.polynomial.multi_polynomial_libsingular.MPolynomial_libsingular#
Bases:
MPolynomial
A multivariate polynomial implemented using libSINGULAR.
- add_m_mul_q(m, q)#
Return
self + m*q
, wherem
must be a monomial andq
a polynomial.INPUT:
m
- a monomialq
- a polynomial
EXAMPLES:
sage: P.<x,y,z>=PolynomialRing(QQ,3) sage: x.add_m_mul_q(y,z) y*z + x
- coefficient(degrees)#
Return the coefficient of the variables with the degrees specified in the python dictionary
degrees
. Mathematically, this is the coefficient in the base ring adjoined by the variables of this ring not listed indegrees
. However, the result has the same parent as this polynomial.This function contrasts with the function
monomial_coefficient
which returns the coefficient in the base ring of a monomial.INPUT:
degrees
- Can be any of:a dictionary of degree restrictions
a list of degree restrictions (with None in the unrestricted variables)
a monomial (very fast, but not as flexible)
- OUTPUT:
element of the parent of this element.
Note
For coefficients of specific monomials, look at
monomial_coefficient()
.EXAMPLES:
sage: R.<x,y> = QQ[] sage: f=x*y+y+5 sage: f.coefficient({x:0,y:1}) 1 sage: f.coefficient({x:0}) y + 5 sage: f=(1+y+y^2)*(1+x+x^2) sage: f.coefficient({x:0}) y^2 + y + 1 sage: f.coefficient([0,None]) y^2 + y + 1 sage: f.coefficient(x) y^2 + y + 1
Note that exponents have all variables specified:
sage: x.coefficient(x.exponents()[0]) 1 sage: f.coefficient([1,0]) 1 sage: f.coefficient({x:1,y:0}) 1
Be aware that this may not be what you think! The physical appearance of the variable x is deceiving – particularly if the exponent would be a variable.
sage: f.coefficient(x^0) # outputs the full polynomial x^2*y^2 + x^2*y + x*y^2 + x^2 + x*y + y^2 + x + y + 1 sage: R.<x,y> = GF(389)[] sage: f=x*y+5 sage: c=f.coefficient({x:0,y:0}); c 5 sage: parent(c) Multivariate Polynomial Ring in x, y over Finite Field of size 389
AUTHOR:
Joel B. Mohler (2007.10.31)
- coefficients()#
Return the nonzero coefficients of this polynomial in a list. The returned list is decreasingly ordered by the term ordering of the parent.
EXAMPLES:
sage: R.<x,y,z> = PolynomialRing(QQ, order='degrevlex') sage: f=23*x^6*y^7 + x^3*y+6*x^7*z sage: f.coefficients() [23, 6, 1] sage: R.<x,y,z> = PolynomialRing(QQ, order='lex') sage: f=23*x^6*y^7 + x^3*y+6*x^7*z sage: f.coefficients() [6, 23, 1]
AUTHOR:
Didier Deshommes
- constant_coefficient()#
Return the constant coefficient of this multivariate polynomial.
EXAMPLES:
sage: P.<x, y> = QQ[] sage: f = 3*x^2 - 2*y + 7*x^2*y^2 + 5 sage: f.constant_coefficient() 5 sage: f = 3*x^2 sage: f.constant_coefficient() 0
- degree(x=None, std_grading=False)#
Return the degree of this polynomial.
INPUT:
x
– (default:None
) a generator of the parent ring
OUTPUT:
If
x
is not given, return the maximum degree of the monomials of the polynomial. Note that the degree of a monomial is affected by the gradings given to the generators of the parent ring. Ifx
is given, it is (or coercible to) a generator of the parent ring and the output is the maximum degree inx
. This is not affected by the gradings of the generators.EXAMPLES:
sage: R.<x, y> = QQ[] sage: f = y^2 - x^9 - x sage: f.degree(x) 9 sage: f.degree(y) 2 sage: (y^10*x - 7*x^2*y^5 + 5*x^3).degree(x) 3 sage: (y^10*x - 7*x^2*y^5 + 5*x^3).degree(y) 10
The term ordering of the parent ring determines the grading of the generators.
sage: T = TermOrder('wdegrevlex', (1,2,3,4)) sage: R = PolynomialRing(QQ, 'x', 12, order=T+T+T) sage: [x.degree() for x in R.gens()] [1, 2, 3, 4, 1, 2, 3, 4, 1, 2, 3, 4]
A matrix term ordering determines the grading of the generators by the first row of the matrix.
sage: m = matrix(3, [3,2,1,1,1,0,1,0,0]) sage: m [3 2 1] [1 1 0] [1 0 0] sage: R.<x,y,z> = PolynomialRing(QQ, order=TermOrder(m)) sage: x.degree(), y.degree(), z.degree() (3, 2, 1) sage: f = x^3*y + x*z^4 sage: f.degree() 11
If the first row contains zero, the grading becomes the standard one.
sage: m = matrix(3, [3,0,1,1,1,0,1,0,0]) sage: m [3 0 1] [1 1 0] [1 0 0] sage: R.<x,y,z> = PolynomialRing(QQ, order=TermOrder(m)) sage: x.degree(), y.degree(), z.degree() (1, 1, 1) sage: f = x^3*y + x*z^4 sage: f.degree() 5
To get the degree with the standard grading regardless of the term ordering of the parent ring, use
std_grading=True
.sage: f.degree(std_grading=True) 5
- degrees()#
Returns a tuple with the maximal degree of each variable in this polynomial. The list of degrees is ordered by the order of the generators.
EXAMPLES:
sage: R.<y0,y1,y2> = PolynomialRing(QQ,3) sage: q = 3*y0*y1*y1*y2; q 3*y0*y1^2*y2 sage: q.degrees() (1, 2, 1) sage: (q + y0^5).degrees() (5, 2, 1)
- dict()#
Return a dictionary representing self. This dictionary is in the same format as the generic MPolynomial: The dictionary consists of
ETuple:coefficient
pairs.EXAMPLES:
sage: R.<x,y,z> = QQ[] sage: f=2*x*y^3*z^2 + 1/7*x^2 + 2/3 sage: f.dict() {(0, 0, 0): 2/3, (1, 3, 2): 2, (2, 0, 0): 1/7}
- divides(other)#
Return
True
if this polynomial dividesother
.EXAMPLES:
sage: R.<x,y,z> = QQ[] sage: p = 3*x*y + 2*y*z + x*z sage: q = x + y + z + 1 sage: r = p * q sage: p.divides(r) True sage: q.divides(p) False sage: r.divides(0) True sage: R.zero().divides(r) False sage: R.zero().divides(0) True
- exponents(as_ETuples=True)#
Return the exponents of the monomials appearing in this polynomial.
INPUT:
as_ETuples
– (default:True
) ifTrue
returns the result as an list of ETuples, otherwise returns a list of tuples
EXAMPLES:
sage: R.<a,b,c> = QQ[] sage: f = a^3 + b + 2*b^2 sage: f.exponents() [(3, 0, 0), (0, 2, 0), (0, 1, 0)] sage: f.exponents(as_ETuples=False) [(3, 0, 0), (0, 2, 0), (0, 1, 0)]
- factor(proof=None)#
Return the factorization of this polynomial.
INPUT:
proof
- ignored.
EXAMPLES:
sage: R.<x, y> = QQ[] sage: f = (x^3 + 2*y^2*x) * (x^2 + x + 1); f x^5 + 2*x^3*y^2 + x^4 + 2*x^2*y^2 + x^3 + 2*x*y^2 sage: F = f.factor() sage: F x * (x^2 + x + 1) * (x^2 + 2*y^2)
Next we factor the same polynomial, but over the finite field of order 3.:
sage: R.<x, y> = GF(3)[] sage: f = (x^3 + 2*y^2*x) * (x^2 + x + 1); f x^5 - x^3*y^2 + x^4 - x^2*y^2 + x^3 - x*y^2 sage: F = f.factor() sage: F # order is somewhat random (-1) * x * (-x + y) * (x + y) * (x - 1)^2
Next we factor a polynomial, but over a finite field of order 9.:
sage: K.<a> = GF(3^2) sage: R.<x, y> = K[] sage: f = (x^3 + 2*a*y^2*x) * (x^2 + x + 1); f x^5 + (-a)*x^3*y^2 + x^4 + (-a)*x^2*y^2 + x^3 + (-a)*x*y^2 sage: F = f.factor() sage: F ((-a)) * x * (x - 1)^2 * ((-a + 1)*x^2 + y^2) sage: f - F 0
Next we factor a polynomial over a number field.:
sage: p = var('p') sage: K.<s> = NumberField(p^3-2) sage: KXY.<x,y> = K[] sage: factor(x^3 - 2*y^3) (x + (-s)*y) * (x^2 + s*x*y + (s^2)*y^2) sage: k = (x^3-2*y^3)^5*(x+s*y)^2*(2/3 + s^2) sage: k.factor() ((s^2 + 2/3)) * (x + s*y)^2 * (x + (-s)*y)^5 * (x^2 + s*x*y + (s^2)*y^2)^5
This shows that ticket trac ticket #2780 is fixed, i.e. that the unit part of the factorization is set correctly:
sage: x = var('x') sage: K.<a> = NumberField(x^2 + 1) sage: R.<y, z> = PolynomialRing(K) sage: f = 2*y^2 + 2*z^2 sage: F = f.factor(); F.unit() 2
Another example:
sage: R.<x,y,z> = GF(32003)[] sage: f = 9*(x-1)^2*(y+z) sage: f.factor() (9) * (y + z) * (x - 1)^2 sage: R.<x,w,v,u> = QQ['x','w','v','u'] sage: p = (4*v^4*u^2 - 16*v^2*u^4 + 16*u^6 - 4*v^4*u + 8*v^2*u^3 + v^4) sage: p.factor() (-2*v^2*u + 4*u^3 + v^2)^2 sage: R.<a,b,c,d> = QQ[] sage: f = (-2) * (a - d) * (-a + b) * (b - d) * (a - c) * (b - c) * (c - d) sage: F = f.factor(); F (-2) * (c - d) * (-b + c) * (b - d) * (-a + c) * (-a + b) * (a - d) sage: F[0][0] c - d sage: F.unit() -2
Constant elements are factorized in the base rings.
sage: P.<x,y> = ZZ[] sage: P(2^3*7).factor() 2^3 * 7 sage: P.<x,y> = GF(2)[] sage: P(1).factor() 1
Factorization for finite prime fields with characteristic \(> 2^{29}\) is not supported
sage: q = 1073741789 sage: T.<aa, bb> = PolynomialRing(GF(q)) sage: f = aa^2 + 12124343*bb*aa + 32434598*bb^2 sage: f.factor() Traceback (most recent call last): ... NotImplementedError: Factorization of multivariate polynomials over prime fields with characteristic > 2^29 is not implemented.
Factorization over the integers is now supported, see trac ticket #17840:
sage: P.<x,y> = PolynomialRing(ZZ) sage: f = 12 * (3*x*y + 4) * (5*x - 2) * (2*y + 7)^2 sage: f.factor() 2^2 * 3 * (2*y + 7)^2 * (5*x - 2) * (3*x*y + 4) sage: g = -12 * (x^2 - y^2) sage: g.factor() (-1) * 2^2 * 3 * (x - y) * (x + y) sage: factor(-4*x*y - 2*x + 2*y + 1) (-1) * (2*y + 1) * (2*x - 1)
Factorization over non-integral domains is not supported
sage: R.<x,y> = PolynomialRing(Zmod(4)) sage: f = (2*x + 1) * (x^2 + x + 1) sage: f.factor() Traceback (most recent call last): ... NotImplementedError: Factorization of multivariate polynomials over Ring of integers modulo 4 is not implemented.
- gcd(right, algorithm=None, **kwds)#
Return the greatest common divisor of self and right.
INPUT:
right
- polynomialalgorithm
-ezgcd
- EZGCD algorithm -modular
- multi-modular algorithm (default)**kwds
- ignored
EXAMPLES:
sage: P.<x,y,z> = QQ[] sage: f = (x*y*z)^6 - 1 sage: g = (x*y*z)^4 - 1 sage: f.gcd(g) x^2*y^2*z^2 - 1 sage: GCD([x^3 - 3*x + 2, x^4 - 1, x^6 -1]) x - 1 sage: R.<x,y> = QQ[] sage: f = (x^3 + 2*y^2*x)^2 sage: g = x^2*y^2 sage: f.gcd(g) x^2
We compute a gcd over a finite field:
sage: F.<u> = GF(31^2) sage: R.<x,y,z> = F[] sage: p = x^3 + (1+u)*y^3 + z^3 sage: q = p^3 * (x - y + z*u) sage: gcd(p,q) x^3 + (u + 1)*y^3 + z^3 sage: gcd(p,q) # yes, twice -- tests that singular ring is properly set. x^3 + (u + 1)*y^3 + z^3
We compute a gcd over a number field:
sage: x = polygen(QQ) sage: F.<u> = NumberField(x^3 - 2) sage: R.<x,y,z> = F[] sage: p = x^3 + (1+u)*y^3 + z^3 sage: q = p^3 * (x - y + z*u) sage: gcd(p,q) x^3 + (u + 1)*y^3 + z^3
- global_height(prec=None)#
Return the (projective) global height of the polynomial.
This returns the absolute logarithmic height of the coefficients thought of as a projective point.
INPUT:
prec
– desired floating point precision (default: default RealField precision).
OUTPUT:
a real number.
EXAMPLES:
sage: R.<x,y> = PolynomialRing(QQ) sage: f = 3*x^3 + 2*x*y^2 sage: exp(f.global_height()) 3.00000000000000
sage: K.<k> = CyclotomicField(3) sage: R.<x,y> = PolynomialRing(K, sparse=True) sage: f = k*x*y + 1 sage: exp(f.global_height()) 1.00000000000000
Scaling should not change the result:
sage: R.<x,y> = PolynomialRing(QQ) sage: f = 1/25*x^2 + 25/3*x*y + y^2 sage: f.global_height() 6.43775164973640 sage: g = 100 * f sage: g.global_height() 6.43775164973640
sage: R.<x> = PolynomialRing(QQ) sage: K.<k> = NumberField(x^2 + 5) sage: T.<t,w> = PolynomialRing(K) sage: f = 1/1331 * t^2 + 5 * w + 7 sage: f.global_height() 9.13959596745043
sage: R.<x,y> = QQ[] sage: f = 1/123*x*y + 12 sage: f.global_height(prec=2) 8.0
sage: R.<x,y> = QQ[] sage: f = 0*x*y sage: f.global_height() 0.000000000000000
- gradient()#
Return a list of partial derivatives of this polynomial, ordered by the variables of the parent.
EXAMPLES:
sage: P.<x,y,z> = PolynomialRing(QQ,3) sage: f= x*y + 1 sage: f.gradient() [y, x, 0]
- hamming_weight()#
Return the number of non-zero coefficients of this polynomial.
This is also called weight,
hamming_weight()
or sparsity.EXAMPLES:
sage: R.<x, y> = ZZ[] sage: f = x^3 - y sage: f.number_of_terms() 2 sage: R(0).number_of_terms() 0 sage: f = (x+y)^100 sage: f.number_of_terms() 101
The method
hamming_weight()
is an alias:sage: f.hamming_weight() 101
- integral(var)#
Integrate this polynomial with respect to the provided variable.
One requires that \(\QQ\) is contained in the ring.
INPUT:
variable
– the integral is taken with respect to variable
EXAMPLES:
sage: R.<x, y> = PolynomialRing(QQ, 2) sage: f = 3*x^3*y^2 + 5*y^2 + 3*x + 2 sage: f.integral(x) 3/4*x^4*y^2 + 5*x*y^2 + 3/2*x^2 + 2*x sage: f.integral(y) x^3*y^3 + 5/3*y^3 + 3*x*y + 2*y
Check that trac ticket #15896 is solved:
sage: s = x+y sage: s.integral(x)+x 1/2*x^2 + x*y + x sage: s.integral(x)*s 1/2*x^3 + 3/2*x^2*y + x*y^2
- inverse_of_unit()#
Return the inverse of this polynomial if it is a unit.
EXAMPLES:
sage: R.<x,y> = QQ[] sage: x.inverse_of_unit() Traceback (most recent call last): ... ArithmeticError: Element is not a unit. sage: R(1/2).inverse_of_unit() 2
- is_constant()#
Return
True
if this polynomial is constant.EXAMPLES:
sage: P.<x,y,z> = PolynomialRing(GF(127)) sage: x.is_constant() False sage: P(1).is_constant() True
- is_homogeneous()#
Return
True
if this polynomial is homogeneous.EXAMPLES:
sage: P.<x,y> = PolynomialRing(RationalField(), 2) sage: (x+y).is_homogeneous() True sage: (x.parent()(0)).is_homogeneous() True sage: (x+y^2).is_homogeneous() False sage: (x^2 + y^2).is_homogeneous() True sage: (x^2 + y^2*x).is_homogeneous() False sage: (x^2*y + y^2*x).is_homogeneous() True
- is_monomial()#
Return
True
if this polynomial is a monomial. A monomial is defined to be a product of generators with coefficient 1.EXAMPLES:
sage: P.<x,y,z> = PolynomialRing(QQ) sage: x.is_monomial() True sage: (2*x).is_monomial() False sage: (x*y).is_monomial() True sage: (x*y + x).is_monomial() False sage: P(2).is_monomial() False sage: P.zero().is_monomial() False
- is_squarefree()#
Return
True
if this polynomial is square free.EXAMPLES:
sage: P.<x,y,z> = PolynomialRing(QQ) sage: f= x^2 + 2*x*y + 1/2*z sage: f.is_squarefree() True sage: h = f^2 sage: h.is_squarefree() False
- is_term()#
Return
True
ifself
is a term, which we define to be a product of generators times some coefficient, which need not be 1.Use
is_monomial()
to require that the coefficient be 1.EXAMPLES:
sage: P.<x,y,z> = PolynomialRing(QQ) sage: x.is_term() True sage: (2*x).is_term() True sage: (x*y).is_term() True sage: (x*y + x).is_term() False sage: P(2).is_term() True sage: P.zero().is_term() False
- is_univariate()#
Return
True
if self is a univariate polynomial, that is if self contains only one variable.EXAMPLES:
sage: P.<x,y,z> = GF(2)[] sage: f = x^2 + 1 sage: f.is_univariate() True sage: f = y*x^2 + 1 sage: f.is_univariate() False sage: f = P(0) sage: f.is_univariate() True
- is_zero()#
Return
True
if this polynomial is zero.EXAMPLES:
sage: P.<x,y> = PolynomialRing(QQ) sage: x.is_zero() False sage: (x-x).is_zero() True
- iterator_exp_coeff(as_ETuples=True)#
Iterate over
self
as pairs of ((E)Tuple, coefficient).INPUT:
as_ETuples
– (default:True
) ifTrue
iterate over pairs whose first element is an ETuple, otherwise as a tuples
EXAMPLES:
sage: R.<a,b,c> = QQ[] sage: f = a*c^3 + a^2*b + 2*b^4 sage: list(f.iterator_exp_coeff()) [((0, 4, 0), 2), ((1, 0, 3), 1), ((2, 1, 0), 1)] sage: list(f.iterator_exp_coeff(as_ETuples=False)) [((0, 4, 0), 2), ((1, 0, 3), 1), ((2, 1, 0), 1)] sage: R.<a,b,c> = PolynomialRing(QQ, 3, order='lex') sage: f = a*c^3 + a^2*b + 2*b^4 sage: list(f.iterator_exp_coeff()) [((2, 1, 0), 1), ((1, 0, 3), 1), ((0, 4, 0), 2)]
- lc()#
Leading coefficient of this polynomial with respect to the term order of
self.parent()
.EXAMPLES:
sage: R.<x,y,z>=PolynomialRing(GF(7),3,order='lex') sage: f = 3*x^1*y^2 + 2*y^3*z^4 sage: f.lc() 3 sage: f = 5*x^3*y^2*z^4 + 4*x^3*y^2*z^1 sage: f.lc() 5
- lcm(g)#
Return the least common multiple of
self
and \(g\).EXAMPLES:
sage: P.<x,y,z> = QQ[] sage: p = (x+y)*(y+z) sage: q = (z^4+2)*(y+z) sage: lcm(p,q) x*y*z^4 + y^2*z^4 + x*z^5 + y*z^5 + 2*x*y + 2*y^2 + 2*x*z + 2*y*z sage: P.<x,y,z> = ZZ[] sage: p = 2*(x+y)*(y+z) sage: q = 3*(z^4+2)*(y+z) sage: lcm(p,q) 6*x*y*z^4 + 6*y^2*z^4 + 6*x*z^5 + 6*y*z^5 + 12*x*y + 12*y^2 + 12*x*z + 12*y*z sage: r.<x,y> = PolynomialRing(GF(2**8, 'a'), 2) sage: a = r.base_ring().0 sage: f = (a^2+a)*x^2*y + (a^4+a^3+a)*y + a^5 sage: f.lcm(x^4) (a^2 + a)*x^6*y + (a^4 + a^3 + a)*x^4*y + (a^5)*x^4 sage: w = var('w') sage: r.<x,y> = PolynomialRing(NumberField(w^4 + 1, 'a'), 2) sage: a = r.base_ring().0 sage: f = (a^2+a)*x^2*y + (a^4+a^3+a)*y + a^5 sage: f.lcm(x^4) (a^2 + a)*x^6*y + (a^3 + a - 1)*x^4*y + (-a)*x^4
- lift(I)#
given an ideal
I = (f_1,...,f_r)
and someg (== self)
inI
, finds_1,...,s_r
such thatg = s_1 f_1 + ... + s_r f_r
.A
ValueError
exception is raised ifg (== self)
does not belong toI
.EXAMPLES:
sage: A.<x,y> = PolynomialRing(QQ,2,order='degrevlex') sage: I = A.ideal([x^10 + x^9*y^2, y^8 - x^2*y^7 ]) sage: f = x*y^13 + y^12 sage: M = f.lift(I) sage: M [y^7, x^7*y^2 + x^8 + x^5*y^3 + x^6*y + x^3*y^4 + x^4*y^2 + x*y^5 + x^2*y^3 + y^4] sage: sum( map( mul , zip( M, I.gens() ) ) ) == f True
Check that trac ticket #13671 is fixed:
sage: R.<x1,x2> = QQ[] sage: I = R.ideal(x2**2 + x1 - 2, x1**2 - 1) sage: f = I.gen(0) + x2*I.gen(1) sage: f.lift(I) [1, x2] sage: (f+1).lift(I) Traceback (most recent call last): ... ValueError: polynomial is not in the ideal sage: f.lift(I) [1, x2]
- lm()#
Returns the lead monomial of self with respect to the term order of
self.parent()
. In Sage a monomial is a product of variables in some power without a coefficient.EXAMPLES:
sage: R.<x,y,z>=PolynomialRing(GF(7),3,order='lex') sage: f = x^1*y^2 + y^3*z^4 sage: f.lm() x*y^2 sage: f = x^3*y^2*z^4 + x^3*y^2*z^1 sage: f.lm() x^3*y^2*z^4 sage: R.<x,y,z>=PolynomialRing(QQ,3,order='deglex') sage: f = x^1*y^2*z^3 + x^3*y^2*z^0 sage: f.lm() x*y^2*z^3 sage: f = x^1*y^2*z^4 + x^1*y^1*z^5 sage: f.lm() x*y^2*z^4 sage: R.<x,y,z>=PolynomialRing(GF(127),3,order='degrevlex') sage: f = x^1*y^5*z^2 + x^4*y^1*z^3 sage: f.lm() x*y^5*z^2 sage: f = x^4*y^7*z^1 + x^4*y^2*z^3 sage: f.lm() x^4*y^7*z
- local_height(v, prec=None)#
Return the maximum of the local height of the coefficients of this polynomial.
INPUT:
v
– a prime or prime ideal of the base ring.prec
– desired floating point precision (default: default RealField precision).
OUTPUT:
a real number.
EXAMPLES:
sage: R.<x,y> = PolynomialRing(QQ) sage: f = 1/1331*x^2 + 1/4000*y^2 sage: f.local_height(1331) 7.19368581839511
sage: R.<x> = QQ[] sage: K.<k> = NumberField(x^2 - 5) sage: T.<t,w> = K[] sage: I = K.ideal(3) sage: f = 1/3*t*w + 3 sage: f.local_height(I) 1.09861228866811
sage: R.<x,y> = QQ[] sage: f = 1/2*x*y + 2 sage: f.local_height(2, prec=2) 0.75
- local_height_arch(i, prec=None)#
Return the maximum of the local height at the
i
-th infinite place of the coefficients of this polynomial.INPUT:
i
– an integer.prec
– desired floating point precision (default: default RealField precision).
OUTPUT:
a real number.
EXAMPLES:
sage: R.<x,y> = PolynomialRing(QQ) sage: f = 210*x*y sage: f.local_height_arch(0) 5.34710753071747
sage: R.<x> = QQ[] sage: K.<k> = NumberField(x^2 - 5) sage: T.<t,w> = K[] sage: f = 1/2*t*w + 3 sage: f.local_height_arch(1, prec=52) 1.09861228866811
sage: R.<x,y> = QQ[] sage: f = 1/2*x*y + 3 sage: f.local_height_arch(0, prec=2) 1.0
- lt()#
Leading term of this polynomial. In Sage a term is a product of variables in some power and a coefficient.
EXAMPLES:
sage: R.<x,y,z>=PolynomialRing(GF(7),3,order='lex') sage: f = 3*x^1*y^2 + 2*y^3*z^4 sage: f.lt() 3*x*y^2 sage: f = 5*x^3*y^2*z^4 + 4*x^3*y^2*z^1 sage: f.lt() -2*x^3*y^2*z^4
- monomial_coefficient(mon)#
Return the coefficient in the base ring of the monomial mon in
self
, where mon must have the same parent as self.This function contrasts with the function
coefficient
which returns the coefficient of a monomial viewing this polynomial in a polynomial ring over a base ring having fewer variables.INPUT:
mon
- a monomial
OUTPUT:
coefficient in base ring
See also
For coefficients in a base ring of fewer variables, look at
coefficient
.EXAMPLES:
sage: P.<x,y> = QQ[] The parent of the return is a member of the base ring. sage: f = 2 * x * y sage: c = f.monomial_coefficient(x*y); c 2 sage: c.parent() Rational Field sage: f = y^2 + y^2*x - x^9 - 7*x + 5*x*y sage: f.monomial_coefficient(y^2) 1 sage: f.monomial_coefficient(x*y) 5 sage: f.monomial_coefficient(x^9) -1 sage: f.monomial_coefficient(x^10) 0
- monomials()#
Return the list of monomials in self. The returned list is decreasingly ordered by the term ordering of
self.parent()
.EXAMPLES:
sage: P.<x,y,z> = QQ[] sage: f = x + 3/2*y*z^2 + 2/3 sage: f.monomials() [y*z^2, x, 1] sage: f = P(3/2) sage: f.monomials() [1]
- number_of_terms()#
Return the number of non-zero coefficients of this polynomial.
This is also called weight,
hamming_weight()
or sparsity.EXAMPLES:
sage: R.<x, y> = ZZ[] sage: f = x^3 - y sage: f.number_of_terms() 2 sage: R(0).number_of_terms() 0 sage: f = (x+y)^100 sage: f.number_of_terms() 101
The method
hamming_weight()
is an alias:sage: f.hamming_weight() 101
- numerator()#
Return a numerator of self computed as self * self.denominator()
If the base_field of self is the Rational Field then the numerator is a polynomial whose base_ring is the Integer Ring, this is done for compatibility to the univariate case.
Warning
This is not the numerator of the rational function defined by self, which would always be self since self is a polynomial.
EXAMPLES:
First we compute the numerator of a polynomial with integer coefficients, which is of course self.
sage: R.<x, y> = ZZ[] sage: f = x^3 + 17*y + 1 sage: f.numerator() x^3 + 17*y + 1 sage: f == f.numerator() True
Next we compute the numerator of a polynomial with rational coefficients.
sage: R.<x,y> = PolynomialRing(QQ) sage: f = (1/17)*x^19 - (2/3)*y + 1/3; f 1/17*x^19 - 2/3*y + 1/3 sage: f.numerator() 3*x^19 - 34*y + 17 sage: f == f.numerator() False sage: f.numerator().base_ring() Integer Ring
We check that the computation of numerator and denominator is valid.
sage: K=QQ['x,y'] sage: f=K.random_element() sage: f.numerator() / f.denominator() == f True
The following tests against a bug fixed in trac ticket #11780:
sage: P.<foo,bar> = ZZ[] sage: Q.<foo,bar> = QQ[] sage: f = Q.random_element() sage: f.numerator().parent() is P True
- nvariables()#
Return the number variables in this polynomial.
EXAMPLES:
sage: P.<x,y,z> = PolynomialRing(GF(127)) sage: f = x*y + z sage: f.nvariables() 3 sage: f = x + y sage: f.nvariables() 2
- quo_rem(right)#
Returns quotient and remainder of self and right.
EXAMPLES:
sage: R.<x,y> = QQ[] sage: f = y*x^2 + x + 1 sage: f.quo_rem(x) (x*y + 1, 1) sage: f.quo_rem(y) (x^2, x + 1) sage: R.<x,y> = ZZ[] sage: f = 2*y*x^2 + x + 1 sage: f.quo_rem(x) (2*x*y + 1, 1) sage: f.quo_rem(y) (2*x^2, x + 1) sage: f.quo_rem(3*x) (0, 2*x^2*y + x + 1)
- reduce(I)#
Return a remainder of this polynomial modulo the polynomials in
I
.INPUT:
I
- an ideal or a list/set/iterable of polynomials.
OUTPUT:
A polynomial
r
such that:self
-r
is in the ideal generated byI
.No term in
r
is divisible by any of the leading monomials ofI
.
The result
r
is canonical if:I
is an ideal, and Sage can compute a Groebner basis of it.I
is a list/set/iterable that is a (strong) Groebner basis for the term order ofself
. (A strong Groebner basis is such that for every leading termt
of the ideal generated byI
, there exists an elementg
ofI
such that the leading term ofg
dividest
.)
The result
r
is implementation-dependent (and possibly order-dependent) otherwise. IfI
is an ideal and no Groebner basis can be computed, its list of generatorsI.gens()
is used for the reduction.EXAMPLES:
sage: P.<x,y,z> = QQ[] sage: f1 = -2 * x^2 + x^3 sage: f2 = -2 * y + x* y sage: f3 = -x^2 + y^2 sage: F = Ideal([f1,f2,f3]) sage: g = x*y - 3*x*y^2 sage: g.reduce(F) -6*y^2 + 2*y sage: g.reduce(F.gens()) -6*y^2 + 2*y
\(\ZZ\) is also supported.
sage: P.<x,y,z> = ZZ[] sage: f1 = -2 * x^2 + x^3 sage: f2 = -2 * y + x* y sage: f3 = -x^2 + y^2 sage: F = Ideal([f1,f2,f3]) sage: g = x*y - 3*x*y^2 sage: g.reduce(F) -6*y^2 + 2*y sage: g.reduce(F.gens()) -6*y^2 + 2*y sage: f = 3*x sage: f.reduce([2*x,y]) 3*x
The reduction is not canonical when
I
is not a Groebner basis:sage: A.<x,y> = QQ[] sage: (x+y).reduce([x+y, x-y]) 2*y sage: (x+y).reduce([x-y, x+y]) 0
- resultant(other, variable=None)#
Compute the resultant of this polynomial and the first argument with respect to the variable given as the second argument.
If a second argument is not provide the first variable of the parent is chosen.
INPUT:
other
- polynomialvariable
- optional variable (default:None
)
EXAMPLES:
sage: P.<x,y> = PolynomialRing(QQ,2) sage: a = x+y sage: b = x^3-y^3 sage: c = a.resultant(b); c -2*y^3 sage: d = a.resultant(b,y); d 2*x^3
The SINGULAR example:
sage: R.<x,y,z> = PolynomialRing(GF(32003),3) sage: f = 3 * (x+2)^3 + y sage: g = x+y+z sage: f.resultant(g,x) 3*y^3 + 9*y^2*z + 9*y*z^2 + 3*z^3 - 18*y^2 - 36*y*z - 18*z^2 + 35*y + 36*z - 24
Resultants are also supported over the Integers:
sage: R.<x,y,a,b,u>=PolynomialRing(ZZ, 5, order='lex') sage: r = (x^4*y^2+x^2*y-y).resultant(x*y-y*a-x*b+a*b+u,x) sage: r y^6*a^4 - 4*y^5*a^4*b - 4*y^5*a^3*u + y^5*a^2 - y^5 + 6*y^4*a^4*b^2 + 12*y^4*a^3*b*u - 4*y^4*a^2*b + 6*y^4*a^2*u^2 - 2*y^4*a*u + 4*y^4*b - 4*y^3*a^4*b^3 - 12*y^3*a^3*b^2*u + 6*y^3*a^2*b^2 - 12*y^3*a^2*b*u^2 + 6*y^3*a*b*u - 4*y^3*a*u^3 - 6*y^3*b^2 + y^3*u^2 + y^2*a^4*b^4 + 4*y^2*a^3*b^3*u - 4*y^2*a^2*b^3 + 6*y^2*a^2*b^2*u^2 - 6*y^2*a*b^2*u + 4*y^2*a*b*u^3 + 4*y^2*b^3 - 2*y^2*b*u^2 + y^2*u^4 + y*a^2*b^4 + 2*y*a*b^3*u - y*b^4 + y*b^2*u^2
- sub_m_mul_q(m, q)#
Return
self - m*q
, wherem
must be a monomial andq
a polynomial.INPUT:
m
- a monomialq
- a polynomial
EXAMPLES:
sage: P.<x,y,z>=PolynomialRing(QQ,3) sage: x.sub_m_mul_q(y,z) -y*z + x
- subs(fixed=None, **kw)#
Fixes some given variables in a given multivariate polynomial and returns the changed multivariate polynomials. The polynomial itself is not affected. The variable,value pairs for fixing are to be provided as dictionary of the form
{variable:value}
.This is a special case of evaluating the polynomial with some of the variables constants and the others the original variables, but should be much faster if only few variables are to be fixed.
INPUT:
fixed
- (optional) dict with variable:value pairs**kw
- names parameters
- OUTPUT:
a new multivariate polynomial
EXAMPLES:
sage: R.<x,y> = QQ[] sage: f = x^2 + y + x^2*y^2 + 5 sage: f(5,y) 25*y^2 + y + 30 sage: f.subs({x:5}) 25*y^2 + y + 30 sage: f.subs(x=5) 25*y^2 + y + 30 sage: P.<x,y,z> = PolynomialRing(GF(2),3) sage: f = x + y + 1 sage: f.subs({x:y+1}) 0 sage: f.subs(x=y) 1 sage: f.subs(x=x) x + y + 1 sage: f.subs({x:z}) y + z + 1 sage: f.subs(x=z+1) y + z sage: f.subs(x=1/y) (y^2 + y + 1)/y sage: f.subs({x:1/y}) (y^2 + y + 1)/y
The parameters are substituted in order and without side effects:
sage: R.<x,y>=QQ[] sage: g=x+y sage: g.subs({x:x+1,y:x*y}) x*y + x + 1 sage: g.subs({x:x+1}).subs({y:x*y}) x*y + x + 1 sage: g.subs({y:x*y}).subs({x:x+1}) x*y + x + y + 1
sage: R.<x,y> = QQ[] sage: f = x + 2*y sage: f.subs(x=y,y=x) 2*x + y
- total_degree(std_grading=False)#
Return the total degree of
self
, which is the maximum degree of all monomials inself
.EXAMPLES:
sage: R.<x,y,z> = QQ[] sage: f = 2*x*y^3*z^2 sage: f.total_degree() 6 sage: f = 4*x^2*y^2*z^3 sage: f.total_degree() 7 sage: f = 99*x^6*y^3*z^9 sage: f.total_degree() 18 sage: f = x*y^3*z^6+3*x^2 sage: f.total_degree() 10 sage: f = z^3+8*x^4*y^5*z sage: f.total_degree() 10 sage: f = z^9+10*x^4+y^8*x^2 sage: f.total_degree() 10
A matrix term ordering changes the grading. To get the total degree using the standard grading, use
std_grading=True
:sage: tord = TermOrder(matrix(3, [3,2,1,1,1,0,1,0,0])) sage: tord Matrix term order with matrix [3 2 1] [1 1 0] [1 0 0] sage: R.<x,y,z> = PolynomialRing(QQ, order=tord) sage: f = x^2*y sage: f.total_degree() 8 sage: f.total_degree(std_grading=True) 3
- univariate_polynomial(R=None)#
Returns a univariate polynomial associated to this multivariate polynomial.
INPUT:
R
- (default:None
) PolynomialRing
If this polynomial is not in at most one variable, then a
ValueError
exception is raised. This is checked using theis_univariate()
method. The new Polynomial is over the same base ring as the givenMPolynomial
and in the variablex
if no ringR
is provided.EXAMPLES:
sage: R.<x, y> = QQ[] sage: f = 3*x^2 - 2*y + 7*x^2*y^2 + 5 sage: f.univariate_polynomial() Traceback (most recent call last): ... TypeError: polynomial must involve at most one variable sage: g = f.subs({x:10}); g 700*y^2 - 2*y + 305 sage: g.univariate_polynomial () 700*y^2 - 2*y + 305 sage: g.univariate_polynomial(PolynomialRing(QQ,'z')) 700*z^2 - 2*z + 305
Here’s an example with a constant multivariate polynomial:
sage: g = R(1) sage: h = g.univariate_polynomial(); h 1 sage: h.parent() Univariate Polynomial Ring in x over Rational Field
- variable(i=0)#
Return the i-th variable occurring in self. The index i is the index in
self.variables()
.EXAMPLES:
sage: P.<x,y,z> = GF(2)[] sage: f = x*z^2 + z + 1 sage: f.variables() (x, z) sage: f.variable(1) z
- variables()#
Return a tuple of all variables occurring in self.
EXAMPLES:
sage: P.<x,y,z> = GF(2)[] sage: f = x*z^2 + z + 1 sage: f.variables() (x, z)
- sage.rings.polynomial.multi_polynomial_libsingular.unpickle_MPolynomialRing_libsingular(base_ring, names, term_order)#
inverse function for
MPolynomialRing_libsingular.__reduce__
EXAMPLES:
sage: P.<x,y> = PolynomialRing(QQ) sage: loads(dumps(P)) is P # indirect doctest True
- sage.rings.polynomial.multi_polynomial_libsingular.unpickle_MPolynomial_libsingular(R, d)#
Deserialize an
MPolynomial_libsingular
objectINPUT:
R
- the base ringd
- a Python dictionary as returned byMPolynomial_libsingular.dict()
EXAMPLES:
sage: P.<x,y> = PolynomialRing(QQ) sage: loads(dumps(x)) == x # indirect doctest True