Symmetric Reduction of Infinite Polynomials

SymmetricReductionStrategy provides a framework for efficient symmetric reduction of Infinite Polynomials, see infinite_polynomial_element.

AUTHORS:

THEORY:

According to M. Aschenbrenner and C. Hillar [AB2007], Symmetric Reduction of an element \(p\) of an Infinite Polynomial Ring \(X\) by some other element \(q\) means the following:

  1. Let \(M\) and \(N\) be the leading terms of \(p\) and \(q\).

  2. Test whether there is a permutation \(P\) that does not diminish the variable indices occurring in \(N\) and preserves their order, so that there is some term \(T\in X\) with \(T N^P = M\). If there is no such permutation, return \(p\).

  3. Replace \(p\) by \(p-T q^P\) and continue with step 1.

When reducing one polynomial \(p\) with respect to a list \(L\) of other polynomials, there usually is a choice of order on which the efficiency crucially depends. Also it helps to modify the polynomials on the list in order to simplify the basic reduction steps.

The preparation of \(L\) may be expensive. Hence, if the same list is used many times then it is reasonable to perform the preparation only once. This is the background of SymmetricReductionStrategy.

Our current strategy is to keep the number of terms in the polynomials as small as possible. For this, we sort \(L\) by increasing number of terms. If several elements of \(L\) allow for a reduction of \(p\), we choose the one with the smallest number of terms. Later on, it should be possible to implement further strategies for choice.

When adding a new polynomial \(q\) to \(L\), we first reduce \(q\) with respect to \(L\). Then, we test heuristically whether it is possible to reduce the number of terms of the elements of \(L\) by reduction modulo \(q\). That way, we see best chances to keep the number of terms in intermediate reduction steps relatively small.

EXAMPLES:

First, we create an infinite polynomial ring and one of its elements:

sage: X.<x,y> = InfinitePolynomialRing(QQ)
sage: p = y[1]*y[3] + y[1]^2*x[3]
>>> from sage.all import *
>>> X = InfinitePolynomialRing(QQ, names=('x', 'y',)); (x, y,) = X._first_ngens(2)
>>> p = y[Integer(1)]*y[Integer(3)] + y[Integer(1)]**Integer(2)*x[Integer(3)]
X.<x,y> = InfinitePolynomialRing(QQ)
p = y[1]*y[3] + y[1]^2*x[3]

We want to symmetrically reduce it by another polynomial. So, we put this other polynomial into a list and create a Symmetric Reduction Strategy object:

sage: from sage.rings.polynomial.symmetric_reduction import SymmetricReductionStrategy
sage: S = SymmetricReductionStrategy(X, [y[2]^2*x[1]])
sage: S
Symmetric Reduction Strategy in
 Infinite polynomial ring in x, y over Rational Field, modulo
    x_1*y_2^2
sage: S.reduce(p)
x_3*y_1^2 + y_3*y_1
>>> from sage.all import *
>>> from sage.rings.polynomial.symmetric_reduction import SymmetricReductionStrategy
>>> S = SymmetricReductionStrategy(X, [y[Integer(2)]**Integer(2)*x[Integer(1)]])
>>> S
Symmetric Reduction Strategy in
 Infinite polynomial ring in x, y over Rational Field, modulo
    x_1*y_2^2
>>> S.reduce(p)
x_3*y_1^2 + y_3*y_1
from sage.rings.polynomial.symmetric_reduction import SymmetricReductionStrategy
S = SymmetricReductionStrategy(X, [y[2]^2*x[1]])
S
S.reduce(p)

The preceding is correct, since any permutation that turns y[2]^2*x[1] into a factor of y[1]^2*x[3] interchanges the variable indices 1 and 2 – which is not allowed in a symmetric reduction. However, reduction by y[1]^2*x[2] works, since one can change variable index 1 into 2 and 2 into 3. So, we add this to S:

sage: S.add_generator(y[1]^2*x[2])
sage: S
Symmetric Reduction Strategy in
 Infinite polynomial ring in x, y over Rational Field, modulo
    x_2*y_1^2,
    x_1*y_2^2
sage: S.reduce(p)                                                                   # needs sage.combinat
y_3*y_1
>>> from sage.all import *
>>> S.add_generator(y[Integer(1)]**Integer(2)*x[Integer(2)])
>>> S
Symmetric Reduction Strategy in
 Infinite polynomial ring in x, y over Rational Field, modulo
    x_2*y_1^2,
    x_1*y_2^2
>>> S.reduce(p)                                                                   # needs sage.combinat
y_3*y_1
S.add_generator(y[1]^2*x[2])
S
S.reduce(p)                                                                   # needs sage.combinat

The next example shows that tail reduction is not done, unless it is explicitly advised:

sage: S.reduce(x[3] + 2*x[2]*y[1]^2 + 3*y[2]^2*x[1])                                # needs sage.combinat
x_3 + 2*x_2*y_1^2 + 3*x_1*y_2^2
sage: S.tailreduce(x[3] + 2*x[2]*y[1]^2 + 3*y[2]^2*x[1])                            # needs sage.combinat
x_3
>>> from sage.all import *
>>> S.reduce(x[Integer(3)] + Integer(2)*x[Integer(2)]*y[Integer(1)]**Integer(2) + Integer(3)*y[Integer(2)]**Integer(2)*x[Integer(1)])                                # needs sage.combinat
x_3 + 2*x_2*y_1^2 + 3*x_1*y_2^2
>>> S.tailreduce(x[Integer(3)] + Integer(2)*x[Integer(2)]*y[Integer(1)]**Integer(2) + Integer(3)*y[Integer(2)]**Integer(2)*x[Integer(1)])                            # needs sage.combinat
x_3
S.reduce(x[3] + 2*x[2]*y[1]^2 + 3*y[2]^2*x[1])                                # needs sage.combinat
S.tailreduce(x[3] + 2*x[2]*y[1]^2 + 3*y[2]^2*x[1])                            # needs sage.combinat

However, it is possible to ask for tailreduction already when the Symmetric Reduction Strategy is created:

sage: S2 = SymmetricReductionStrategy(X, [y[2]^2*x[1],y[1]^2*x[2]], tailreduce=True)
sage: S2
Symmetric Reduction Strategy in
 Infinite polynomial ring in x, y over Rational Field, modulo
    x_2*y_1^2,
    x_1*y_2^2
with tailreduction
sage: S2.reduce(x[3] + 2*x[2]*y[1]^2 + 3*y[2]^2*x[1])                               # needs sage.combinat
x_3
>>> from sage.all import *
>>> S2 = SymmetricReductionStrategy(X, [y[Integer(2)]**Integer(2)*x[Integer(1)],y[Integer(1)]**Integer(2)*x[Integer(2)]], tailreduce=True)
>>> S2
Symmetric Reduction Strategy in
 Infinite polynomial ring in x, y over Rational Field, modulo
    x_2*y_1^2,
    x_1*y_2^2
with tailreduction
>>> S2.reduce(x[Integer(3)] + Integer(2)*x[Integer(2)]*y[Integer(1)]**Integer(2) + Integer(3)*y[Integer(2)]**Integer(2)*x[Integer(1)])                               # needs sage.combinat
x_3
S2 = SymmetricReductionStrategy(X, [y[2]^2*x[1],y[1]^2*x[2]], tailreduce=True)
S2
S2.reduce(x[3] + 2*x[2]*y[1]^2 + 3*y[2]^2*x[1])                               # needs sage.combinat
class sage.rings.polynomial.symmetric_reduction.SymmetricReductionStrategy[source]

Bases: object

A framework for efficient symmetric reduction of InfinitePolynomial, see infinite_polynomial_element.

INPUT:

  • Parent – an Infinite Polynomial Ring, see infinite_polynomial_element

  • L – list (default: the empty list); list of elements of Parent with respect to which will be reduced

  • good_input – boolean (default: None); if this optional parameter is true, it is assumed that each element of L is symmetrically reduced with respect to the previous elements of L

