Introduction

Nemo allow the creation of dense, univariate polynomials over any computable ring . There are two different kinds of implementation: a generic one for the case where no specific implementation exists, and efficient implementations of polynomials over numerous specific rings, usually provided by C/C++ libraries.

The following table shows each of the polynomial types available in Nemo, the base ring , and the Julia/Nemo types for that kind of polynomial (the type information is mainly of concern to developers).

Base ring Library Element type Parent type
Generic ring Nemo GenPoly{T} GenPolyRing{T}
Flint fmpz_poly FmpzPolyRing
(small ) Flint nmod_poly NmodPolyRing
(large ) Flint fmpz_mod_poly FmpzModPolyRing
Flint fmpq_poly FmpqPolyRing
(small ) Flint fq_nmod_poly FqNmodPolyRing
(large ) Flint fq_poly FqPolyRing
Arb arb_poly ArbPolyRing
Arb acb_poly AcbPolyRing

The string representation of the variable and the base ring of a generic polynomial is stored in its parent object.

All polynomial element types belong to the abstract type PolyElem and all of the polynomial ring types belong to the abstract type PolyRing. This enables one to write generic functions that can accept any Nemo polynomial type.

Polynomial ring constructors

In order to construct polynomials in Nemo, one must first construct the polynomial ring itself. This is accomplished with the following constructor.

PolynomialRing(::Ring, ::AbstractString, ::Bool)

A shorthand version of this function is provided: given a base ring R, we abbreviate the constructor as follows.

R["x"]

Here are some examples of creating polynomial rings and making use of the resulting parent objects to coerce various elements into the polynomial ring.

R, x = PolynomialRing(ZZ, "x")
S, y = PolynomialRing(R, "y")
T, z = QQ["z"]

f = R()
g = R(123)
h = S(ZZ(1234))
k = S(x + 1)
m = T(z + 1)

Polynomial element constructors

Once a polynomial ring is constructed, there are various ways to construct polynomials in that ring.

The easiest way is simply using the generator returned by the PolynomialRing constructor and and build up the polynomial using basic arithmetic. Julia has quite flexible notation for the construction of polynomials in this way.

In addition we provide the following functions for constructing certain useful polynomials.

# Base.zeroMethod.

zero(R::PolyRing)

Return the zero polynomial in the given polynomial ring.

source

# Base.oneMethod.

one(R::PolyRing)

Return the constant polynomial in the given polynomial ring.

source

# Nemo.genMethod.

gen(R::PolyRing)

Return the generator of the given polynomial ring.

source

Here are some examples of constructing polynomials.

R, x = PolynomialRing(ZZ, "x")
S, y = PolynomialRing(R, "y")

f = x^3 + 3x + 21
g = (x + 1)*y^2 + 2x + 1

h = zero(S)
k = one(R)
m = gen(S)

Basic functionality

All univariate polynomial modules in Nemo must provide the functionality listed in this section. (Note that only some of these functions are useful to a user.)

Developers who are writing their own polynomial module, whether as an interface to a C library, or as some kind of generic module, must provide all of these functions for custom univariate polynomial types in Nemo.

We write U for the type of the polynomials in the polynomial ring and T for the type of elements of the coefficient ring.

All of these functions are provided for all existing polynomial types in Nemo.

parent_type{U <: PolyElem}(::Type{U})

Given the type of polynomial elements, should return the type of the corresponding parent object.

elem_type(R::PolyRing)

Given a parent object for the polynomial ring, return the type of elements of the polynomial ring.

Base.hash(a::PolyElem, h::UInt)

Return a UInt hexadecimal hash of the polynomial . This should be xor'd with a fixed random hexadecimal specific to the polynomial type. The hash of each coefficient should be xor'd with the supplied parameter h as part of computing the hash.

fit!(a::PolyElem, n::Int)

By reallocating if necessary, ensure that the given polynomial has space for at least coefficients. This function does not change the length of the polynomial and will only ever increase the number of allocated coefficients. Any coefficients added by this function are initialised to zero.

normalise(a::PolyElem, n::Int)

Return the normalised length of the given polynomial, assuming its current length is . Its normalised length is such that it either has nonzero leading term or is the zero polynomial. Note that this function doesn't normalise the polynomial. That can be done with a subsequent call to set_length! using the length returned by normalise.

set_length!(a::PolyElem, n::Int)

Set the length of an existing polynomial that has sufficient space allocated, i.e. a polynomial for which no reallocation is needed. Note that if the Julia type definition for a custom polynomial type has a field, length, which corresponds to the current length of the polynomial, then the developer doesn't need to supply this function, as the supplied generic implementation will work. Note that it can change the length to any value from zero to the number of coefficients currently allocated and initialised.

length(a::PolyElem)

Return the current length (not the number of allocated coefficients), of the given polynomial. Note that this function only needs to be provided by a developer for a custom polynomial type if the Julia type definition for polynomial elements doesn't contain a field length corresponding to the current length of the polynomial. Otherwise the supplied generic implementation will work.

coeff(a::PolyElem, n::Int)

Return the degree n coefficient of the given polynomial. Note coefficients are numbered from n = 0 for the constant coefficient. If is bigger then the degree of the polynomial, the function returns a zero coefficient. We require .

