PolyDict engine for generic multivariate polynomial rings

This module provides an implementation of the underlying arithmetic for multi-variate polynomial rings using Python dicts.

This class is not meant for end users, but instead for implementing multivariate polynomial rings over a completely general base. It does not do strong type checking or have parents, etc. For speed, it has been implemented in Cython.

The functions in this file use the ‘dictionary representation’ of multivariate polynomials

{(e1,...,er):c1,...} <-> c1*x1^e1*...*xr^er+...,

which we call a polydict. The exponent tuple (e1,...,er) in this representation is an instance of the class ETuple. This class behaves like a normal Python tuple but also offers advanced access methods for sparse monomials like positions of non-zero exponents etc.

AUTHORS:

  • William Stein
  • David Joyner
  • Martin Albrecht (ETuple)
  • Joel B. Mohler (2008-03-17) – ETuple rewrite as sparse C array
class sage.rings.polynomial.polydict.ETuple

Bases: object

Representation of the exponents of a polydict monomial. If (0,0,3,0,5) is the exponent tuple of x_2^3*x_4^5 then this class only stores {2:3,4:5} instead of the full tuple. This sparse information may be obtained by provided methods.

The index/value data is all stored in the _data C int array member variable. For the example above, the C array would contain 2,3,4,5. The indices are interlaced with the values.

This data structure is very nice to work with for some functions implemented in this class, but tricky for others. One reason that I really like the format is that it requires a single memory allocation for all of the values. A hash table would require more allocations and presumably be slower. I didn’t benchmark this question (although, there is no question that this is much faster than the prior use of python dicts).

combine_to_positives(other)

Given a pair of ETuples (self, other), returns a triple of ETuples (a, b, c) so that self = a + b,other = a + c and b and c have all positive entries.

EXAMPLES:

sage: from sage.rings.polynomial.polydict import ETuple
sage: e = ETuple([-2,1,-5, 3, 1,0])
sage: f = ETuple([1,-3,-3,4,0,2])
sage: e.combine_to_positives(f)
((-2, -3, -5, 3, 0, 0), (0, 4, 0, 0, 1, 0), (3, 0, 2, 1, 0, 2))
common_nonzero_positions(other, sort=False)

Returns an optionally sorted list of non zero positions either in self or other, i.e. the only positions that need to be considered for any vector operation.

EXAMPLES:

sage: from sage.rings.polynomial.polydict import ETuple
sage: e = ETuple([1,0,2])
sage: f = ETuple([0,0,1])
sage: e.common_nonzero_positions(f)
set([0, 2])
sage: e.common_nonzero_positions(f,sort=True)
[0, 2]
eadd(other)

Vector addition of self with other.

EXAMPLES:

sage: from sage.rings.polynomial.polydict import ETuple
sage: e = ETuple([1,0,2])
sage: f = ETuple([0,1,1])
sage: e.eadd(f)
(1, 1, 3)

Verify that trac 6428 has been addressed:

sage: R.<y,z> = Frac(QQ['x'])[]
sage: type(y)                 
<class 'sage.rings.polynomial.multi_polynomial_element.MPolynomial_polydict'>
sage: y^(2^32)
...
OverflowError: Exponent overflow (2147483648).
eadd_p(other, pos)

Adds other to self at position pos.

EXAMPLES:

sage: from sage.rings.polynomial.polydict import ETuple
sage: e = ETuple([1,0,2])
sage: e.eadd_p(5, 1)
(1, 5, 2)
sage: e = ETuple([0]*7)
sage: e.eadd_p(5,4)
(0, 0, 0, 0, 5, 0, 0)
emax(other)

Vector subtraction of self with other.

EXAMPLES:

sage: from sage.rings.polynomial.polydict import ETuple
sage: e = ETuple([1,0,2])
sage: f = ETuple([0,1,1])
sage: e.emax(f)
(1, 1, 2)
sage: e=ETuple((1,2,3,4))
sage: f=ETuple((4,0,2,1))
sage: f.emax(e)
(4, 2, 3, 4)
sage: e=ETuple((1,-2,-2,4))
sage: f=ETuple((4,0,0,0))
sage: f.emax(e)
(4, 0, 0, 4)
sage: f.emax(e).nonzero_positions()
[0, 3]
emin(other)

Vector subtraction of self with other.

EXAMPLES:

sage: from sage.rings.polynomial.polydict import ETuple
sage: e = ETuple([1,0,2])
sage: f = ETuple([0,1,1])
sage: e.emin(f)
(0, 0, 1)
sage: e = ETuple([1,0,-1])
sage: f = ETuple([0,-2,1])
sage: e.emin(f)
(0, -2, -1)
emul(factor)

Scalar Vector multiplication of self.

EXAMPLES:

sage: from sage.rings.polynomial.polydict import ETuple
sage: e = ETuple([1,0,2])
sage: e.emul(2)
(2, 0, 4)
esub(other)