EXAMPLES:

sage: X.<y> = InfinitePolynomialRing(QQ)
sage: from sage.rings.polynomial.symmetric_reduction import SymmetricReductionStrategy
sage: S = SymmetricReductionStrategy(X, [y[2]^2*y[1],y[1]^2*y[2]], good_input=True)
sage: S.reduce(y[3] + 2*y[2]*y[1]^2 + 3*y[2]^2*y[1])
y_3 + 3*y_2^2*y_1 + 2*y_2*y_1^2
sage: S.tailreduce(y[3] + 2*y[2]*y[1]^2 + 3*y[2]^2*y[1])                        # needs sage.combinat
y_3
>>> from sage.all import *
>>> X = InfinitePolynomialRing(QQ, names=('y',)); (y,) = X._first_ngens(1)
>>> from sage.rings.polynomial.symmetric_reduction import SymmetricReductionStrategy
>>> S = SymmetricReductionStrategy(X, [y[Integer(2)]**Integer(2)*y[Integer(1)],y[Integer(1)]**Integer(2)*y[Integer(2)]], good_input=True)
>>> S.reduce(y[Integer(3)] + Integer(2)*y[Integer(2)]*y[Integer(1)]**Integer(2) + Integer(3)*y[Integer(2)]**Integer(2)*y[Integer(1)])
y_3 + 3*y_2^2*y_1 + 2*y_2*y_1^2
>>> S.tailreduce(y[Integer(3)] + Integer(2)*y[Integer(2)]*y[Integer(1)]**Integer(2) + Integer(3)*y[Integer(2)]**Integer(2)*y[Integer(1)])                        # needs sage.combinat
y_3
X.<y> = InfinitePolynomialRing(QQ)
from sage.rings.polynomial.symmetric_reduction import SymmetricReductionStrategy
S = SymmetricReductionStrategy(X, [y[2]^2*y[1],y[1]^2*y[2]], good_input=True)
S.reduce(y[3] + 2*y[2]*y[1]^2 + 3*y[2]^2*y[1])
S.tailreduce(y[3] + 2*y[2]*y[1]^2 + 3*y[2]^2*y[1])                        # needs sage.combinat
add_generator(p, good_input=None)[source]

Add another polynomial to self.

INPUT:

  • p – an element of the underlying infinite polynomial ring

  • good_input – boolean (default: None); if True, it is assumed that p is reduced with respect to self. Otherwise, this reduction will be done first (which may cost some time).

Note

Previously added polynomials may be modified. All input is prepared in view of an efficient symmetric reduction.

EXAMPLES:

sage: from sage.rings.polynomial.symmetric_reduction import SymmetricReductionStrategy
sage: X.<x,y> = InfinitePolynomialRing(QQ)
sage: S = SymmetricReductionStrategy(X)
sage: S
Symmetric Reduction Strategy in
 Infinite polynomial ring in x, y over Rational Field
sage: S.add_generator(y[3] + y[1]*(x[3]+x[1]))
sage: S
Symmetric Reduction Strategy in
 Infinite polynomial ring in x, y over Rational Field, modulo
    x_3*y_1 + x_1*y_1 + y_3
>>> from sage.all import *
>>> from sage.rings.polynomial.symmetric_reduction import SymmetricReductionStrategy
>>> X = InfinitePolynomialRing(QQ, names=('x', 'y',)); (x, y,) = X._first_ngens(2)
>>> S = SymmetricReductionStrategy(X)
>>> S
Symmetric Reduction Strategy in
 Infinite polynomial ring in x, y over Rational Field
>>> S.add_generator(y[Integer(3)] + y[Integer(1)]*(x[Integer(3)]+x[Integer(1)]))
>>> S
Symmetric Reduction Strategy in
 Infinite polynomial ring in x, y over Rational Field, modulo
    x_3*y_1 + x_1*y_1 + y_3