setcoeff!{T <: RingElem}(a::PolyElem{T}, n::Int, c::T)

Set the coefficient of the degree term of the given polynomial to the given value a. The polynomial is not normalised automatically after this operation, however the polynomial is automatically resized if there is not sufficient allocated space.

deepcopy(a::PolyElem)

Construct a copy of the given polynomial and return it. This function must recursively construct copies of all of the internal data in the given polynomial. Nemo polynomials are mutable and so returning shallow copies is not sufficient.

mul!(c::PolyElem, a::PolyElem, b::PolyElem)

Multiply by and set the existing polynomial to the result. This function is provided for performance reasons as it saves allocating a new object for the result and eliminates associated garbage collection.

addeq!(c::PolyElem, a::PolyElem)

In-place addition. Adds to and sets to the result. This function is provided for performance reasons as it saves allocating a new object for the result and eliminates associated garbage collection.

Given a parent object S for a polynomial ring, the following coercion functions are provided to coerce various elements into the polynomial ring. Developers provide these by overloading the call operator for the polynomial parent objects.

S()

Coerce zero into the ring .

S(n::Integer)
S(n::fmpz)

Coerce an integer value or Flint integer into the polynomial ring .

S(n::T)

Coerces an element of the base ring, of type T into .

S(A::Array{T, 1})

Take an array of elements in the base ring, of type T and construct the polynomial with those coefficients, starting with the constant coefficient.

S(f::PolyElem)

Take a polynomial that is already in the ring and simply return it. A copy of the original is not made.

S(c::RingElem)

Try to coerce the given ring element into the polynomial ring. This only succeeds if can be coerced into the base ring.

In addition to the above, developers of custom polynomials must ensure the parent object of a polynomial type constains a field base_ring specifying the base ring, a field S containing a symbol (not a string) representing the variable name of the polynomial ring. They must also ensure that each polynomial element contains a field parent specifying the parent object of the polynomial.

Typically a developer will also overload the PolynomialRing generic function to create polynomials of the custom type they are implementing.

Basic manipulation

Numerous functions are provided to manipulate polynomials and to set and retrieve coefficients and other basic data associated with the polynomials. Also see the section on basic functionality above.

# Nemo.base_ringMethod.

base_ring(R::PolyRing)

Return the base ring of the given polynomial ring.

source

# Nemo.base_ringMethod.

base_ring(a::PolyElem)

Return the base ring of the polynomial ring of the given polynomial.

source

# Base.parentMethod.

parent(a::PolyElem)

Return the parent of the given polynomial.

source

# Base.varMethod.

var(a::PolyRing)

Return the internal name of the generator of the polynomial ring. Note that this is returned as a Symbol not a String.

source

# Nemo.varsMethod.

vars(a::PolyRing)

Return an array of the variable names for the polynomial ring. Note that this is returned as an array of Symbol not String.

source

# Nemo.degreeMethod.

degree(a::PolyElem)

Return the degree of the given polynomial. This is defined to be one less than the length, even for constant polynomials.

source

# Nemo.modulusMethod.

modulus{T <: ResElem}(a::PolyElem{T})

Return the modulus of the coefficients of the given polynomial.

source

# Nemo.leadMethod.

lead(x::PolyElem)

Return the leading coefficient of the given polynomial. This will be the nonzero coefficient of the term with highest degree unless the polynomial in the zero polynomial, in which case a zero coefficient is returned.

source

# Nemo.trailMethod.

trail(x::PolyElem)

Return the trailing coefficient of the given polynomial. This will be the nonzero coefficient of the term with lowest degree unless the polynomial in the zero polynomial, in which case a zero coefficient is returned.

source

# Nemo.iszeroMethod.

iszero(a::PolyElem)

Return true if the given polynomial is zero, otherwise return false.

source

# Nemo.isoneMethod.

isone(a::PolyElem)

Return true if the given polynomial is the constant polynomial , otherwise return false.

source

# Nemo.isgenMethod.

isgen(a::PolyElem)

Return true if the given polynomial is the constant generator of its polynomial ring, otherwise return false.

source

# Nemo.isunitMethod.

isunit(a::PolyElem)

Return true if the given polynomial is a unit in its polynomial ring, otherwise return false.

source

# Nemo.ismonomialMethod.

ismonomial(a::PolyElem)

Return true if the given polynomial is a monomial.

source

# Nemo.istermMethod.

isterm(a::PolyElem)

Return true if the given polynomial is has one term. This function is recursive, with all scalar types returning true.

source

# Base.denMethod.

den(a::fmpq_poly)

Return the least common denominator of the coefficients of the polynomial .

source

Here are some examples of basic manipulation of polynomials.

R, x = PolynomialRing(ZZ, "x")
S, y = PolynomialRing(R, "y")
T, z = PolynomialRing(QQ, "z")

a = zero(S)
b = one(S)

c = ZZ(1)//2*z^2 + ZZ(1)//3
d = x*y^2 + (x + 1)*y + 3

U = base_ring(S)
V = base_ring(y + 1)
v = var(S)
T = parent(y + 1)

f = lead(d)

g = isgen(y)
h = isone(b)
k = iszero(a)
m = isunit(b)
n = degree(d)
p = length(b)
q = den(c)