Vector subtraction of self with other.

EXAMPLES:

sage: from sage.rings.polynomial.polydict import ETuple
sage: e = ETuple([1,0,2])
sage: f = ETuple([0,1,1])
sage: e.esub(f)
(1, -1, 1)
nonzero_positions(sort=False)

Returns the positions of non-zero exponents in the tuple.

INPUT:

  • sort – if True a sorted list is returned. If False an unsorted list is returned. (default: False)

EXAMPLES:

sage: from sage.rings.polynomial.polydict import ETuple
sage: e = ETuple([1,0,2])
sage: e.nonzero_positions()
[0, 2]
nonzero_values(sort=True)

Returns the non-zero values of the tuple.

INPUT:

  • sort – if True the values are sorted by their indices. Otherwise the values are returned unsorted. (default: True)

EXAMPLES:

sage: from sage.rings.polynomial.polydict import ETuple        
sage: e = ETuple([2,0,1])
sage: e.nonzero_values()
[2, 1]
sage: f = ETuple([0,-1,1])
sage: f.nonzero_values(sort=True)
[-1, 1]
reversed()

Returns the reversed ETuple of self.

EXAMPLES:

sage: from sage.rings.polynomial.polydict import ETuple
sage: e = ETuple([1,2,3])
sage: e.reversed()
(3, 2, 1)
sparse_iter()

Iterator over the elements of self where the elements are returned as (i,e) where i is the position of e in the tuple.

EXAMPLES:

sage: from sage.rings.polynomial.polydict import ETuple
sage: e = ETuple([1,0,2,0,3])
sage: list(e.sparse_iter())
[(0, 1), (2, 2), (4, 3)]
class sage.rings.polynomial.polydict.ETupleIter

Bases: object

next()
x.next() -> the next value, or raise StopIteration
class sage.rings.polynomial.polydict.PolyDict

Bases: object

coefficient(mon)

Return a polydict that defines a polynomial in 1 less number of variables that gives the coefficient of mon in this polynomial.

The coefficient is defined as follows. If f is this polynomial, then the coefficient is the sum T/mon where the sum is over terms T in f that are exactly divisible by mon.

coefficients()

Return the coefficients of self.

EXAMPLES:

sage: from sage.rings.polynomial.polydict import PolyDict        
sage: f = PolyDict({(2,3):2, (1,2):3, (2,1):4})
sage: f.coefficients()
[3, 2, 4]
compare(other, fn=None)
degree(x=None)
dict()

Return a copy of the dict that defines self. It is safe to change this. For a reference, use dictref.

EXAMPLES:

sage: from sage.rings.polynomial.polydict import PolyDict        
sage: f = PolyDict({(2,3):2, (1,2):3, (2,1):4})
sage: f.dict()
{(1, 2): 3, (2, 3): 2, (2, 1): 4}
exponents()

Return the exponents of self.

EXAMPLES:

sage: from sage.rings.polynomial.polydict import PolyDict        
sage: f = PolyDict({(2,3):2, (1,2):3, (2,1):4})
sage: f.exponents()
[(1, 2), (2, 3), (2, 1)]
homogenize(var)
is_homogeneous()
latex(vars, atomic_exponents=True, atomic_coefficients=True, cmpfn=None)

Return a nice polynomial latex representation of this PolyDict, where the vars are substituted in.

INPUT:

  • vars – list
  • atomic_exponents – bool (default: True)
  • atomic_coefficients – bool (default: True)

EXAMPLES:

sage: from sage.rings.polynomial.polydict import PolyDict        
sage: f = PolyDict({(2,3):2, (1,2):3, (2,1):4})
sage: f.latex(['a','WW'])
'2 a^{2}WW^{3} + 4 a^{2}WW + 3 aWW^{2}'

When atomic_exponents is False, the exponents are surrounded in parenthesis, since ^ has such high precedence:

# I've removed fractional exponent support in ETuple when moving to a sparse C integer array
#sage: f = PolyDict({(2/3,3,5):2, (1,2,1):3, (2,1,1):4}, force_int_exponents=False)        
#sage: f.latex(['a','b','c'], atomic_exponents=False)
#'4 a^{2}bc + 3 ab^{2}c + 2 a^{2/3}b^{3}c^{5}'
lcmt(greater_etuple)

Provides functionality of lc, lm, and lt by calling the tuple compare function on the provided term order T.

INPUT:

  • greater_etuple – a term order
list()

Return a list that defines self. It is safe to change this.

EXAMPLES:

sage: from sage.rings.polynomial.polydict import PolyDict
sage: f = PolyDict({(2,3):2, (1,2):3, (2,1):4})
sage: f.list()
[[3, [1, 2]], [2, [2, 3]], [4, [2, 1]]]
max_exp()

