Base class for elements of multivariate polynomial rings
Bases: sage.structure.element.CommutativeRingElement
INPUT:
Returns the named of the arguments of self, in the order they are accepted from call.
EXAMPLES:
sage: R.<x,y> = ZZ[]
sage: x.args()
(x, y)
Return a copy of this polynomial but with coefficients in R, if at all possible.
INPUT:
EXAMPLES:
sage: R.<x,y> = QQ[]
sage: f = x^3 + 3/5*y + 1
sage: f.change_ring(GF(7))
x^3 + 2*y + 1
sage: R.<x,y> = GF(9,'a')[]
sage: (x+2*y).change_ring(GF(3))
x - y
Return the nonzero coefficients of this polynomial in a list. The returned list is decreasingly ordered by the term ordering of self.parent(), i.e. the list of coefficients matches the list of monomials returned by sage.rings.polynomial.multi_polynomial_libsingular.MPolynomial_libsingular.monomials().
EXAMPLES:
sage: R.<x,y,z> = PolynomialRing(QQ,3,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,3,order='lex')
sage: f=23*x^6*y^7 + x^3*y+6*x^7*z
sage: f.coefficients()
[6, 23, 1]
Test the same stuff with base ring – different implementation:
sage: R.<x,y,z> = PolynomialRing(ZZ,3,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(ZZ,3,order='lex')
sage: f=23*x^6*y^7 + x^3*y+6*x^7*z
sage: f.coefficients()
[6, 23, 1]
AUTHOR:
Returns the content of this polynomial. Here, we define content as the gcd of the coefficients in the base ring.
EXAMPLES:
sage: R.<x,y> = ZZ[]
sage: f = 4*x+6*y
sage: f.content()
2
sage: f.content().parent()
Integer Ring
TESTS:
Since trac ticket #10771, the gcd in QQ restricts to the gcd in ZZ:
sage: R.<x,y> = QQ[]
sage: f = 4*x+6*y
sage: f.content(); f.content().parent()
2
Rational Field
Return a denominator of self.
First, the lcm of the denominators of the entries of self is computed and returned. If this computation fails, the unit of the parent of self is returned.
Note that some subclases may implement its own denominator function.
Warning
This is not the denominator of the rational function defined by self, which would always be 1 since self is a polynomial.
EXAMPLES:
First we compute the denominator of a polynomial with integer coefficients, which is of course 1.
sage: R.<x,y> = ZZ[]
sage: f = x^3 + 17*y + x + y
sage: f.denominator()
1
Next we compute the denominator of a polynomial over a number field.
sage: R.<x,y> = NumberField(symbolic_expression(x^2+3) ,'a')['x,y']
sage: f = (1/17)*x^19 + (1/6)*y - (2/3)*x + 1/3; f
1/17*x^19 - 2/3*x + 1/6*y + 1/3
sage: f.denominator()
102
Finally, we try to compute the denominator of a polynomial with coefficients in the real numbers, which is a ring whose elements do not have a denominator method.
sage: R.<a,b,c> = RR[]
sage: f = a + b + RR('0.3'); f
a + b + 0.300000000000000
sage: f.denominator()
1.00000000000000
Check that the denominator is an element over the base whenever the base has no denominator function. This closes #9063
sage: R.<a,b,c> = GF(5)[]
sage: x = R(0)
sage: x.denominator()
1
sage: type(x.denominator())
<type 'sage.rings.finite_rings.integer_mod.IntegerMod_int'>
sage: type(a.denominator())
<type 'sage.rings.finite_rings.integer_mod.IntegerMod_int'>
sage: from sage.rings.polynomial.multi_polynomial_element import MPolynomial
sage: isinstance(a / b, MPolynomial)
False
sage: isinstance(a.numerator() / a.denominator(), MPolynomial)
True
The formal derivative of this polynomial, with respect to variables supplied in args.
Multiple variables and iteration counts may be supplied; see documentation for the global derivative() function for more details.
See also
_derivative()
EXAMPLES:
Polynomials implemented via Singular:
sage: R.<x, y> = PolynomialRing(FiniteField(5))
sage: f = x^3*y^5 + x^7*y
sage: type(f)
<type 'sage.rings.polynomial.multi_polynomial_libsingular.MPolynomial_libsingular'>
sage: f.derivative(x)
2*x^6*y - 2*x^2*y^5
sage: f.derivative(y)
x^7
Generic multivariate polynomials:
sage: R.<t> = PowerSeriesRing(QQ)
sage: S.<x, y> = PolynomialRing(R)
sage: f = (t^2 + O(t^3))*x^2*y^3 + (37*t^4 + O(t^5))*x^3
sage: type(f)
<class 'sage.rings.polynomial.multi_polynomial_element.MPolynomial_polydict'>
sage: f.derivative(x) # with respect to x
(2*t^2 + O(t^3))*x*y^3 + (111*t^4 + O(t^5))*x^2
sage: f.derivative(y) # with respect to y
(3*t^2 + O(t^3))*x^2*y^2
sage: f.derivative(t) # with respect to t (recurses into base ring)
(2*t + O(t^2))*x^2*y^3 + (148*t^3 + O(t^4))*x^3
sage: f.derivative(x, y) # with respect to x and then y
(6*t^2 + O(t^3))*x*y^2
sage: f.derivative(y, 3) # with respect to y three times
(6*t^2 + O(t^3))*x^2
sage: f.derivative() # can't figure out the variable
Traceback (most recent call last):
...
ValueError: must specify which variable to differentiate with respect to
Polynomials over the symbolic ring (just for fun....):
sage: x = var("x")
sage: S.<u, v> = PolynomialRing(SR)
sage: f = u*v*x
sage: f.derivative(x) == u*v
True
sage: f.derivative(u) == v*x
True
Return a list of partial derivatives of this polynomial, ordered by the variables of self.parent().
EXAMPLES:
sage: P.<x,y,z> = PolynomialRing(ZZ,3)
sage: f = x*y + 1
sage: f.gradient()
[y, x, 0]
Return the homogenization of this polynomial.
The polynomial itself is returned if it is homogeneous already. Otherwise, the monomials are multiplied with the smallest powers of var such that they all have the same total degree.
INPUT:
OUTPUT:
If var specifies a variable in the polynomial ring, then a homogeneous element in that ring is returned. Otherwise, a homogeneous element is returned in a polynomial ring with an extra last variable var.
EXAMPLES:
sage: R.<x,y> = QQ[]
sage: f = x^2 + y + 1 + 5*x*y^10
sage: f.homogenize()
5*x*y^10 + x^2*h^9 + y*h^10 + h^11
The parameter var can be used to specify the name of the variable:
sage: g = f.homogenize('z'); g
5*x*y^10 + x^2*z^9 + y*z^10 + z^11
sage: g.parent()
Multivariate Polynomial Ring in x, y, z over Rational Field
However, if the polynomial is homogeneous already, then that parameter is ignored and no extra variable is added to the polynomial ring:
sage: f = x^2 + y^2
sage: g = f.homogenize('z'); g
x^2 + y^2
sage: g.parent()
Multivariate Polynomial Ring in x, y over Rational Field
If you want the ring of the result to be independent of whether the polynomial is homogenized, you can use var to use an existing variable to homogenize:
sage: R.<x,y,z> = QQ[]
sage: f = x^2 + y^2
sage: g = f.homogenize(z); g
x^2 + y^2
sage: g.parent()
Multivariate Polynomial Ring in x, y, z over Rational Field
sage: f = x^2 - y
sage: g = f.homogenize(z); g
x^2 - y*z
sage: g.parent()
Multivariate Polynomial Ring in x, y, z over Rational Field
The parameter var can also be given as a zero-based index in the list of variables:
sage: g = f.homogenize(2); g
x^2 - y*z
If the variable specified by var is not present in the polynomial, then setting it to 1 yields the original polynomial:
sage: g(x,y,1)
x^2 - y
If it is present already, this might not be the case:
sage: g = f.homogenize(x); g
x^2 - x*y
sage: g(1,y,z)
-y + 1
In particular, this can be surprising in positive characteristic:
sage: R.<x,y> = GF(2)[]
sage: f = x + 1
sage: f.homogenize(x)
0
TESTS:
sage: R = PolynomialRing(QQ, 'x', 5)
sage: p = R.random_element()
sage: q1 = p.homogenize()
sage: q2 = p.homogenize()
sage: q1.parent() is q2.parent()
True
Returns an inverse of self modulo the polynomial ideal ,
namely a multivariate polynomial
such that
self * f - 1 belongs to
.
OUTPUT:
- a multivariate polynomial representing the inverse of f modulo I
EXAMPLES:
sage: R.<x1,x2> = QQ[]
sage: I = R.ideal(x2**2 + x1 - 2, x1**2 - 1)
sage: f = x1 + 3*x2^2; g = f.inverse_mod(I); g
1/16*x1 + 3/16
sage: (f*g).reduce(I)
1
Test a non-invertible element:
sage: R.<x1,x2> = QQ[]
sage: I = R.ideal(x2**2 + x1 - 2, x1**2 - 1)
sage: f = x1 + x2
sage: f.inverse_mod(I)
Traceback (most recent call last):
...
ArithmeticError: element is non-invertible
Returns True if this polynomial is a generator of its parent.
EXAMPLES:
sage: R.<x,y>=ZZ[]
sage: x.is_generator()
True
sage: (x+y-y).is_generator()
True
sage: (x*y).is_generator()
False
sage: R.<x,y>=QQ[]
sage: x.is_generator()
True
sage: (x+y-y).is_generator()
True
sage: (x*y).is_generator()
False
Return True if self is a homogeneous polynomial.
TESTS:
sage: from sage.rings.polynomial.multi_polynomial import MPolynomial
sage: P.<x, y> = PolynomialRing(QQ, 2)
sage: MPolynomial.is_homogeneous(x+y)
True
sage: MPolynomial.is_homogeneous(P(0))
True
sage: MPolynomial.is_homogeneous(x+y^2)
False
sage: MPolynomial.is_homogeneous(x^2 + y^2)
True
sage: MPolynomial.is_homogeneous(x^2 + y^2*x)
False
sage: MPolynomial.is_homogeneous(x^2*y + y^2*x)
True
Note
This is a generic implementation which is likely overridden by subclasses.
Return the Jacobian ideal of the polynomial self.
EXAMPLES:
sage: R.<x,y,z> = QQ[]
sage: f = x^3 + y^3 + z^3
sage: f.jacobian_ideal()
Ideal (3*x^2, 3*y^2, 3*z^2) of Multivariate Polynomial Ring in x, y, z over Rational Field
given an ideal I = (f_1,...,f_r) and some g (== self) in I, find s_1,...,s_r such that g = s_1 f_1 + ... + s_r f_r.
EXAMPLE:
sage: A.<x,y> = PolynomialRing(CC,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
This is an implementation of the Macaulay Resultant. It computes the resultant of universal polynomials as well as polynomials with constant coefficients. This is a project done in sage days 55. It’s based on the implementation in Maple by Manfred Minimair, which in turn is based on the references [CLO], [Can], [Mac]. It calculates the Macaulay resultant for a list of Polynomials, up to sign!
AUTHORS:
INPUT:
works when args[0] is the list of polynomials, or args is itself the list of polynomials
OUTPUT:
EXAMPLES:
The number of polynomials has to match the number of variables:
sage: R.<x,y,z> = PolynomialRing(QQ,3)
sage: y.macaulay_resultant(x+z)
Traceback (most recent call last):
...
TypeError: number of polynomials(= 2) must equal number of variables (= 3)
The polynomials need to be all homogeneous:
sage: R.<x,y,z> = PolynomialRing(QQ,3)
sage: y.macaulay_resultant([x+z, z+x^3])
Traceback (most recent call last):
...
TypeError: resultant for non-homogeneous polynomials is not supported
All polynomials must be in the same ring:
sage: R.<x,y,z> = PolynomialRing(QQ,3)
sage: S.<x,y> = PolynomialRing(QQ, 2)
sage: y.macaulay_resultant(z+x,z)
Traceback (most recent call last):
...
TypeError: not all inputs are polynomials in the calling ring
The following example recreates Proposition 2.10 in Ch.3 of Using Algebraic Geometry:
sage: K.<x,y> = PolynomialRing(ZZ, 2)
sage: flist,R = K._macaulay_resultant_universal_polynomials([1,1,2])
sage: flist[0].macaulay_resultant(flist[1:])
u2^2*u4^2*u6 - 2*u1*u2*u4*u5*u6 + u1^2*u5^2*u6 - u2^2*u3*u4*u7 + u1*u2*u3*u5*u7 + u0*u2*u4*u5*u7 - u0*u1*u5^2*u7 + u1*u2*u3*u4*u8 - u0*u2*u4^2*u8 - u1^2*u3*u5*u8 + u0*u1*u4*u5*u8 + u2^2*u3^2*u9 - 2*u0*u2*u3*u5*u9 + u0^2*u5^2*u9 - u1*u2*u3^2*u10 + u0*u2*u3*u4*u10 + u0*u1*u3*u5*u10 - u0^2*u4*u5*u10 + u1^2*u3^2*u11 - 2*u0*u1*u3*u4*u11 + u0^2*u4^2*u11
The following example degenerates into the determinant of a matrix:
sage: K.<x,y> = PolynomialRing(ZZ, 2)
sage: flist,R = K._macaulay_resultant_universal_polynomials([1,1,1])
sage: flist[0].macaulay_resultant(flist[1:])
-u2*u4*u6 + u1*u5*u6 + u2*u3*u7 - u0*u5*u7 - u1*u3*u8 + u0*u4*u8
The following example is by Patrick Ingram(arxiv:1310.4114):
sage: U = PolynomialRing(ZZ,'y',2); y0,y1 = U.gens()
sage: R = PolynomialRing(U,'x',3); x0,x1,x2 = R.gens()
sage: f0 = y0*x2^2 - x0^2 + 2*x1*x2
sage: f1 = y1*x2^2 - x1^2 + 2*x0*x2
sage: f2 = x0*x1 - x2^2
sage: f0.macaulay_resultant(f1,f2)
y0^2*y1^2 - 4*y0^3 - 4*y1^3 + 18*y0*y1 - 27
a simple example with constant rational coefficients:
sage: R.<x,y,z,w> = PolynomialRing(QQ,4)
sage: w.macaulay_resultant([z,y,x])
1
an example where the resultant vanishes:
sage: R.<x,y,z> = PolynomialRing(QQ,3)
sage: (x+y).macaulay_resultant([y^2,x])
0
an example of bad reduction at a prime p = 5:
sage: R.<x,y,z> = PolynomialRing(QQ,3)
sage: y.macaulay_resultant([x^3+25*y^2*x,5*z])
125
The input can given as an unpacked list of polynomials:
sage: R.<x,y,z> = PolynomialRing(QQ,3)
sage: y.macaulay_resultant(x^3+25*y^2*x,5*z)
125
an example when the coefficients live in a finite field:
sage: F = FiniteField(11)
sage: R.<x,y,z,w> = PolynomialRing(F,4)
sage: z.macaulay_resultant([x^3,5*y,w])
4
example when the denominator in the algorithm vanishes(in this case the resultant is the constant term of the quotient of char polynomials of numerator/denominator):
sage: R.<x,y,z> = PolynomialRing(QQ,3)
sage: y.macaulay_resultant([x+z, z^2])
-1
when there are only 2 polynomials, macaulay resultant degenerates to the traditional resultant:
sage: R.<x> = PolynomialRing(QQ,1)
sage: f = x^2+1; g = x^5+1
sage: fh = f.homogenize()
sage: gh = g.homogenize()
sage: RH = fh.parent()
sage: f.resultant(g) == fh.macaulay_resultant(gh)
True
Returns the polynomial obtained by applying f to the non-zero coefficients of self.
If f is a sage.categories.map.Map, then the resulting polynomial will be defined over the codomain of f. Otherwise, the resulting polynomial will be over the same ring as self. Set new_base_ring to override this behaviour.
INPUT:
EXAMPLES:
sage: k.<a> = GF(9); R.<x,y> = k[]; f = x*a + 2*x^3*y*a + a
sage: f.map_coefficients(lambda a : a + 1)
(-a + 1)*x^3*y + (a + 1)*x + (a + 1)
Examples with different base ring:
sage: R.<r> = GF(9); S.<s> = GF(81)
sage: h = Hom(R,S)[0]; h
Ring morphism:
From: Finite Field in r of size 3^2
To: Finite Field in s of size 3^4
Defn: r |--> 2*s^3 + 2*s^2 + 1
sage: T.<X,Y> = R[]
sage: f = r*X+Y
sage: g = f.map_coefficients(h); g
(-s^3 - s^2 + 1)*X + Y
sage: g.parent()
Multivariate Polynomial Ring in X, Y over Finite Field in s of size 3^4
sage: h = lambda x: x.trace()
sage: g = f.map_coefficients(h); g
X - Y
sage: g.parent()
Multivariate Polynomial Ring in X, Y over Finite Field in r of size 3^2
sage: g = f.map_coefficients(h, new_base_ring=GF(3)); g
X - Y
sage: g.parent()
Multivariate Polynomial Ring in X, Y over Finite Field of size 3
Return the Newton polytope of this polynomial.
EXAMPLES:
sage: R.<x,y> = QQ[]
sage: f = 1 + x*y + x^3 + y^3
sage: P = f.newton_polytope()
sage: P
A 2-dimensional polyhedron in ZZ^2 defined as the convex hull of 3 vertices
sage: P.is_simple()
True
TESTS:
sage: R.<x,y> = QQ[]
sage: R(0).newton_polytope()
The empty polyhedron in ZZ^0
sage: R(1).newton_polytope()
A 0-dimensional polyhedron in ZZ^2 defined as the convex hull of 1 vertex
sage: R(x^2+y^2).newton_polytope().integral_points()
((0, 2), (1, 1), (2, 0))
Return a numerator of self computed as self * self.denominator()
Note that some subclases may implement its own numerator function.
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*x + y + 1
sage: f.numerator()
x^3 + 17*x + y + 1
sage: f == f.numerator()
True
Next we compute the numerator of a polynomial over a number field.
sage: R.<x,y> = NumberField(symbolic_expression(x^2+3) ,'a')['x,y']
sage: f = (1/17)*y^19 - (2/3)*x + 1/3; f
1/17*y^19 - 2/3*x + 1/3
sage: f.numerator()
3*y^19 - 34*x + 17
sage: f == f.numerator()
False
We try to compute the numerator of a polynomial with coefficients in the finite field of 3 elements.
sage: K.<x,y,z> = GF(3)['x, y, z']
sage: f = 2*x*z + 2*z^2 + 2*y + 1; f
-x*z - z^2 - y + 1
sage: f.numerator()
-x*z - z^2 - y + 1
We check that the computation the numerator and denominator are valid
sage: K=NumberField(symbolic_expression('x^3+2'),'a')['x']['s,t']
sage: f=K.random_element()
sage: f.numerator() / f.denominator() == f
True
sage: R=RR['x,y,z']
sage: f=R.random_element()
sage: f.numerator() / f.denominator() == f
True
Let var be one of the variables of the parent of self. This returns self viewed as a univariate polynomial in var over the polynomial ring generated by all the other variables of the parent.
EXAMPLES:
sage: R.<x,w,z> = QQ[]
sage: f = x^3 + 3*w*x + w^5 + (17*w^3)*x + z^5
sage: f.polynomial(x)
x^3 + (17*w^3 + 3*w)*x + w^5 + z^5
sage: parent(f.polynomial(x))
Univariate Polynomial Ring in x over Multivariate Polynomial Ring in w, z over Rational Field
sage: f.polynomial(w)
w^5 + 17*x*w^3 + 3*x*w + z^5 + x^3
sage: f.polynomial(z)
z^5 + w^5 + 17*x*w^3 + x^3 + 3*x*w
sage: R.<x,w,z,k> = ZZ[]
sage: f = x^3 + 3*w*x + w^5 + (17*w^3)*x + z^5 +x*w*z*k + 5
sage: f.polynomial(x)
x^3 + (17*w^3 + w*z*k + 3*w)*x + w^5 + z^5 + 5
sage: f.polynomial(w)
w^5 + 17*x*w^3 + (x*z*k + 3*x)*w + z^5 + x^3 + 5
sage: f.polynomial(z)
z^5 + x*w*k*z + w^5 + 17*x*w^3 + x^3 + 3*x*w + 5
sage: f.polynomial(k)
x*w*z*k + w^5 + z^5 + 17*x*w^3 + x^3 + 3*x*w + 5
sage: R.<x,y>=GF(5)[]
sage: f=x^2+x+y
sage: f.polynomial(x)
x^2 + x + y
sage: f.polynomial(y)
y + x^2 + x
Given two nonzero polynomials self and right, returns the Sylvester matrix of the polynomials with respect to a given variable.
Note that the Sylvester matrix is not defined if one of the polynomials is zero.
INPUT:
OUTPUT:
EXAMPLES:
sage: R.<x, y> = PolynomialRing(ZZ)
sage: f = (y + 1)*x + 3*x**2
sage: g = (y + 2)*x + 4*x**2
sage: M = f.sylvester_matrix(g, x)
sage: print M
[ 3 y + 1 0 0]
[ 0 3 y + 1 0]
[ 4 y + 2 0 0]
[ 0 4 y + 2 0]
If the polynomials share a non-constant common factor then the determinant of the Sylvester matrix will be zero:
sage: M.determinant()
0
sage: f.sylvester_matrix(1 + g, x).determinant()
y^2 - y + 7
If both polynomials are of positive degree with respect to variable, the determinant of the Sylvester matrix is the resultant:
sage: f = R.random_element(4)
sage: g = R.random_element(4)
sage: f.sylvester_matrix(g, x).determinant() == f.resultant(g, x)
True
TEST:
The variable is optional:
sage: f = x + y
sage: g = x + y
sage: f.sylvester_matrix(g)
[1 y]
[1 y]
Polynomials must be defined over compatible base rings:
sage: K.<x, y> = QQ[]
sage: f = x + y
sage: L.<x, y> = ZZ[]
sage: g = x + y
sage: R.<x, y> = GF(25, 'a')[]
sage: h = x + y
sage: f.sylvester_matrix(g, 'x')
[1 y]
[1 y]
sage: g.sylvester_matrix(h, 'x')
[1 y]
[1 y]
sage: f.sylvester_matrix(h, 'x')
Traceback (most recent call last):
...
TypeError: no common canonical parent for objects with parents: 'Multivariate Polynomial Ring in x, y over Rational Field' and 'Multivariate Polynomial Ring in x, y over Finite Field in a of size 5^2'
sage: K.<x, y, z> = QQ[]
sage: f = x + y
sage: L.<x, z> = QQ[]
sage: g = x + z
sage: f.sylvester_matrix(g)
[1 y]
[1 z]
Corner cases:
sage: K.<x ,y>=QQ[]
sage: f = x^2+1
sage: g = K(0)
sage: f.sylvester_matrix(g)
Traceback (most recent call last):
...
ValueError: The Sylvester matrix is not defined for zero polynomials
sage: g.sylvester_matrix(f)
Traceback (most recent call last):
...
ValueError: The Sylvester matrix is not defined for zero polynomials
sage: g.sylvester_matrix(g)
Traceback (most recent call last):
...
ValueError: The Sylvester matrix is not defined for zero polynomials
sage: K(3).sylvester_matrix(x^2)
[3 0]
[0 3]
sage: K(3).sylvester_matrix(K(4))
[]
Returns a new multivariate polynomial obtained from self by deleting all terms that involve the given variable to a power at least n.
Return the weighted degree of self, which is the maximum weighted degree of all monomials in self; the weighted degree of a monomial is the sum of all powers of the variables in the monomial, each power multiplied with its respective weight in weights.
This method is given for convenience. It is faster to use polynomial rings with weighted term orders and the standard degree function.
INPUT:
EXAMPLES:
sage: R.<x,y,z> = GF(7)[]
sage: p = x^3 + y + x*z^2
sage: p.weighted_degree({z:0, x:1, y:2})
3
sage: p.weighted_degree(1, 2, 0)
3
sage: p.weighted_degree((1, 4, 2))
5
sage: p.weighted_degree((1, 4, 1))
4
sage: p.weighted_degree(2**64, 2**50, 2**128)
680564733841876926945195958937245974528
sage: q = R.random_element(100, 20) #random
sage: q.weighted_degree(1, 1, 1) == q.total_degree()
True
You may also work with negative weights
sage: p.weighted_degree(-1, -2, -1)
-2
Note that only integer weights are allowed
sage: p.weighted_degree(x,1,1)
Traceback (most recent call last):
...
TypeError
sage: p.weighted_degree(2/1,1,1)
6
The weighted_degree coincides with the degree of a weighted polynomial ring, but the later is faster.
sage: K = PolynomialRing(QQ, 'x,y', order=TermOrder('wdegrevlex', (2,3)))
sage: p = K.random_element(10)
sage: p.degree() == p.weighted_degree(2,3)
True
TESTS:
sage: R = PolynomialRing(QQ, 'a', 5)
sage: f = R.random_element(terms=20)
sage: w = random_vector(ZZ,5)
sage: d1 = f.weighted_degree(w)
sage: d2 = (f*1.0).weighted_degree(w)
sage: d1 == d2
True