Arithmetic operators

All the usual arithmetic operators are overloaded for Nemo polynomials. Note that Julia uses the single slash for floating point division. Therefore to perform exact division in a ring we use divexact. To construct an element of a fraction field one can use the double slash operator //.

The following standard operators and functions are provided.

Function Operation
-(a::PolyElem) unary minus
+{T <: RingElem}(a::PolyElem{T}, b::PolyElem{T}) addition
-{T <: RingElem}(a::PolyElem{T}, b::PolyElem{T}) subtraction
*{T <: RingElem}(a::PolyElem{T}, b::PolyElem{T}) multiplication
divexact{T <: RingElem}(a::PolyElem{T}, b::PolyElem{T}) exact division

The following ad hoc operators and functions are also provided.

Function Operation
+(a::Integer, b::PolyElem) addition
+(a::PolyElem, b::Integer) addition
+(a::fmpz, b::PolyElem) addition
+(a::PolyElem, b::fmpz) addition
+{T <: RingElem}(a::T, b::PolyElem{T}) addition
+{T <: RingElem}(a::PolyElem{T}, b::T) addition
-(a::Integer, b::PolyElem) subtraction
-(a::PolyElem, b::Integer) subtraction
-(a::fmpz, b::PolyElem) subtraction
-(a::PolyElem, b::fmpz) subtraction
-{T <: RingElem}(a::T, b::PolyElem{T}) subtraction
-{T <: RingElem}(a::PolyElem{T}, b::T) subtraction
*(a::Integer, b::PolyElem) multiplication
*(a::PolyElem, b::Integer) multiplication
*(a::fmpz, b::PolyElem) multiplication
*(a::PolyElem, b::fmpz) multiplication
*{T <: RingElem}(a::T, b::PolyElem{T}) multiplication
*{T <: RingElem}(a::PolyElem{T}, b::T) multiplication
divexact(a::PolyElem, b::Integer) exact division
divexact(a::PolyElem, b::fmpz) exact division
divexact{T <: RingElem}(a::PolyElem{T}, b::T) exact division
^(a::PolyElem, n::Int) powering

If the appropriate promote_rule and coercion exists, these operators can also be used with elements of other rings. Nemo will try to coerce the operands to the dominating type and then apply the operator.

Here are some examples of arithmetic operations on polynomials.

R, x = PolynomialRing(ZZ, "x")
S, y = PolynomialRing(R, "y")

f = x*y^2 + (x + 1)*y + 3
g = (x + 1)*y + (x^3 + 2x + 2)

h = f - g
k = f*g
m = f + g
n = g - 4
p = fmpz(5) - g
q = f*7
r = divexact(f, -1)
s = divexact(g*(x + 1), x + 1)
t = f^3

Comparison operators

The following comparison operators are implemented for polynomials in Nemo. Julia provides the corresponding != operator automatically.

Function

isequal{T <: RingElem}(a::PolyElem{T}, b::PolyElem{T}) =={T <: RingElem}(a::PolyElem{T}, b::PolyElem{T})

The isequal operation returns true if and only if all the coefficients of the polynomial are precisely equal as compared by isequal. This is a stronger form of equality, used for comparing inexact coefficients, such as elements of a power series ring, the -adics, or the reals or complex numbers. Two elements are precisely equal only if they have the same precision or bounds in addition to being arithmetically equal.

We also have the following ad hoc comparison operators.

Function

=={T <: RingElem}(a::PolyElem{T}, b::T) =={T <: RingElem}(a::T, b::PolyElem{T}) ==(a::PolyElem, b::Integer) ==(a::Integer, b::PolyElem) ==(a::PolyElem, b::fmpz) ==(a::fmpz, b::PolyElem)

Here are some examples of comparisons.

R, x = PolynomialRing(ZZ, "x")
S, y = PolynomialRing(R, "y")

f = x*y^2 + (x + 1)*y + 3
g = x*y^2 + (x + 1)*y + 3
h = S(3)

f == g
isequal(f, g)
f != 3
g != x
h == fmpz(3)

Truncation

# Base.truncateMethod.

truncate(a::PolyElem, n::Int)

Return truncated to terms.

source

# Nemo.mullowMethod.

mullow{T <: RingElem}(a::PolyElem{T}, b::PolyElem{T}, n::Int)

Return truncated to terms.

source

Here are some examples of truncated operations.

R, x = PolynomialRing(ZZ, "x")
S, y = PolynomialRing(R, "y")

f = x*y^2 + (x + 1)*y + 3
g = (x + 1)*y + (x^3 + 2x + 2)

h = truncate(f, 1)
k = mullow(f, g, 4)

Reversal

# Base.reverseMethod.

reverse(x::PolyElem, len::Int)

Return the reverse of the polynomial , thought of as a polynomial of the given length (the polynomial will be notionally truncated or padded with zeroes before the leading term if necessary to match the specified length). The resulting polynomial is normalised. If len is negative we throw a DomainError().

source

# Base.reverseMethod.

reverse(x::PolyElem)

Return the reverse of the polynomial , i.e. the leading coefficient of becomes the constant coefficient of the result, etc. The resulting polynomial is normalised.