from sage.rings.polynomial.symmetric_reduction import SymmetricReductionStrategy
X.<x,y> = InfinitePolynomialRing(QQ)
S = SymmetricReductionStrategy(X)
S
S.add_generator(y[3] + y[1]*(x[3]+x[1]))
S

Note that the first added polynomial will be simplified when adding a suitable second polynomial:

sage: S.add_generator(x[2] + x[1])                                          # needs sage.combinat
sage: S                                                                     # needs sage.combinat
Symmetric Reduction Strategy in
 Infinite polynomial ring in x, y over Rational Field, modulo
    y_3,
    x_2 + x_1
>>> from sage.all import *
>>> S.add_generator(x[Integer(2)] + x[Integer(1)])                                          # needs sage.combinat
>>> S                                                                     # needs sage.combinat
Symmetric Reduction Strategy in
 Infinite polynomial ring in x, y over Rational Field, modulo
    y_3,
    x_2 + x_1
S.add_generator(x[2] + x[1])                                          # needs sage.combinat
S                                                                     # needs sage.combinat

By default, reduction is applied to any newly added polynomial. This can be avoided by specifying the optional parameter ‘good_input’:

sage: # needs sage.combinat
sage: S.add_generator(y[2] + y[1]*x[2])
sage: S
Symmetric Reduction Strategy in
 Infinite polynomial ring in x, y over Rational Field, modulo
    y_3,
    x_1*y_1 - y_2,
    x_2 + x_1
sage: S.reduce(x[3] + x[2])
-2*x_1
sage: S.add_generator(x[3] + x[2], good_input=True)
sage: S
Symmetric Reduction Strategy in
 Infinite polynomial ring in x, y over Rational Field, modulo
    y_3,
    x_3 + x_2,
    x_1*y_1 - y_2,
    x_2 + x_1
>>> from sage.all import *
>>> # needs sage.combinat
>>> S.add_generator(y[Integer(2)] + y[Integer(1)]*x[Integer(2)])
>>> S
Symmetric Reduction Strategy in
 Infinite polynomial ring in x, y over Rational Field, modulo
    y_3,
    x_1*y_1 - y_2,
    x_2 + x_1
>>> S.reduce(x[Integer(3)] + x[Integer(2)])
-2*x_1
>>> S.add_generator(x[Integer(3)] + x[Integer(2)], good_input=True)
>>> S
Symmetric Reduction Strategy in
 Infinite polynomial ring in x, y over Rational Field, modulo
    y_3,
    x_3 + x_2,
    x_1*y_1 - y_2,
    x_2 + x_1
# needs sage.combinat
S.add_generator(y[2] + y[1]*x[2])
S
S.reduce(x[3] + x[2])
S.add_generator(x[3] + x[2], good_input=True)
S

In the previous example, x[3] + x[2] is added without being reduced to zero.

gens()[source]

Return the list of Infinite Polynomials modulo which self reduces.

EXAMPLES:

sage: X.<y> = InfinitePolynomialRing(QQ)
sage: from sage.rings.polynomial.symmetric_reduction import SymmetricReductionStrategy
sage: S = SymmetricReductionStrategy(X, [y[2]^2*y[1],y[1]^2*y[2]])
sage: S
Symmetric Reduction Strategy in
 Infinite polynomial ring in y over Rational Field, modulo
    y_2*y_1^2,
    y_2^2*y_1
sage: S.gens()
[y_2*y_1^2, y_2^2*y_1]
>>> from sage.all import *
>>> X = InfinitePolynomialRing(QQ, names=('y',)); (y,) = X._first_ngens(1)
>>> from sage.rings.polynomial.symmetric_reduction import SymmetricReductionStrategy
>>> S = SymmetricReductionStrategy(X, [y[Integer(2)]**Integer(2)*y[Integer(1)],y[Integer(1)]**Integer(2)*y[Integer(2)]])
>>> S
Symmetric Reduction Strategy in
 Infinite polynomial ring in y over Rational Field, modulo
    y_2*y_1^2,
    y_2^2*y_1
