Base class for elements of multivariate polynomial rings¶
- class sage.rings.polynomial.multi_polynomial.MPolynomial¶
Bases:
sage.structure.element.CommutativeRingElement
- args()¶
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)
- change_ring(R)¶
Return a copy of this polynomial but with coefficients in
R
, if at all possible.INPUT:
R
– a ring or morphism.
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
sage: K.<z> = CyclotomicField(3) sage: R.<x,y> = K[] sage: f = x^2 + z*y sage: f.change_ring(K.embeddings(CC)[1]) x^2 + (-0.500000000000000 - 0.866025403784438*I)*y
- coefficients()¶
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 bysage.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 \(\ZZ\) – 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:
Didier Deshommes
- content()¶
Returns the content of this polynomial. Here, we define content as the gcd of the coefficients in the base ring.
See also
EXAMPLES:
sage: R.<x,y> = ZZ[] sage: f = 4*x+6*y sage: f.content() 2 sage: f.content().parent() Integer Ring
- content_ideal()¶
Return the content ideal of this polynomial, defined as the ideal generated by its coefficients.
See also
EXAMPLES:
sage: R.<x,y> = ZZ[] sage: f = 2*x*y + 6*x - 4*y + 2 sage: f.content_ideal() Principal ideal (2) of Integer Ring sage: S.<z,t> = R[] sage: g = x*z + y*t sage: g.content_ideal() Ideal (x, y) of Multivariate Polynomial Ring in x, y over Integer Ring
- denominator()¶
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 subclasses 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 trac ticket #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
- derivative(*args)¶
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
- discriminant(variable)¶
Returns the discriminant of self with respect to the given variable.
INPUT:
variable
- The variable with respect to which we computethe discriminant
OUTPUT:
An element of the base ring of the polynomial ring.
EXAMPLES:
sage: R.<x,y,z>=QQ[] sage: f=4*x*y^2 + 1/4*x*y*z + 3/2*x*z^2 - 1/2*z^2 sage: f.discriminant(x) 1 sage: f.discriminant(y) -383/16*x^2*z^2 + 8*x*z^2 sage: f.discriminant(z) -383/16*x^2*y^2 + 8*x*y^2
Note that, unlike the univariate case, the result lives in the same ring as the polynomial:
sage: R.<x,y>=QQ[] sage: f=x^5*y+3*x^2*y^2-2*x+y-1 sage: f.discriminant(y) x^10 + 2*x^5 + 24*x^3 + 12*x^2 + 1 sage: f.polynomial(y).discriminant() x^10 + 2*x^5 + 24*x^3 + 12*x^2 + 1 sage: f.discriminant(y).parent()==f.polynomial(y).discriminant().parent() False
- AUTHOR:
Miguel Marco
- gcd(other)¶
Return a greatest common divisor of this polynomial and
other
.INPUT:
other
– a polynomial with the same parent as this polynomial
EXAMPLES:
sage: Q.<z> = Frac(QQ['z']) sage: R.<x,y> = Q[] sage: r = x*y - (2*z-1)/(z^2+z+1) * x + y/z sage: p = r * (x + z*y - 1/z^2) sage: q = r * (x*y*z + 1) sage: gcd(p,q) (z^3 + z^2 + z)*x*y + (-2*z^2 + z)*x + (z^2 + z + 1)*y
Polynomials over polynomial rings are converted to a simpler polynomial ring with all variables to compute the gcd:
sage: A.<z,t> = ZZ[] sage: B.<x,y> = A[] sage: r = x*y*z*t+1 sage: p = r * (x - y + z - t + 1) sage: q = r * (x*z - y*t) sage: gcd(p,q) z*t*x*y + 1 sage: _.parent() Multivariate Polynomial Ring in x, y over Multivariate Polynomial Ring in z, t over Integer Ring
Some multivariate polynomial rings have no gcd implementation:
sage: R.<x,y> =GaussianIntegers()[] sage: x.gcd(x) Traceback (most recent call last): ... NotImplementedError: GCD is not implemented for multivariate polynomials over Gaussian Integers in Number Field in I with defining polynomial x^2 + 1 with I = 1*I
- gradient()¶
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]
- homogeneous_components()¶
Return the homogeneous components of this polynomial.
OUTPUT:
A dictionary mapping degrees to homogeneous polynomials.
EXAMPLES:
sage: R.<x,y> = QQ[] sage: (x^3 + 2*x*y^3 + 4*y^3 + y).homogeneous_components() {1: y, 3: x^3 + 4*y^3, 4: 2*x*y^3} sage: R.zero().homogeneous_components() {}
In case of weighted term orders, the polynomials are homogeneous with respect to the weights:
sage: S.<a,b,c> = PolynomialRing(ZZ, order=TermOrder('wdegrevlex', (1,2,3))) sage: (a^6 + b^3 + b*c + a^2*c + c + a + 1).homogeneous_components() {0: 1, 1: a, 3: c, 5: a^2*c + b*c, 6: a^6 + b^3}
- homogenize(var='h')¶
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:
var
– a variable in the polynomial ring (as a string, an element of the ring, or a zero-based index in the list of variables) or a name for a new variable (default:'h'
)
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 variablevar
.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
- inverse_mod(I)¶
Returns an inverse of self modulo the polynomial ideal \(I\), namely a multivariate polynomial \(f\) such that
self * f - 1
belongs to \(I\).- INPUT:
I
– an ideal of the polynomial ring in which self lives
OUTPUT:
a multivariate polynomial representing the inverse of
f
moduloI
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
- is_generator()¶
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
- is_homogeneous()¶
Return
True
if self is a homogeneous polynomial.Note
This is a generic implementation which is likely overridden by subclasses.
- is_nilpotent()¶
Return
True
ifself
is nilpotent, i.e., some power ofself
is 0.EXAMPLES:
sage: R.<x,y> = QQbar[] sage: (x+y).is_nilpotent() False sage: R(0).is_nilpotent() True sage: _.<x,y> = Zmod(4)[] sage: (2*x).is_nilpotent() True sage: (2+y*x).is_nilpotent() False sage: _.<x,y> = Zmod(36)[] sage: (4+6*x).is_nilpotent() False sage: (6*x + 12*y + 18*x*y + 24*(x^2+y^2)).is_nilpotent() True
- is_square(root=False)¶
Test whether this polynomial is a square root.
INPUT:
root
- if set toTrue
return a pair(True, root)
whereroot
is a square root or(False, None)
if it is not a square.
EXAMPLES:
sage: R.<a,b> = QQ[] sage: a.is_square() False sage: ((1+a*b^2)^2).is_square() True sage: ((1+a*b^2)^2).is_square(root=True) (True, a*b^2 + 1)
- is_symmetric(group=None)¶
Return whether this polynomial is symmetric.
INPUT:
group
(default: symmetric group) – if set, test whether the polynomial is invariant with respect to the given permutation group
EXAMPLES:
sage: R.<x,y,z> = QQ[] sage: p = (x+y+z)**2 - 3 * (x+y)*(x+z)*(y+z) sage: p.is_symmetric() True sage: (x + y - z).is_symmetric() False sage: R.one().is_symmetric() True sage: p = (x-y)*(y-z)*(z-x) sage: p.is_symmetric() False sage: p.is_symmetric(AlternatingGroup(3)) True sage: R.<x,y> = QQ[] sage: ((x + y)**2).is_symmetric() True sage: R.one().is_symmetric() True sage: (x + 2*y).is_symmetric() False
An example with a GAP permutation group (here the quaternions):
sage: R = PolynomialRing(QQ, 'x', 8) sage: x = R.gens() sage: p = sum(prod(x[i] for i in e) for e in [(0,1,2), (0,1,7), (0,2,7), (1,2,7), (3,4,5), (3,4,6), (3,5,6), (4,5,6)]) sage: p.is_symmetric(libgap.TransitiveGroup(8, 5)) True sage: p = sum(prod(x[i] for i in e) for e in [(0,1,2), (0,1,7), (0,2,7), (1,2,7), (3,4,5), (3,4,6), (3,5,6)]) sage: p.is_symmetric(libgap.TransitiveGroup(8, 5)) False
- is_unit()¶
Return
True
ifself
is a unit, that is, has a multiplicative inverse.EXAMPLES:
sage: R.<x,y> = QQbar[] sage: (x+y).is_unit() False sage: R(0).is_unit() False sage: R(-1).is_unit() True sage: R(-1 + x).is_unit() False sage: R(2).is_unit() True
Check that trac ticket #22454 is fixed:
sage: _.<x,y> = Zmod(4)[] sage: (1 + 2*x).is_unit() True sage: (x*y).is_unit() False sage: _.<x,y> = Zmod(36)[] sage: (7+ 6*x + 12*y - 18*x*y).is_unit() 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)]
- jacobian_ideal()¶
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
- 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
.EXAMPLES:
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
- macaulay_resultant(*args)¶
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:
Hao Chen, Solomon Vishkautsan (7-2014)
INPUT:
args
– a list of \(n-1\) homogeneous polynomials in \(n\) variables.works when
args[0]
is the list of polynomials, orargs
is itself the list of polynomials
OUTPUT:
the macaulay resultant
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 \(3*3\) 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
- map_coefficients(f, new_base_ring=None)¶
Returns the polynomial obtained by applying
f
to the non-zero coefficients of self.If
f
is asage.categories.map.Map
, then the resulting polynomial will be defined over the codomain off
. Otherwise, the resulting polynomial will be over the same ring as self. Setnew_base_ring
to override this behaviour.INPUT:
f
– a callable that will be applied to the coefficients of self.new_base_ring
(optional) – if given, the resulting polynomial will be defined over this ring.
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
- newton_polytope()¶
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
- nth_root(n)¶
Return a \(n\)-th root of this element.
If there is no such root, a
ValueError
is raised.EXAMPLES:
sage: R.<x,y,z> = QQ[] sage: a = 32 * (x*y + 1)^5 * (x+y+z)^5 sage: a.nth_root(5) 2*x^2*y + 2*x*y^2 + 2*x*y*z + 2*x + 2*y + 2*z sage: b = x + 2*y + 3*z sage: b.nth_root(42) Traceback (most recent call last): ... ValueError: not a 42nd power sage: R.<x,y> = QQ[] sage: S.<z,t> = R[] sage: T.<u,v> = S[] sage: p = (1 + x*u + y + v) * (1 + z*t) sage: (p**3).nth_root(3) (x*z*t + x)*u + (z*t + 1)*v + (y + 1)*z*t + y + 1 sage: (p**3).nth_root(3).parent() is p.parent() True sage: ((1+x+z+t)**2).nth_root(3) Traceback (most recent call last): ... ValueError: not a 3rd power
- numerator()¶
Return a numerator of self computed as self * self.denominator()
Note that some subclasses 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
- polynomial(var)¶
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
- reduced_form(**kwds)¶
Return a reduced form of this polynomial.
The algorithm is from Stoll and Cremona’s “On the Reduction Theory of Binary Forms” [CS2003]. This takes a two variable homogeneous polynomial and finds a reduced form. This is a \(SL(2,\ZZ)\)-equivalent binary form whose covariant in the upper half plane is in the fundamental domain. If the polynomial has multiple roots, they are removed and the algorithm is applied to the portion without multiple roots.
This reduction should also minimize the sum of the squares of the coefficients, but this is not always the case. By default the coefficient minimizing algorithm in [HS2018] is applied. The coefficients can be minimized either with respect to the sum of their squares or the maximum of their global heights.
A portion of the algorithm uses Newton’s method to find a solution to a system of equations. If Newton’s method fails to converge to a point in the upper half plane, the function will use the less precise \(z_0\) covariant from the \(Q_0\) form as defined on page 7 of [CS2003]. Additionally, if this polynomial has a root with multiplicity at least half the total degree of the polynomial, then we must also use the \(z_0\) covariant. See [CS2003] for details.
Note that, if the covariant is within
error_limit
of the boundary but outside the fundamental domain, our function will erroneously move it to within the fundamental domain, hence our conjugation will be off by 1. If you don’t want this to happen, decrease yourerror_limit
and increase your precision.Implemented by Rebecca Lauren Miller as part of GSOC 2016. Smallest coefficients added by Ben Hutz July 2018.
INPUT:
keywords:
prec
– integer, sets the precision (default:300)return_conjugation
– boolean. Returns element of \(SL(2, \ZZ)\) (default:True)error_limit
– sets the error tolerance (default:0.000001)smallest_coeffs
– (default: True), boolean, whether to find the model with smallest coefficientsnorm_type
– either'norm'
or'height'
. What type of norm to use for smallest coefficientsemb
– (optional) embedding of based field into CC
OUTPUT:
a polynomial (reduced binary form)
a matrix (element of \(SL(2, \ZZ)\))
- TODO: When Newton’s Method doesn’t converge to a root in the upper half plane.
Now we just return z0. It would be better to modify and find the unique root in the upper half plane.
EXAMPLES:
sage: R.<x,h> = PolynomialRing(QQ) sage: f = 19*x^8 - 262*x^7*h + 1507*x^6*h^2 - 4784*x^5*h^3 + 9202*x^4*h^4\ -10962*x^3*h^5 + 7844*x^2*h^6 - 3040*x*h^7 + 475*h^8 sage: f.reduced_form(prec=200, smallest_coeffs=False) ( -x^8 - 2*x^7*h + 7*x^6*h^2 + 16*x^5*h^3 + 2*x^4*h^4 - 2*x^3*h^5 + 4*x^2*h^6 - 5*h^8, [ 1 -2] [ 1 -1] )
An example where the multiplicity is too high:
sage: R.<x,y> = PolynomialRing(QQ) sage: f = x^3 + 378666*x^2*y - 12444444*x*y^2 + 1234567890*y^3 sage: j = f * (x-545*y)^9 sage: j.reduced_form(prec=200, smallest_coeffs=False) Traceback (most recent call last): ... ValueError: cannot have a root with multiplicity >= 12/2
An example where Newton’s Method does not find the right root:
sage: R.<x,y> = PolynomialRing(QQ) sage: F = x^6 + 3*x^5*y - 8*x^4*y^2 - 2*x^3*y^3 - 44*x^2*y^4 - 8*x*y^5 sage: F.reduced_form(smallest_coeffs=False, prec=400) Traceback (most recent call last): ... ArithmeticError: Newton's method converged to z not in the upper half plane
An example with covariant on the boundary, therefore a non-unique form:
sage: R.<x,y> = PolynomialRing(QQ) sage: F = 5*x^2*y - 5*x*y^2 - 30*y^3 sage: F.reduced_form(smallest_coeffs=False) ( [1 1] 5*x^2*y + 5*x*y^2 - 30*y^3, [0 1] )
An example where precision needs to be increased:
sage: R.<x,y> = PolynomialRing(QQ) sage: F=-16*x^7 - 114*x^6*y - 345*x^5*y^2 - 599*x^4*y^3 - 666*x^3*y^4 - 481*x^2*y^5 - 207*x*y^6 - 40*y^7 sage: F.reduced_form(prec=50, smallest_coeffs=False) Traceback (most recent call last): ... ValueError: accuracy of Newton's root not within tolerance(0.0000124... > 1e-06), increase precision sage: F.reduced_form(prec=100, smallest_coeffs=False) ( [-1 -1] -x^5*y^2 - 24*x^3*y^4 - 3*x^2*y^5 - 2*x*y^6 + 16*y^7, [ 1 0] )
sage: R.<x,y> = PolynomialRing(QQ) sage: F = - 8*x^4 - 3933*x^3*y - 725085*x^2*y^2 - 59411592*x*y^3 - 1825511633*y^4 sage: F.reduced_form(return_conjugation=False) x^4 + 9*x^3*y - 3*x*y^3 - 8*y^4
sage: R.<x,y> = QQ[] sage: F = -2*x^3 + 2*x^2*y + 3*x*y^2 + 127*y^3 sage: F.reduced_form() ( [1 4] -2*x^3 - 22*x^2*y - 77*x*y^2 + 43*y^3, [0 1] )
sage: R.<x,y> = QQ[] sage: F = -2*x^3 + 2*x^2*y + 3*x*y^2 + 127*y^3 sage: F.reduced_form(norm_type='height') ( [5 4] -58*x^3 - 47*x^2*y + 52*x*y^2 + 43*y^3, [1 1] )
sage: R.<x,y,z> = PolynomialRing(QQ) sage: F = x^4 + x^3*y*z + y^2*z sage: F.reduced_form() Traceback (most recent call last): ... ValueError: (=x^3*y*z + x^4 + y^2*z) must have two variables
sage: R.<x,y> = PolynomialRing(ZZ) sage: F = - 8*x^6 - 3933*x^3*y - 725085*x^2*y^2 - 59411592*x*y^3 - 99*y^6 sage: F.reduced_form(return_conjugation=False) Traceback (most recent call last): ... ValueError: (=-8*x^6 - 99*y^6 - 3933*x^3*y - 725085*x^2*y^2 - 59411592*x*y^3) must be homogeneous
sage: R.<x,y> = PolynomialRing(RR) sage: F = 217.992172373276*x^3 + 96023.1505442490*x^2*y + 1.40987971253579e7*x*y^2\ + 6.90016027113216e8*y^3 sage: F.reduced_form(smallest_coeffs=False) # tol 1e-8 ( -39.5673942565918*x^3 + 111.874026298523*x^2*y + 231.052762985229*x*y^2 - 138.380829811096*y^3, [-147 -148] [ 1 1] )
sage: R.<x,y> = PolynomialRing(CC) sage: F = (0.759099196558145 + 0.845425869641446*CC.0)*x^3 + (84.8317207268542 + 93.8840848648033*CC.0)*x^2*y\ + (3159.07040755858 + 3475.33037377779*CC.0)*x*y^2 + (39202.5965389079 + 42882.5139724962*CC.0)*y^3 sage: F.reduced_form(smallest_coeffs=False) # tol 1e-11 ( (-0.759099196558145 - 0.845425869641446*I)*x^3 + (-0.571709908900118 - 0.0418133346027929*I)*x^2*y + (0.856525964330103 - 0.0721403997649759*I)*x*y^2 + (-0.965531044130330 + 0.754252314465703*I)*y^3, [-1 37] [ 0 -1] )
- specialization(D=None, phi=None)¶
Specialization of this polynomial.
Given a family of polynomials defined over a polynomial ring. A specialization is a particular member of that family. The specialization can be specified either by a dictionary or a
SpecializationMorphism
.INPUT:
D
– dictionary (optional)phi
– SpecializationMorphism (optional)
OUTPUT: a new polynomial
EXAMPLES:
sage: R.<c> = PolynomialRing(QQ) sage: S.<x,y> = PolynomialRing(R) sage: F = x^2 + c*y^2 sage: F.specialization({c:2}) x^2 + 2*y^2
sage: S.<a,b> = PolynomialRing(QQ) sage: P.<x,y,z> = PolynomialRing(S) sage: RR.<c,d> = PolynomialRing(P) sage: f = a*x^2 + b*y^3 + c*y^2 - b*a*d + d^2 - a*c*b*z^2 sage: f.specialization({a:2, z:4, d:2}) (y^2 - 32*b)*c + b*y^3 + 2*x^2 - 4*b + 4
Check that we preserve multi- versus uni-variate:
sage: R.<l> = PolynomialRing(QQ, 1) sage: S.<k> = PolynomialRing(R) sage: K.<a, b, c> = PolynomialRing(S) sage: F = a*k^2 + b*l + c^2 sage: F.specialization({b:56, c:5}).parent() Univariate Polynomial Ring in a over Univariate Polynomial Ring in k over Multivariate Polynomial Ring in l over Rational Field
- subresultants(other, variable=None)¶
Return the nonzero subresultant polynomials of
self
andother
.INPUT:
other
– a polynomial
OUTPUT: a list of polynomials in the same ring as
self
EXAMPLES:
sage: R.<x,y> = QQ[] sage: p = (y^2 + 6)*(x - 1) - y*(x^2 + 1) sage: q = (x^2 + 6)*(y - 1) - x*(y^2 + 1) sage: p.subresultants(q, y) [2*x^6 - 22*x^5 + 102*x^4 - 274*x^3 + 488*x^2 - 552*x + 288, -x^3 - x^2*y + 6*x^2 + 5*x*y - 11*x - 6*y + 6] sage: p.subresultants(q, x) [2*y^6 - 22*y^5 + 102*y^4 - 274*y^3 + 488*y^2 - 552*y + 288, x*y^2 + y^3 - 5*x*y - 6*y^2 + 6*x + 11*y - 6]
- sylvester_matrix(right, variable=None)¶
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:
self , right: multivariate polynomials
variable: optional, compute the Sylvester matrix with respect to this variable. If variable is not provided, the first variable of the polynomial ring is used.
OUTPUT:
The Sylvester matrix of self and right.
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: 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
- truncate(var, n)¶
Returns a new multivariate polynomial obtained from self by deleting all terms that involve the given variable to a power at least n.
- weighted_degree(*weights)¶
Return the weighted degree of
self
, which is the maximum weighted degree of all monomials inself
; 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 inweights
.This method is given for convenience. It is faster to use polynomial rings with weighted term orders and the standard
degree
function.INPUT:
weights
- Either individual numbers, an iterable or a dictionary, specifying the weights of each variable. If it is a dictionary, it maps each variable ofself
to its weight. If it is a sequence of individual numbers or a tuple, the weights are specified in the order of the generators as given byself.parent().gens()
:
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: unable to convert non-constant polynomial x to Integer Ring sage: p.weighted_degree(2/1,1,1) 6
The
weighted_degree
coincides with thedegree
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
- sage.rings.polynomial.multi_polynomial.is_MPolynomial(x)¶