source

Here are some examples of reversal.

R, x = PolynomialRing(ZZ, "x")
S, y = PolynomialRing(R, "y")

f = x*y^2 + (x + 1)*y + 3

g = reverse(f, 7)
h = reverse(f)

Shifting

# Nemo.shift_leftMethod.

shift_left(x::PolyElem, n::Int)

Return the polynomial shifted left by terms, i.e. multiplied by .

source

# Nemo.shift_rightMethod.

shift_right(f::PolyElem, n::Int)

Return the polynomial shifted right by terms, i.e. divided by .

source

Here are some examples of shifting.

R, x = PolynomialRing(ZZ, "x")
S, y = PolynomialRing(R, "y")

f = x*y^2 + (x + 1)*y + 3

g = shift_left(f, 7)
h = shift_right(f, 2)

Modulo arithmetic

For polynomials over a field or residue ring, we can reduce modulo a given polynomial. This isn't always well-defined in the case of a residue ring, but when it is well-defined, we obtain the correct result. If Nemo encounters an impossible inverse, an exception will be raised.

# Nemo.mulmodMethod.

mulmod{T <: Union{ResElem, FieldElem}}(a::PolyElem{T}, b::PolyElem{T}, d::PolyElem{T})

Return .

source

# Nemo.powmodMethod.

powmod{T <: Union{ResElem, FieldElem}}(a::PolyElem{T}, b::Int, d::PolyElem{T})

Return . There are no restrictions on .

source

# Nemo.powmodMethod.

powmod(x::fmpz_mod_poly, e::fmpz, y::fmpz_mod_poly)

Return .

source

# Base.invmodMethod.

invmod{T <: Union{ResElem, FieldElem}}(a::PolyElem{T}, b::PolyElem{T})

Return .

source

Here are some examples of modular arithmetic.

R, x = PolynomialRing(QQ, "x")
S = ResidueRing(R, x^3 + 3x + 1)
T, y = PolynomialRing(S, "y")

f = (3*x^2 + x + 2)*y + x^2 + 1
g = (5*x^2 + 2*x + 1)*y^2 + 2x*y + x + 1
h = (3*x^3 + 2*x^2 + x + 7)*y^5 + 2x*y + 1

invmod(f, g)
mulmod(f, g, h)
powmod(f, 3, h)

Euclidean division

For polynomials over a field, we have a euclidean domain, and in many cases for polynomials over a residue ring things behave as though we had a euclidean domain so long as we don't hit an impossible inverse. For such rings we define euclidean division of polynomials. If an impossible inverse is hit, we raise an exception.

# Base.modMethod.

mod{T <: Union{ResElem, FieldElem}}(f::PolyElem{T}, g::PolyElem{T})

Return .

source

# Base.divremMethod.

divrem{T <: Union{ResElem, FieldElem}}(f::PolyElem{T}, g::PolyElem{T})

Return a tuple such that where is the euclidean quotient of by .

source

Here are some examples of euclidean division.

R = ResidueRing(ZZ, 7)
S, x = PolynomialRing(R, "x")
T = ResidueRing(S, x^3 + 3x + 1)
U, y = PolynomialRing(T, "y")

f = y^3 + x*y^2 + (x + 1)*y + 3
g = (x + 1)*y^2 + (x^3 + 2x + 2)

h = mod(f, g)
q, r = divrem(f, g)

Pseudodivision

Given two polynomials , pseudodivision computes polynomials and with length length such that where length length and is the leading coefficient of .

We call the pseudoquotient and the pseudoremainder.

# Nemo.pseudoremMethod.

pseudorem{T <: RingElem}(f::PolyElem{T}, g::PolyElem{T})

Return the pseudoremainder of divided by . If we throw a DivideError().

source

# Nemo.pseudodivremMethod.

pseudodivrem{T <: RingElem}(f::PolyElem{T}, g::PolyElem{T})

Return a tuple consisting of the pseudoquotient and pseudoremainder of divided by . If we throw a DivideError().

source

Here are some examples of pseudodivision.

R, x = PolynomialRing(ZZ, "x")
S, y = PolynomialRing(R, "y")

f = x*y^2 + (x + 1)*y + 3
g = (x + 1)*y + (x^3 + 2x + 2)

h = pseudorem(f, g)
q, r = pseudodivrem(f, g)

Content, primitive part, GCD and LCM

In Nemo, we allow computation of the greatest common divisor of polynomials over any ring. This is enabled by making use of pseudoremainders when we aren't working over a euclidean domain or something mimicking such a domain. In certain cases this allows us to return a greatest common divisor when it otherwise wouldn't be possible. However, a greatest common divisor is not necessarily unique, or even well-defined.

If an impossible inverse is encountered whilst computing the greatest common divisor, an exception is thrown.

# Base.gcdMethod.

gcd{T <: RingElem}(a::PolyElem{T}, b::PolyElem{T})

Return a greatest common divisor of and if it exists.

source

# Base.lcmMethod.

lcm{T <: RingElem}(a::PolyElem{T}, b::PolyElem{T})

Return a least common multiple of and if it exists.

source

# Nemo.contentMethod.

content(a::PolyElem)

Return the content of , i.e. the greatest common divisor of its coefficients.

