Generic sparse distributed multivariate polynomials
AbstractAlgebra.jl provides a module, implemented in src/generic/MPoly.jl
for generic sparse distributed multivariate polynomials over any commutative ring belonging to the AbstractAlgebra abstract type hierarchy.
This modules implements the Multivariate Polynomial interface, including the sparse distributed, random access part of the interface.
All of the generic functionality is part of a submodule of AbstractAlgebra called Generic
. This is exported by default so that it is not necessary to qualify the function names with the submodule name.
Multivariates are implemented in this module using a Julia array of coefficients and a 2-dimensional Julia array of UInt
s for the exponent vectors. Note that exponent $n$ is represented by the $n$-th column of the exponent array, not the $n$-th row. This is because Julia uses a column major representation.
Types and parent objects
Multivariate polynomials implemented in AbstractAlgebra.jl have type Generic.MPoly{T}
where T
is the type of elements of the coefficient ring.
The polynomials are implemented using a Julia array of coefficients and a 2-dimensional Julia array of UInt
s for the exponent vectors. Note that exponent $n$ is represented by the $n$-th column of the exponent array, not the $n$-th row. This is because Julia uses a column major representation. See the file src/generic/GenericTypes.jl
for details.
The top bit of each UInt
is reserved for overflow detection.
Parent objects of such polynomials have type Generic.MPolyRing{T}
.
The string representation of the variables of the polynomial ring and the base/coefficient ring $R$ and the ordering are stored in the parent object.
The polynomial element types belong to the abstract type AbstractAlgebra.MPolyElem{T}
and the polynomial ring types belong to the abstract type AbstractAlgebra.MPolyRing{T}
.
Note that both the generic polynomial ring type Generic.MPolyRing{T}
and the abstract type it belongs to, AbstractAlgebra.MPolyRing{T}
are both called MPolyRing
. The former is a (parameterised) concrete type for a polynomial ring over a given base ring whose elements have type T
. The latter is an abstract type representing all multivariate polynomial ring types in AbstractAlgebra.jl, whether generic or very specialised (e.g. supplied by a C library).
Polynomial ring constructors
In order to construct multivariate polynomials in AbstractAlgebra.jl, one must first construct the polynomial ring itself. This is accomplished with the following constructor.
PolynomialRing(R::AbstractAlgebra.Ring, S::Array{String, 1}; cached::Bool = true, ordering::Symbol=:lex)
Given a base ring R
and and array S
of strings specifying how the generators (variables) should be printed, return a tuple S, (x, ...)
representing the new polynomial ring $S = R[x, \ldots]$ and a tuple of the generators $(x, ...)$ of the ring. By default the parent object S
will depend only on R
and (x, ...)
and will be cached. Setting the optional argument cached
to false
will prevent the parent object S
from being cached.
The optional named argument ordering
can be used to specify an ordering. The currently supported options are :lex
, :deglex
and :degrevlex
.
Like for univariate polynomials, a shorthand version of this function is provided when the number of generators is greater than 1
: given a base ring R
, we abbreviate the constructor as follows:
R["x", "y", ...]
Here are some examples of creating multivariate polynomial rings and making use of the resulting parent objects to coerce various elements into the polynomial ring.
Examples
julia> R, (x, y) = PolynomialRing(ZZ, ["x", "y"]; ordering=:deglex)
(Multivariate Polynomial Ring in x, y over Integers, AbstractAlgebra.Generic.MPoly{BigInt}[x, y])
julia> T, (z, t) = QQ["z", "t"]
(Multivariate Polynomial Ring in z, t over Rationals, AbstractAlgebra.Generic.MPoly{Rational{BigInt}}[z, t])
julia> f = R()
0
julia> g = R(123)
123
julia> h = R(BigInt(1234))
1234
julia> k = R(x + 1)
x + 1
julia> m = R(x + y + 1)
x + y + 1
julia> derivative(k, 1)
1
julia> derivative(k, 2)
0
All of the examples here are generic polynomial rings, but specialised implementations of polynomial rings provided by external modules will also usually provide a PolynomialRing
constructor to allow creation of their polynomial rings.
Polynomial functionality provided by AbstractAlgebra.jl
Basic manipulation
AbstractAlgebra.Generic.vars
— Method.vars(p::AbstractAlgebra.MPolyElem{T}) where {T <: RingElement}
Return the variables actually occuring in $p$.
AbstractAlgebra.Generic.var_index
— Method.var_index(p::AbstractAlgebra.MPolyElem{T}) where {T <: RingElement}
Return the index of the given variable $x$. If $x$ is not a variable in a multivariate polynomial ring, an exception is raised.
AbstractAlgebra.Generic.degree
— Method.degree(f::AbstractAlgebra.MPolyElem{T}, i::Int) where T <: RingElement
Return the degree of the polynomial $f$ in terms of the i-th variable.
AbstractAlgebra.Generic.degree
— Method.degree(f::AbstractAlgebra.MPolyElem{T}, x::AbstractAlgebra.MPolyElem{T}) where T <: RingElement
Return the degree of the polynomial $f$ in terms of the variable $x$.
AbstractAlgebra.Generic.degrees
— Method.degrees(f::AbstractAlgebra.MPolyElem{T}) where T <: RingElement
Return an array of the degrees of the polynomial $f$ in terms of each variable.
AbstractAlgebra.Generic.isconstant
— Method.isconstant(x::AbstractAlgebra.MPolyElem{T}) where T <: RingElement
Return
true
ifx
is a degree zero polynomial or the zero polynomial, i.e. a constant polynomial.
AbstractAlgebra.Generic.isterm
— Method.isterm(x::MPoly)
Return
true
if the given polynomial has precisely one term.
AbstractAlgebra.ismonomial
— Method.ismonomial(x::AbstractAlgebra.MPolyElem)
Return
true
if the given polynomial has precisely one term whose coefficient is one.
AbstractAlgebra.Generic.coeff
— Method.coeff(f::AbstractAlgebra.MPolyElem{T}, m::AbstractAlgebra.MPolyElem{T}) where T <: RingElement
Return the coefficient of the monomial $m$ of the polynomial $f$. If there is no such monomial, zero is returned.
Examples
julia> R, (x, y) = PolynomialRing(ZZ, ["x", "y"])
(Multivariate Polynomial Ring in x, y over Integers, AbstractAlgebra.Generic.MPoly{BigInt}[x, y])
julia> f = x^2 + 2x + 1
x^2 + 2*x + 1
julia> V = vars(f)
1-element Array{AbstractAlgebra.Generic.MPoly{BigInt},1}:
x
julia> var_index(y) == 2
true
julia> degree(f, x) == 2
true
julia> degree(f, 2) == 0
true
julia> d = degrees(f)
2-element Array{Int64,1}:
2
0
julia> isconstant(R(1))
true
julia> isterm(2x)
true
julia> ismonomial(y)
true
julia> isunit(R(1))
true
julia> c = coeff(f, x^2)
1
Changing base (coefficient) rings
In order to substitute the variables of a polynomial $f$ over a ring $T$ by elements in a $T$-algebra $S$, you first have to change the base ring of $f$ using the following function, where $g$ is a function representing the structure homomorphism of the $T$-algebra $S$.
AbstractAlgebra.change_base_ring
— Method.change_base_ring(R::Ring, p::MPolyElem{<: RingElement}; parent::MPolyRing, cached::Bool)
Return the polynomial obtained by coercing the non-zero coefficients of
p
intoR
.If the optional
parent
keyword is provided, the polynomial will be an element ofparent
. The caching of the parent object can be controlled via thecached
keyword argument.
AbstractAlgebra.Generic.map_coeffs
— Method.map_coeffs(f, p::MPolyElem{<: RingElement}; parent::MPolyRing)
Transform the polynomial
p
by applyingf
on each non-zero coefficient.If the optional
parent
keyword is provided, the polynomial will be an element ofparent
. The caching of the parent object can be controlled via thecached
keyword argument.
Examples
julia> R, (x, y) = PolynomialRing(ZZ, ["x", "y"])
(Multivariate Polynomial Ring in x, y over Integers, AbstractAlgebra.Generic.MPoly{BigInt}[x, y])
julia> fz = x^2*y^2 + x + 1
x^2*y^2 + x + 1
julia> fq = change_base_ring(QQ, fz)
x^2*y^2 + x + 1
In case a specific parent ring is constructed, it can also be passed to the function.
Examples
julia> R, (x, y) = PolynomialRing(ZZ, ["x", "y"])
(Multivariate Polynomial Ring in x, y over Integers, AbstractAlgebra.Generic.MPoly{BigInt}[x, y])
julia> S, = PolynomialRing(QQ, ["x", "y"])
(Multivariate Polynomial Ring in x, y over Rationals, AbstractAlgebra.Generic.MPoly{Rational{BigInt}}[x, y])
julia> fz = x^5 + y^3 + 1
x^5 + y^3 + 1
julia> fq = change_base_ring(QQ, fz, parent=S)
x^5 + y^3 + 1
Multivariate coefficients
In order to return the "coefficient" (as a multivariate polynomial in the same ring), of a given monomial (in which some of the variables may not appear and others may be required to appear to exponent zero), we can use the following function.
AbstractAlgebra.Generic.coeff
— Method.coeff(a::AbstractAlgebra.MPolyElem{T}, vars::Vector{Int}, exps::Vector{Int}) where T <: RingElement
Return the "coefficient" of $a$ (as a multivariate polynomial in the same ring) of the monomial consisting of the product of the variables of the given indices raised to the given exponents (note that not all variables need to appear and the exponents can be zero). E.g.
coeff(f, [1, 3], [0, 2])
returns the coefficient of $x^0*z^2$ in the polynomial $f$ (assuming variables $x, y, z$ in that order).
AbstractAlgebra.Generic.coeff
— Method.coeff(a::T, vars::Vector{T}, exps::Vector{Int}) where T <: AbstractAlgebra.MPolyElem
Return the "coefficient" of $a$ (as a multivariate polynomial in the same ring) of the monomial consisting of the product of the given variables to the given exponents (note that not all variables need to appear and the exponents can be zero). E.g.
coeff(f, [x, z], [0, 2])
returns the coefficient of $x^0*z^2$ in the polynomial $f$.
Examples
julia> R, (x, y, z) = PolynomialRing(ZZ, ["x", "y", "z"])
(Multivariate Polynomial Ring in x, y, z over Integers, AbstractAlgebra.Generic.MPoly{BigInt}[x, y, z])
julia> f = x^4*y^2*z^2 - 2x^4*y*z^2 + 4x^4*z^2 + 2x^2*y^2 + x + 1
x^4*y^2*z^2 - 2*x^4*y*z^2 + 4*x^4*z^2 + 2*x^2*y^2 + x + 1
julia> coeff(f, [1, 3], [4, 2]) == coeff(f, [x, z], [4, 2])
true
Inflation/deflation
AbstractAlgebra.Generic.deflation
— Method.deflation(f::AbstractAlgebra.MPolyElem{T}) where T <: RingElement
Compute deflation parameters for the exponents of the polynomial $f$. This is a pair of arrays of integers, the first array of which (the shift) gives the minimum exponent for each variable of the polynomial, and the second of which (the deflation) gives the gcds of all the exponents after subtracting the shift, again per variable. This functionality is used by gcd (and can be used by factorisation algorithms).
AbstractAlgebra.Generic.deflate
— Method.deflate(f::AbstractAlgebra.MPolyElem{T}, v::Vector{Int}) where T <: RingElement
Return a polynomial with the same coefficients as $f$ but whose exponents have been shifted down by the given shifts (supplied as an array of shifts, one for each variable, then deflated (divided) by the given exponents (again supplied as an array of deflation factors, one for each variable). The algorithm automatically replaces a deflation of $0$ by $1$, to avoid division by $0$.
AbstractAlgebra.Generic.inflate
— Method.inflate(f::AbstractAlgebra.MPolyElem{T}, v::Vector{Int}) where T <: RingElement
Return a polynomial with the same coefficients as $f$ but whose exponents have been inflated (multiplied) by the given deflation exponents (supplied as an array of inflation factors, one for each variable) and then shifted by the given shifts (again supplied as an array of shifts, one for each variable).
Examples
julia> R, (x, y) = PolynomialRing(ZZ, ["x", "y"])
(Multivariate Polynomial Ring in x, y over Integers, AbstractAlgebra.Generic.MPoly{BigInt}[x, y])
julia> f = x^7*y^8 + 3*x^4*y^8 - x^4*y^2 + 5x*y^5 - x*y^2
x^7*y^8 + 3*x^4*y^8 - x^4*y^2 + 5*x*y^5 - x*y^2
julia> def, shift = deflation(f)
([1, 2], [3, 3])
julia> f1 = deflate(f, def, shift)
x^2*y^2 + 3*x*y^2 - x + 5*y - 1
julia> f2 = inflate(f1, def, shift)
x^7*y^8 + 3*x^4*y^8 - x^4*y^2 + 5*x*y^5 - x*y^2
julia> f2 == f
true
Conversions
AbstractAlgebra.Generic.to_univariate
— Method.to_univariate(R::AbstractAlgebra.PolyRing{T}, p::AbstractAlgebra.MPolyElem{T}) where T <: AbstractAlgebra.RingElement
Assuming the polynomial $p$ is actually a univariate polynomial, convert the polynomial to a univariate polynomial in the given univariate polynomial ring $R$. An exception is raised if the polynomial $p$ involves more than one variable.
Examples
julia> R, (x, y) = PolynomialRing(ZZ, ["x", "y"])
(Multivariate Polynomial Ring in x, y over Integers, AbstractAlgebra.Generic.MPoly{BigInt}[x, y])
julia> S, z = PolynomialRing(ZZ, "z")
(Univariate Polynomial Ring in z over Integers, z)
julia> f = 2x^5 + 3x^4 - 2x^2 - 1
2*x^5 + 3*x^4 - 2*x^2 - 1
julia> g = to_univariate(S, f)
2*z^5 + 3*z^4 - 2*z^2 - 1
Evaluation
The following function allows evaluation of a polynomial at all its variables. The result is always in the ring that a product of a coefficient and one of the values belongs to, i.e. if all the values are in the coefficient ring, the result of the evaluation will be too.
AbstractAlgebra.Generic.evaluate
— Method.evaluate(a::AbstractAlgebra.MPolyElem{T}, vals::Vector{U}) where {T <: RingElement, U <: RingElement}
Evaluate the polynomial by substituting in the array of values for each of the variables. The evaluation will succeed if multiplication is defined between elements of the coefficient ring of $a$ and elements of the supplied vector.
The following functions allow evaluation of a polynomial at some of its variables. Note that the result will be a product of values and an element of the polynomial ring, i.e. even if all the values are in the coefficient ring and all variables are given values, the result will be a constant polynomial, not a coefficient.
AbstractAlgebra.Generic.evaluate
— Method.evaluate(a::AbstractAlgebra.MPolyElem{T}, vars::Vector{Int}, vals::Vector{U}) where {T <: RingElement, U <: RingElement}
Evaluate the polynomial by substituting in the supplied values in the array
vals
for the corresponding variables with indices given by the arrayvars
. The evaluation will succeed if multiplication is defined between elements of the coefficient ring of $a$ and elements ofvals
.
AbstractAlgebra.Generic.evaluate
— Method.evaluate(a::S, vars::Vector{S}, vals::Vector{U}) where {S <: AbstractAlgebra.MPolyElem{T}, U <: RingElement} where T <: RingElement
Evaluate the polynomial by substituting in the supplied values in the array
vals
for the corresponding variables (supplied as polynomials) given by the arrayvars
. The evaluation will succeed if multiplication is defined between elements of the coefficient ring of $a$ and elements ofvals
.
The following function allows evaluation of a polynomial at values in a not necessarily commutative ring, e.g. elements of a matrix algebra.
AbstractAlgebra.Generic.evaluate
— Method.evaluate(a::AbstractAlgebra.MPolyElem{T}, vals::Vector{U}) where {T <: RingElement, U <: NCRingElem}
Evaluate the polynomial at the supplied values, which may be any ring elements, commutative or non-commutative, but in the same ring. Evaluation always proceeds in the order of the variables as supplied when creating the polynomial ring to which $a$ belongs. The evaluation will succeed if a product of a coefficient of the polynomial by one of the values is defined.
Examples
julia> R, (x, y) = PolynomialRing(ZZ, ["x", "y"])
(Multivariate Polynomial Ring in x, y over Integers, AbstractAlgebra.Generic.MPoly{BigInt}[x, y])
julia> f = 2x^2*y^2 + 3x + y + 1
2*x^2*y^2 + 3*x + y + 1
julia> evaluate(f, BigInt[1, 2])
14
julia> evaluate(f, [QQ(1), QQ(2)])
14//1
julia> evaluate(f, [1, 2])
14
julia> f(1, 2) == 14
true
julia> evaluate(f, [x + y, 2y - x])
2*x^4 - 4*x^3*y - 6*x^2*y^2 + 8*x*y^3 + 2*x + 8*y^4 + 5*y + 1
julia> f(x + y, 2y - x)
2*x^4 - 4*x^3*y - 6*x^2*y^2 + 8*x*y^3 + 2*x + 8*y^4 + 5*y + 1
julia> R, (x, y, z) = PolynomialRing(ZZ, ["x", "y", "z"])
(Multivariate Polynomial Ring in x, y, z over Integers, AbstractAlgebra.Generic.MPoly{BigInt}[x, y, z])
julia> f = x^2*y^2 + 2x*z + 3y*z + z + 1
x^2*y^2 + 2*x*z + 3*y*z + z + 1
julia> evaluate(f, [1, 3], [3, 4])
9*y^2 + 12*y + 29
julia> evaluate(f, [x, z], [3, 4])
9*y^2 + 12*y + 29
julia> evaluate(f, [1, 2], [x + z, x - z])
x^4 - 2*x^2*z^2 + 5*x*z + z^4 - z^2 + z + 1
julia> S = MatrixAlgebra(ZZ, 2)
Matrix Algebra of degree 2 over Integers
julia> M1 = S([1 2; 3 4])
[1 2]
[3 4]
julia> M2 = S([2 3; 1 -1])
[2 3]
[1 -1]
julia> M3 = S([-1 1; 1 1])
[-1 1]
[ 1 1]
julia> evaluate(f, [M1, M2, M3])
[ 64 83]
[124 149]
Leading coefficients, leading monomials and leading terms
The leading coefficient, leading monomial and leading term of a polynomial p are returned by the following functions:
AbstractAlgebra.Generic.lc
— Method.lc(p::MPolyElem)
Return the leading coefficient of the polynomial p.
AbstractAlgebra.Generic.lm
— Method.lm(p::MPolyElem)
Return the leading monomial of the polynomial p.
AbstractAlgebra.Generic.lt
— Method.lt(p::MPolyElem)
Return the leading term of the polynomial p.
Examples
using AbstractAlgebra
R,(x,y) = PolynomialRing(ZZ, ["x", "y"], ordering=:deglex)
p = 2*x*y + 3*y^3
lt(p)
lm(p)
lc(p)
lt(p) == lc(p) * lm(p)
Least common multiple, greatest common divisor
The greated common divisor of two polynomials a and b is returned by
Base.gcd
— Method.gcd(a::AbstractAlgebra.Generic.MPoly{T}, a::AbstractAlgebra.Generic.MPoly{T}) where {T <: RingElement}
Return the greatest common divisor of a and b in parent(a).
Note that this functionality is currently only provided for AbstractAlgebra generic polynomials. It is not automatically provided for all multivariate rings that implement the multivariate interface.
However, if such a gcd is provided, the least common multiple of two polynomials a and b is returned by
Base.lcm
— Method.lcm(a::AbstractAlgebra.MPolyElem{T}, a::AbstractAlgebra.MPolyElem{T}) where {T <: RingElement}
Return the least common multiple of a and b in parent(a).
Examples
julia> using AbstractAlgebra
julia> R,(x,y) = PolynomialRing(ZZ, ["x", "y"])
(Multivariate Polynomial Ring in x, y over Integers, AbstractAlgebra.Generic.MPoly{BigInt}[x, y])
julia> a = x*y + 2*y
x*y + 2*y
julia> b = x^3*y + y
x^3*y + y
julia> gcd(a,b)
y
julia> lcm(a,b)
x^4*y + 2*x^3*y + x*y + 2*y
julia> lcm(a,b) == a * b // gcd(a,b)
true
Derivations
AbstractAlgebra.Generic.derivative
— Method.derivative(f::AbstractAlgebra.MPolyElem{T}, x::AbstractAlgebra.MPolyElem{T}) where T <: RingElement
Return the partial derivative of
f
with respect tox
. The valuex
must be a generator of the polynomial ring off
.
Examples
julia> R, (x, y) = AbstractAlgebra.PolynomialRing(ZZ, ["x", "y"])
(Multivariate Polynomial Ring in x, y over Integers, AbstractAlgebra.Generic.MPoly{BigInt}[x, y])
julia> f = x*y + x + y + 1
x*y + x + y + 1
julia> derivative(f, x)
y + 1
julia> derivative(f, y)
x + 1
Homogeneous polynomials
It is possible to test whether a polynomial is homogeneous with respect to the standard grading using the function
AbstractAlgebra.Generic.ishomogeneous
— Method.ishomogeneous(x::MPoly{T}) where {T <: RingElement}
Return
true
if the given polynomial is homogeneous with respect to the standard grading andfalse
otherwise.