>>> S.gens()
[y_2*y_1^2, y_2^2*y_1]
X.<y> = InfinitePolynomialRing(QQ)
from sage.rings.polynomial.symmetric_reduction import SymmetricReductionStrategy
S = SymmetricReductionStrategy(X, [y[2]^2*y[1],y[1]^2*y[2]])
S
S.gens()
reduce(p, notail=False, report=None)[source]

Symmetric reduction of an infinite polynomial.

INPUT:

  • p – an element of the underlying infinite polynomial ring

  • notail – boolean (default: False); if True, tail reduction is avoided (but there is no guarantee that there will be no tail reduction at all)

  • report – object (default: None); if not None, print information on the progress of the computation

OUTPUT: reduction of p with respect to self

Note

If tail reduction shall be forced, use tailreduce().

EXAMPLES:

sage: from sage.rings.polynomial.symmetric_reduction import SymmetricReductionStrategy
sage: X.<x,y> = InfinitePolynomialRing(QQ)
sage: S = SymmetricReductionStrategy(X, [y[3]], tailreduce=True)
sage: S.reduce(y[4]*x[1] + y[1]*x[4])
x_4*y_1
sage: S.reduce(y[4]*x[1] + y[1]*x[4], notail=True)
x_4*y_1 + x_1*y_4
>>> from sage.all import *
>>> from sage.rings.polynomial.symmetric_reduction import SymmetricReductionStrategy
>>> X = InfinitePolynomialRing(QQ, names=('x', 'y',)); (x, y,) = X._first_ngens(2)
>>> S = SymmetricReductionStrategy(X, [y[Integer(3)]], tailreduce=True)
>>> S.reduce(y[Integer(4)]*x[Integer(1)] + y[Integer(1)]*x[Integer(4)])
x_4*y_1
>>> S.reduce(y[Integer(4)]*x[Integer(1)] + y[Integer(1)]*x[Integer(4)], notail=True)
x_4*y_1 + x_1*y_4
from sage.rings.polynomial.symmetric_reduction import SymmetricReductionStrategy
X.<x,y> = InfinitePolynomialRing(QQ)
S = SymmetricReductionStrategy(X, [y[3]], tailreduce=True)
S.reduce(y[4]*x[1] + y[1]*x[4])
S.reduce(y[4]*x[1] + y[1]*x[4], notail=True)

Last, we demonstrate the report option:

sage: S = SymmetricReductionStrategy(X, [x[2] + y[1],
....:                                    x[2]*y[3] + x[1]*y[2] + y[4],
....:                                    y[3] + y[2]])
sage: S
Symmetric Reduction Strategy in
 Infinite polynomial ring in x, y over Rational Field, modulo
    y_3 + y_2,
    x_2 + y_1,
    x_1*y_2 + y_4 - y_3*y_1
sage: S.reduce(x[3] + x[1]*y[3] + x[1]*y[1], report=True)
:::>
x_1*y_1 + y_4 - y_3*y_1 - y_1
>>> from sage.all import *
>>> S = SymmetricReductionStrategy(X, [x[Integer(2)] + y[Integer(1)],
...                                    x[Integer(2)]*y[Integer(3)] + x[Integer(1)]*y[Integer(2)] + y[Integer(4)],
...                                    y[Integer(3)] + y[Integer(2)]])
>>> S
Symmetric Reduction Strategy in
 Infinite polynomial ring in x, y over Rational Field, modulo
    y_3 + y_2,
    x_2 + y_1,
    x_1*y_2 + y_4 - y_3*y_1
>>> S.reduce(x[Integer(3)] + x[Integer(1)]*y[Integer(3)] + x[Integer(1)]*y[Integer(1)], report=True)
:::>
x_1*y_1 + y_4 - y_3*y_1 - y_1
S = SymmetricReductionStrategy(X, [x[2] + y[1],
                                   x[2]*y[3] + x[1]*y[2] + y[4],
                                   y[3] + y[2]])