source

# Nemo.primpartMethod.

primpart(a::PolyElem)

Return the primitive part of , i.e. the polynomial divided by its content.

source

# Base.gcdxMethod.

gcdx{T <: RingElem}(a::PolyElem{T}, b::PolyElem{T})

Return a tuple such that is the resultant of and and such that .

source

# Base.gcdxMethod.

gcdx{T <: Union{ResElem, FieldElem}}(a::PolyElem{T}, b::PolyElem{T})

Return a tuple such that is the greatest common divisor of and and such that .

source

# Nemo.gcdinvMethod.

gcdinv{T <: Union{ResElem, FieldElem}}(a::PolyElem{T}, b::PolyElem{T})

Return a tuple such that is the greatest common divisor of and and such that . This function is useful for inverting modulo a polynomial and checking that it really was invertible.

source

Here are some examples of content, primitive part and GCD.

R, x = PolynomialRing(ZZ, "x")
S, y = PolynomialRing(R, "y")

k = x*y^2 + (x + 1)*y + 3
l = (x + 1)*y + (x^3 + 2x + 2)
m = y^2 + x + 1

n = content(k)
p = primpart(k*(x^2 + 1))
q = gcd(k*m, l*m)
r = lcm(k*m, l*m)

R, x = PolynomialRing(QQ, "x")
T = ResidueRing(R, x^3 + 3x + 1)
U, z = PolynomialRing(T, "z")

g = z^3 + 2z + 1
h = z^5 + 1

r, s, t = gcdx(g, h)
u, v = gcdinv(g, h)

Evaluation, composition and substitution

# Nemo.evaluateMethod.

evaluate{T <: RingElem}(a::PolyElem{T}, b::T)

Evaluate the polynomial at the value and return the result.

source

# Nemo.evaluateMethod.

evaluate(a::PolyElem, b::Integer)

Evaluate the polynomial at the value and return the result.

source

# Nemo.evaluateMethod.

evaluate(a::PolyElem, b::fmpz)

Evaluate the polynomial at the value and return the result.

source

Remove and valuation

# Nemo.removeMethod.

remove{T <: RingElem}(z::PolyElem{T}, p::PolyElem{T})

Computes the valuation of at , that is, the largest such that divides . Additionally, is returned as well.

See also valuation, which only returns the valuation.

source

# Nemo.valuationMethod.

valuation{T <: RingElem}(z::PolyElem{T}, p::PolyElem{T})

Computes the valuation of at , that is, the largest such that divides .

See also remove, which also returns .

source

# Nemo.dividesMethod.

divides{T <: RingElem}(f::PolyElem{T}, g::PolyElem{T})

Returns a pair consisting of a flag which is set to true if divides and false otherwise, and a polynomial such that if such a polynomial exists. If not, the value of is undetermined.

source

# Nemo.dividesMethod.

divides{T <: RingElem}(f::PolyElem{T}, g::T)

Returns a pair consisting of a flag which is set to true if divides and false otherwise, and a polynomial such that if such a polynomial exists. If not, the value of is undetermined.

source

# Nemo.evaluate2Method.

evaluate2(x::arb_poly, y::Integer)

Return a tuple consisting of the polynomial evaluated at and its derivative evaluated at .

source

# Nemo.evaluate2Method.

evaluate2(x::arb_poly, y::Float64)

Return a tuple consisting of the polynomial evaluated at and its derivative evaluated at .

source

# Nemo.evaluate2Method.

evaluate2(x::arb_poly, y::fmpz)

Return a tuple consisting of the polynomial evaluated at and its derivative evaluated at .

source

# Nemo.evaluate2Method.

evaluate2(x::arb_poly, y::fmpq)

Return a tuple consisting of the polynomial evaluated at and its derivative evaluated at .

source

# Nemo.evaluate2Method.

evaluate2(x::arb_poly, y::arb)

Return a tuple consisting of the polynomial evaluated at and its derivative evaluated at .

source

# Nemo.evaluate2Method.

evaluate2(x::arb_poly, y::acb)

Return a tuple consisting of the polynomial evaluated at and its derivative evaluated at .

source

# Nemo.evaluate2Method.

evaluate2(x::acb_poly, y::Integer)

Return a tuple consisting of the polynomial evaluated at and its derivative evaluated at .

source

# Nemo.evaluate2Method.

evaluate2(x::acb_poly, y::Float64)

Return a tuple consisting of the polynomial evaluated at and its derivative evaluated at .

source

# Nemo.evaluate2Method.

evaluate2(x::acb_poly, y::fmpq)

Return a tuple consisting of the polynomial evaluated at and its derivative evaluated at .

source

# Nemo.evaluate2Method.

evaluate2(x::acb_poly, y::fmpq)

Return a tuple consisting of the polynomial evaluated at and its derivative evaluated at .

source

# Nemo.evaluate2Method.

evaluate2(x::acb_poly, y::arb)

Return a tuple consisting of the polynomial evaluated at and its derivative evaluated at .

source

# Nemo.evaluate2Method.

evaluate2(x::acb_poly, y::acb)

Return a tuple consisting of the polynomial evaluated at and its derivative evaluated at .

source

# Nemo.composeMethod.