Returns an ETuple containing the maximum exponents appearing. If there are no terms at all in the PolyDict, it returns None.

The nvars parameter is necessary because a PolyDict doesn’t know it from the data it has (and an empty PolyDict offers no clues).

EXAMPLES:

sage: from sage.rings.polynomial.polydict import PolyDict
sage: f = PolyDict({(2,3):2, (1,2):3, (2,1):4})
sage: f.max_exp()
(2, 3)
sage: PolyDict({}).max_exp() # returns None
min_exp()

Returns an ETuple containing the minimum exponents appearing. If there are no terms at all in the PolyDict, it returns None.

The nvars parameter is necessary because a PolyDict doesn’t know it from the data it has (and an empty PolyDict offers no clues).

EXAMPLES:

sage: from sage.rings.polynomial.polydict import PolyDict
sage: f = PolyDict({(2,3):2, (1,2):3, (2,1):4})
sage: f.min_exp()
(1, 1)
sage: PolyDict({}).min_exp() # returns None
monomial_coefficient(mon)

EXAMPLES:

sage: from sage.rings.polynomial.polydict import PolyDict
sage: f = PolyDict({(2,3):2, (1,2):3, (2,1):4})
sage: f.monomial_coefficient(PolyDict({(2,1):1}).dict())
4
poly_repr(vars, atomic_exponents=True, atomic_coefficients=True, cmpfn=None)

Return a nice polynomial string representation of this PolyDict, where the vars are substituted in.

INPUT:

  • vars – list
  • atomic_exponents – bool (default: True)
  • atomic_coefficients – bool (default: True)

EXAMPLES:

sage: from sage.rings.polynomial.polydict import PolyDict
sage: f = PolyDict({(2,3):2, (1,2):3, (2,1):4})
sage: f.poly_repr(['a','WW'])
'2*a^2*WW^3 + 4*a^2*WW + 3*a*WW^2'

When atomic_exponents is False, the exponents are surrounded in parenthesis, since ^ has such high precedence.

# I've removed fractional exponent support in ETuple when moving to a sparse C integer array
#sage: f = PolyDict({(2/3,3,5):2, (1,2,1):3, (2,1,1):4}, force_int_exponents=False)        
#sage: f.poly_repr(['a','b','c'], atomic_exponents=False)
#'4*a^(2)*b*c + 3*a*b^(2)*c + 2*a^(2/3)*b^(3)*c^(5)'

We check to make sure that when we are in characteristic two, we don’t put negative signs on the generators.

sage: Integers(2)['x,y'].gens()
(x, y)
polynomial_coefficient(degrees)

Return a polydict that defines the coefficient in the current polynomial viewed as a tower of polynomial extensions.

INPUT:

  • degrees – a list of degree restrictions; list elements are None if the variable in that position should be unrestricted

EXAMPLES:

sage: from sage.rings.polynomial.polydict import PolyDict
sage: f = PolyDict({(2,3):2, (1,2):3, (2,1):4})
sage: f.polynomial_coefficient([2,None])
PolyDict with representation {(0, 3): 2, (0, 1): 4}
sage: f = PolyDict({(0,3):2, (0,2):3, (2,1):4})
sage: f.polynomial_coefficient([0,None])
PolyDict with representation {(0, 3): 2, (0, 2): 3}
scalar_lmult(s)

Left Scalar Multiplication

EXAMPLES:

sage: from sage.rings.polynomial.polydict import PolyDict
sage: x,y=FreeMonoid(2,'x,y').gens()  # a strange object to live in a polydict, but non-commutative!
sage: f = PolyDict({(2,3):x})
sage: f.scalar_lmult(y)
PolyDict with representation {(2, 3): y*x}
sage: f = PolyDict({(2,3):2, (1,2):3, (2,1):4})
sage: f.scalar_lmult(-2)
PolyDict with representation {(1, 2): -6, (2, 3): -4, (2, 1): -8}
scalar_rmult(s)

Right Scalar Multiplication

EXAMPLES:

sage: from sage.rings.polynomial.polydict import PolyDict
sage: x,y=FreeMonoid(2,'x,y').gens()  # a strange object to live in a polydict, but non-commutative!
sage: f = PolyDict({(2,3):x})
sage: f.scalar_rmult(y)
PolyDict with representation {(2, 3): x*y}
sage: f = PolyDict({(2,3):2, (1,2):3, (2,1):4})
sage: f.scalar_rmult(-2)
PolyDict with representation {(1, 2): -6, (2, 3): -4, (2, 1): -8}
total_degree()
valuation(x=None)
sage.rings.polynomial.polydict.make_ETuple(data, length)
sage.rings.polynomial.polydict.make_PolyDict(data)

Previous topic

Direct low-level access to SINGULAR’s Groebner basis engine via libSINGULAR.

Next topic

Infinite Polynomial Rings

This Page