S
S.reduce(x[3] + x[1]*y[3] + x[1]*y[1], report=True)

Each ‘:’ indicates that one reduction of the leading monomial was performed. Eventually, the ‘>’ indicates that the computation is finished.

reset()[source]

Remove all polynomials from self.

EXAMPLES:

sage: X.<y> = InfinitePolynomialRing(QQ)
sage: from sage.rings.polynomial.symmetric_reduction import SymmetricReductionStrategy
sage: S = SymmetricReductionStrategy(X, [y[2]^2*y[1],y[1]^2*y[2]])
sage: S
Symmetric Reduction Strategy in
 Infinite polynomial ring in y over Rational Field, modulo
    y_2*y_1^2,
    y_2^2*y_1
sage: S.reset()
sage: S
Symmetric Reduction Strategy in Infinite polynomial ring in y over Rational Field
>>> from sage.all import *
>>> X = InfinitePolynomialRing(QQ, names=('y',)); (y,) = X._first_ngens(1)
>>> from sage.rings.polynomial.symmetric_reduction import SymmetricReductionStrategy
>>> S = SymmetricReductionStrategy(X, [y[Integer(2)]**Integer(2)*y[Integer(1)],y[Integer(1)]**Integer(2)*y[Integer(2)]])
>>> S
Symmetric Reduction Strategy in
 Infinite polynomial ring in y over Rational Field, modulo
    y_2*y_1^2,
    y_2^2*y_1
>>> S.reset()
>>> S
Symmetric Reduction Strategy in Infinite polynomial ring in y over Rational Field
X.<y> = InfinitePolynomialRing(QQ)
from sage.rings.polynomial.symmetric_reduction import SymmetricReductionStrategy
S = SymmetricReductionStrategy(X, [y[2]^2*y[1],y[1]^2*y[2]])
S
S.reset()
S
setgens(L)[source]

Define the list of Infinite Polynomials modulo which self reduces.

INPUT:

  • L – list of elements of the underlying infinite polynomial ring

Note

It is not tested if L is a good input. That method simply assigns a copy of L to the generators of self.

EXAMPLES:

sage: from sage.rings.polynomial.symmetric_reduction import SymmetricReductionStrategy
sage: X.<y> = InfinitePolynomialRing(QQ)
sage: S = SymmetricReductionStrategy(X, [y[2]^2*y[1],y[1]^2*y[2]])
sage: R = SymmetricReductionStrategy(X)
sage: R.setgens(S.gens())
sage: R
Symmetric Reduction Strategy in
 Infinite polynomial ring in y over Rational Field, modulo
    y_2*y_1^2,
    y_2^2*y_1
sage: R.gens() is S.gens()
False
sage: R.gens() == S.gens()
True
>>> from sage.all import *
>>> from sage.rings.polynomial.symmetric_reduction import SymmetricReductionStrategy
>>> X = InfinitePolynomialRing(QQ, names=('y',)); (y,) = X._first_ngens(1)
>>> S = SymmetricReductionStrategy(X, [y[Integer(2)]**Integer(2)*y[Integer(1)],y[Integer(1)]**Integer(2)*y[Integer(2)]])
>>> R = SymmetricReductionStrategy(X)
>>> R.setgens(S.gens())
>>> R
Symmetric Reduction Strategy in
 Infinite polynomial ring in y over Rational Field, modulo
    y_2*y_1^2,
    y_2^2*y_1
>>> R.gens() is S.gens()
False
>>> R.gens() == S.gens()
True
from sage.rings.polynomial.symmetric_reduction import SymmetricReductionStrategy
X.<y> = InfinitePolynomialRing(QQ)
S = SymmetricReductionStrategy(X, [y[2]^2*y[1],y[1]^2*y[2]])
R = SymmetricReductionStrategy(X)
R.setgens(S.gens())
R
R.gens() is S.gens()
R.gens() == S.gens()
tailreduce(p, report=None)[source]