compose(a::PolyElem, b::PolyElem)

Compose the polynomial with the polynomial and return the result, i.e. return .

source

# Nemo.substMethod.

subst{T <: RingElem}(f::PolyElem{T}, a::Any)

Evaluate the polynomial at . Note that can be anything, whether a ring element or not.

source

We also overload the functional notation so that the polynomial can be evaluated at by writing . This feature is only available with Julia 0.5 however.

Here are some examples of polynomial evaluation, composition and substitution.

RR = RealField(64)
R, x = PolynomialRing(ZZ, "x")
S, y = PolynomialRing(R, "y")
T, z = PolynomialRing(RR, "z")

f = x*y^2 + (x + 1)*y + 3
g = (x + 1)*y + (x^3 + 2x + 2)
h = z^2 + 2z + 1
M = R[x + 1 2x; x - 3 2x - 1]

k = evaluate(f, 3)
m = evaluate(f, x^2 + 2x + 1)
n = compose(f, g)
p = subst(f, M)
q = f(M)
r = f(23)
s, t = evaluate2(h, RR("2.0 +/- 0.1"))

Derivative and integral

# Nemo.derivativeMethod.

derivative(a::PolyElem)

Return the derivative of the polynomial .

source

# Nemo.integralMethod.

integral{T <: Union{ResElem, FieldElem}}(x::PolyElem{T})

Return the integral of the polynomial .

source

Here are some examples of integral and derivative.

R, x = PolynomialRing(ZZ, "x")
S, y = PolynomialRing(R, "y")
T, z = PolynomialRing(QQ, "z")
U = ResidueRing(T, z^3 + 3z + 1)
V, w = PolynomialRing(U, "w")

f = x*y^2 + (x + 1)*y + 3
g = (z^2 + 2z + 1)*w^2 + (z + 1)*w - 2z + 4

h = derivative(f)
k = integral(g)   

Resultant and discriminant

# Nemo.resultantMethod.

resultant{T <: RingElem}(a::PolyElem{T}, b::PolyElem{T})

Return the resultant of the and .

source

# Nemo.discriminantMethod.

discriminant(a::PolyElem)

Return the discrimnant of the given polynomial.

source

Here are some examples of computing the resultant and discriminant.

R, x = PolynomialRing(ZZ, "x")
S, y = PolynomialRing(R, "y")

f = 3x*y^2 + (x + 1)*y + 3
g = 6(x + 1)*y + (x^3 + 2x + 2)

h = resultant(f, g)
k = discriminant(f)

Newton representation

# Nemo.monomial_to_newton!Method.

monomial_to_newton!{T <: RingElem}(P::Array{T, 1}, roots::Array{T, 1})

Converts a polynomial , given as an array of coefficients, in-place from its coefficients given in the standard monomial basis to the Newton basis for the roots . In other words, this determines output coefficients such that is equal to the input polynomial.

source

# Nemo.newton_to_monomial!Method.

newton_to_monomial!{T <: RingElem}(P::Array{T, 1}, roots::Array{T, 1})

Converts a polynomial , given as an array of coefficients, in-place from its coefficients given in the Newton basis for the roots to the standard monomial basis. In other words, this evaluates where are the input coefficients given by .

source

Here are some examples of conversion to and from Newton representation.

R, x = PolynomialRing(ZZ, "x")
S, y = PolynomialRing(R, "y")

f = 3x*y^2 + (x + 1)*y + 3
g = deepcopy(f)
roots = [R(1), R(2), R(3)]

monomial_to_newton!(g.coeffs, roots)
newton_to_monomial!(g.coeffs, roots)

Multipoint evaluation and interpolation

# Nemo.interpolateMethod.

interpolate{T <: RingElem}(S::PolyRing, x::Array{T, 1}, y::Array{T, 1})

Given two arrays of values and of the same length , find the polynomial in the polynomial ring of length at most such that has the value at the points . The values in the arrays and must belong to the base ring of the polynomial ring . If no such polynomial exists, an exception is raised.

source

Here is an example of interpolation.

R, x = PolynomialRing(ZZ, "x")
S, y = PolynomialRing(R, "y")

xs = [R(1), R(2), R(3), R(4)]
ys = [R(1), R(4), R(9), R(16)]

f = interpolate(S, xs, ys)

Signature

Signature is only available for certain coefficient rings.

# Nemo.signatureMethod.

signature(f::fmpz_poly)

Return the signature of the polynomial , i.e. a tuple such that is the number of real roots of and is half the number of complex roots.

source

# Nemo.signatureMethod.

signature(f::fmpq_poly)

Return the signature of , i.e. a tuple where is the number of real roots of and is half the number of complex roots.

source

Here is an example of signature.

R, x = PolynomialRing(ZZ, "x")

f = x^3 + 3x + 1

(r, s) = signature(f)

Root finding

# Nemo.rootsMethod.

roots(x::acb_poly; target=0, isolate_real=false, initial_prec=0, max_prec=0, max_iter=0)

Attempts to isolate the complex roots of the complex polynomial by iteratively refining balls in which they lie.

This is done by increasing the working precision, starting at initial_prec. The maximal number of iterations can be set using max_iter and the maximal precision can be set using max_prec.

If isolate_real is set and is strictly real, then the real roots will be isolated from the non-real roots. Every root will have either zero, positive or negative real part.

It is assumed that is squarefree.

source

Here is an example of finding complex roots.

CC = ComplexField(64)
C, y = PolynomialRing(CC, "y")

m = y^2 + 2y + 3
n = m + CC("0 +/- 0.0001", "0 +/- 0.0001")

r = roots(n)

p = y^7 - 1

r = roots(n, isolate_real = true)

Construction from roots

# Nemo.from_rootsMethod.

from_roots(R::ArbPolyRing, b::Array{arb, 1})

Construct a polynomial in the given polynomial ring from a list of its roots.

source

# Nemo.from_rootsMethod.

from_roots(R::AcbPolyRing, b::Array{acb, 1})

Construct a polynomial in the given polynomial ring from a list of its roots.

source

Here are some examples of constructing polynomials from their roots.

RR = RealField(64)
R, x = PolynomialRing(RR, "x")

xs = arb[inv(RR(i)) for i=1:5]
f = from_roots(R, xs)

Lifting

When working over a residue ring it is useful to be able to lift to the base ring of the residue ring, e.g. from to .

# Nemo.liftMethod.

function lift(R::FmpzPolyRing, y::nmod_poly)

Lift from a polynomial over to a polynomial over with minimal reduced nonnegative coefficients. The ring R specifies the ring to lift into.

source

# Nemo.liftMethod.

function lift(R::FmpzPolyRing, y::fmpz_mod_poly)

Lift from a polynomial over to a polynomial over with minimal reduced nonnegative coefficients. The ring R specifies the ring to lift into.

source

Here is an example of lifting.

R = ResidueRing(ZZ, 123456789012345678949)
S, x = PolynomialRing(R, "x")
T, y = PolynomialRing(ZZ, "y")

f = x^2 + 2x + 1

a = lift(T, f)

Overlapping and containment

Occasionally it is useful to be able to tell when inexact polynomials overlap or contain other exact or inexact polynomials. The following functions are provided for this purpose.

# Nemo.overlapsMethod.

overlaps(x::arb_poly, y::arb_poly)

Return true if the coefficient balls of overlap the coefficient balls of , otherwise return false.

source

# Nemo.overlapsMethod.

overlaps(x::acb_poly, y::acb_poly)

Return true if the coefficient boxes of overlap the coefficient boxes of , otherwise return false.

source

# Base.containsMethod.

contains(x::arb_poly, y::arb_poly)

Return true if the coefficient balls of contain the corresponding coefficient balls of , otherwise return false.

source

# Base.containsMethod.

contains(x::acb_poly, y::acb_poly)

Return true if the coefficient boxes of contain the corresponding coefficient boxes of , otherwise return false.

source

# Base.containsMethod.

contains(x::arb_poly, y::fmpz_poly)

Return true if the coefficient balls of contain the corresponding exact coefficients of , otherwise return false.

source

# Base.containsMethod.

contains(x::arb_poly, y::fmpq_poly)

Return true if the coefficient balls of contain the corresponding exact coefficients of , otherwise return false.

source

# Base.containsMethod.

contains(x::acb_poly, y::fmpz_poly)

Return true if the coefficient boxes of contain the corresponding exact coefficients of , otherwise return false.

source

# Base.containsMethod.

contains(x::acb_poly, y::fmpq_poly)

Return true if the coefficient boxes of contain the corresponding exact coefficients of , otherwise return false.

source

It is sometimes also useful to be able to determine if there is a unique integer contained in the coefficient of an inexact constant polynomial.

# Nemo.unique_integerMethod.

unique_integer(x::arb_poly)

Return a tuple (t, z) where is true if there is a unique integer contained in each of the coefficients of , otherwise sets to false. In the former case, is set to the integer polynomial.

source

# Nemo.unique_integerMethod.

unique_integer(x::acb_poly)

Return a tuple (t, z) where is true if there is a unique integer contained in the (constant) polynomial , along with that integer in case it is, otherwise sets to false.

source

We also have the following functions.

# Base.isrealMethod.

isreal(x::acb_poly)

Return true if all the coefficients of are real, i.e. have exact zero imaginary parts.

source

Here are some examples of overlapping and containment.

RR = RealField(64)
CC = ComplexField(64)
R, x = PolynomialRing(RR, "x")
C, y = PolynomialRing(CC, "y")
Zx, zx = PolynomialRing(ZZ, "x")
Qx, qx = PolynomialRing(QQ, "x")

f = x^2 + 2x + 1
h = f + RR("0 +/- 0.0001")
k = f + RR("0 +/- 0.0001") * x^4
m = y^2 + 2y + 1
n = m + CC("0 +/- 0.0001", "0 +/- 0.0001")

contains(h, f)
overlaps(f, k)
contains(n, m)
t, z = unique_integer(k)
isreal(n)

Factorisation

Polynomials can only be factorised over certain rings. In general we use the same format for the output as the Julia factorisation function, namely an associative array with polynomial factors as keys and exponents as values.

# Nemo.isirreducibleMethod.

isirreducible(x::nmod_poly)

Return true if is irreducible, otherwise return false.

source

# Nemo.isirreducibleMethod.