Symmetric reduction of an infinite polynomial, with forced tail reduction.

INPUT:

  • p – an element of the underlying infinite polynomial ring

  • report – object (default: None); if not None, print information on the progress of the computation

OUTPUT:

Reduction (including the non-leading elements) of p with respect to self.

EXAMPLES:

sage: from sage.rings.polynomial.symmetric_reduction import SymmetricReductionStrategy
sage: X.<x,y> = InfinitePolynomialRing(QQ)
sage: S = SymmetricReductionStrategy(X, [y[3]])
sage: S.reduce(y[4]*x[1] + y[1]*x[4])
x_4*y_1 + x_1*y_4
sage: S.tailreduce(y[4]*x[1] + y[1]*x[4])                                   # needs sage.combinat
x_4*y_1
>>> from sage.all import *
>>> from sage.rings.polynomial.symmetric_reduction import SymmetricReductionStrategy
>>> X = InfinitePolynomialRing(QQ, names=('x', 'y',)); (x, y,) = X._first_ngens(2)
>>> S = SymmetricReductionStrategy(X, [y[Integer(3)]])
>>> S.reduce(y[Integer(4)]*x[Integer(1)] + y[Integer(1)]*x[Integer(4)])
x_4*y_1 + x_1*y_4
>>> S.tailreduce(y[Integer(4)]*x[Integer(1)] + y[Integer(1)]*x[Integer(4)])                                   # needs sage.combinat
x_4*y_1
from sage.rings.polynomial.symmetric_reduction import SymmetricReductionStrategy
X.<x,y> = InfinitePolynomialRing(QQ)
S = SymmetricReductionStrategy(X, [y[3]])
S.reduce(y[4]*x[1] + y[1]*x[4])
S.tailreduce(y[4]*x[1] + y[1]*x[4])                                   # needs sage.combinat

Last, we demonstrate the ‘report’ option:

sage: S = SymmetricReductionStrategy(X, [x[2] + y[1],
....:                                    x[2]*x[3] + x[1]*y[2] + y[4],
....:                                    y[3] + y[2]])
sage: S
Symmetric Reduction Strategy in
 Infinite polynomial ring in x, y over Rational Field, modulo
    y_3 + y_2,
    x_2 + y_1,
    x_1*y_2 + y_4 + y_1^2
sage: S.tailreduce(x[3] + x[1]*y[3] + x[1]*y[1], report=True)               # needs sage.combinat
T[3]:::>
T[3]:>
x_1*y_1 - y_2 + y_1^2 - y_1
>>> from sage.all import *
>>> S = SymmetricReductionStrategy(X, [x[Integer(2)] + y[Integer(1)],
...                                    x[Integer(2)]*x[Integer(3)] + x[Integer(1)]*y[Integer(2)] + y[Integer(4)],
...                                    y[Integer(3)] + y[Integer(2)]])
>>> S
Symmetric Reduction Strategy in
 Infinite polynomial ring in x, y over Rational Field, modulo
    y_3 + y_2,
    x_2 + y_1,
    x_1*y_2 + y_4 + y_1^2
>>> S.tailreduce(x[Integer(3)] + x[Integer(1)]*y[Integer(3)] + x[Integer(1)]*y[Integer(1)], report=True)               # needs sage.combinat
T[3]:::>
T[3]:>
x_1*y_1 - y_2 + y_1^2 - y_1
S = SymmetricReductionStrategy(X, [x[2] + y[1],
                                   x[2]*x[3] + x[1]*y[2] + y[4],
                                   y[3] + y[2]])
S
S.tailreduce(x[3] + x[1]*y[3] + x[1]*y[1], report=True)               # needs sage.combinat
The protocol means the following.
  • ‘T[3]’ means that we currently do tail reduction for a polynomial with three terms.

  • ‘:::>’ means that there were three reductions of leading terms.

  • The tail of the result of the preceding reduction still has three terms. One reduction of leading terms was possible, and then the final result was obtained.