isirreducible(x::fmpz_mod_poly)

Return true if is irreducible, otherwise return false.

source

# Nemo.isirreducibleMethod.

isirreducible(x::fq_poly)

Return true if is irreducible, otherwise return false.

source

# Nemo.isirreducibleMethod.

isirreducible(x::fq_nmod_poly)

Return true if is irreducible, otherwise return false.

source

# Nemo.issquarefreeMethod.

issquarefree(x::nmod_poly)

Return true if is squarefree, otherwise return false.

source

# Nemo.issquarefreeMethod.

issquarefree(x::fmpz_mod_poly)

Return true if is squarefree, otherwise return false.

source

# Nemo.issquarefreeMethod.

issquarefree(x::fq_poly)

Return true if is squarefree, otherwise return false.

source

# Nemo.issquarefreeMethod.

issquarefree(x::fq_nmod_poly)

Return true if is squarefree, otherwise return false.

source

# Base.factorMethod.

factor(x::fmpz_poly)

Returns the factorization of .

source

# Base.factorMethod.

factor(x::nmod_poly)

Return the factorisation of .

source

# Base.factorMethod.

factor(x::fmpz_mod_poly)

Return the factorisation of .

source

# Base.factorMethod.

factor(x::fq_poly)

Return the factorisation of .

source

# Base.factorMethod.

factor(x::fq_nmod_poly)

Return the factorisation of .

source

# Nemo.factor_squarefreeMethod.

factor_squarefree(x::nmod_poly)

Return the squarefree factorisation of .

source

# Nemo.factor_squarefreeMethod.

factor_squarefree(x::fmpz_mod_poly)

Return the squarefree factorisation of .

source

# Nemo.factor_squarefreeMethod.

factor_squarefree(x::fq_poly)

Return the squarefree factorisation of .

source

# Nemo.factor_squarefreeMethod.

factor_squarefree(x::fq_nmod_poly)

Return the squarefree factorisation of .

source

# Nemo.factor_distinct_degMethod.

factor_distinct_deg(x::nmod_poly)

Return the distinct degree factorisation of a squarefree polynomial .

source

# Nemo.factor_distinct_degMethod.

factor_distinct_deg(x::fmpz_mod_poly)

Return the distinct degree factorisation of a squarefree polynomial .

source

# Nemo.factor_distinct_degMethod.

factor_distinct_deg(x::fq_poly)

Return the distinct degree factorisation of a squarefree polynomial .

source

# Nemo.factor_distinct_degMethod.

factor_distinct_deg(x::fq_nmod_poly)

Return the distinct degree factorisation of a squarefree polynomial .

source

Here are some examples of factorisation.

R = ResidueRing(ZZ, 23)
S, x = PolynomialRing(R, "x")

f = x^2 + 2x + 1
g = x^3 + 3x + 1

R = factor(f*g)
S = factor_squarefree(f*g)
T = factor_distinct_deg((x + 1)*g*(x^5+x^3+x+1))

Special functions

The following special functions can be computed for any polynomial ring. Typically one uses the generator of a polynomial ring to get the respective special polynomials expressed in terms of that generator.

# Nemo.chebyshev_tMethod.

chebyshev_t(n::Int, x::PolyElem)

Return the Chebyshev polynomial of the first kind , defined by .

source

# Nemo.chebyshev_uMethod.

chebyshev_u(n::Int, x::PolyElem)

Return the Chebyshev polynomial of the first kind , defined by .

source

The following special polynomials are only available for certain base rings.

# Nemo.cyclotomicMethod.

cyclotomic(n::Int, x::fmpz_poly)

Return the th cyclotomic polynomial, defined as where runs over all the th primitive roots of unity.

source

# Nemo.swinnerton_dyerMethod.

swinnerton_dyer(n::Int, x::fmpz_poly)

Return the Swinnerton-Dyer polynomial , defined as the integer polynomial where denotes the -th prime number and all combinations of signs are taken. This polynomial has degree and is irreducible over the integers (it is the minimal polynomial of ).

source

# Nemo.cos_minpolyMethod.

cos_minpoly(n::Int, x::fmpz_poly)

Return the minimal polynomial of . For suitable choice of , this gives the minimal polynomial of or for any rational .

source

# Nemo.theta_qexpMethod.

theta_qexp(e::Int, n::Int, x::fmpz_poly)

Return the -expansion to length of the Jacobi theta function raised to the power , i.e. where .

source

# Nemo.eta_qexpMethod.

eta_qexp(e::Int, n::Int, x::fmpz_poly)

Return the -expansion to length of the Dedekind eta function (without the leading factor ) raised to the power , i.e. . In particular, gives the generating function of the partition function , and gives, after multiplication by , the modular discriminant which generates the Ramanujan tau function .

source

Here are some examples of special functions.

R, x = PolynomialRing(ZZ, "x")
S, y = PolynomialRing(R, "y")

f = chebyshev_t(20, y)
g = chebyshev_u(15, y)
h = cyclotomic(120, x)
j = swinnerton_dyer(5, x)
k = cos_minpoly(30, x)
l = theta_qexp(3, 30, x)
m = eta_qexp(24, 30, x)
o = cyclotomic(10, 1 + x + x^2)