Skew Univariate Polynomial Rings

This module provides the SkewPolynomialRing_general, which constructs a general dense univariate skew polynomials over commutative base rings with automorphisms over the base rings. This is usual accessed only indirectly through the constructor sage.rings.polynomial.skew_polynomial_constructor.SkewPolynomialRing().

See SkewPolynomialRing_general for a definition of a univariate skew polynomial ring.

AUTHOR:

  • Xavier Caruso (2012-06-29): initial version
  • Arpit Merchant (2016-08-04): improved docstrings, fixed doctests and refactored classes and methods
  • Johan Rosenkilde (2016-08-03): changes for bug fixes, docstring and doctest errors
class sage.rings.polynomial.skew_polynomial_ring.SkewPolynomialRing_general(base_ring, twist_map, name, sparse, element_class)

Bases: sage.rings.ring.Algebra, sage.structure.unique_representation.UniqueRepresentation

A general implementation of univariate skew polynomialring over a commutative ring.

Let \(R\) be a commutative ring, and let \(\sigma\) be an automorphism of \(R\). The ring of skew polynomials \(R[X, \sigma]\) is the polynomial ring \(R[X]\), where the addition is the usual polynomial addition, but the multiplication operation is defined by the modified rule

\[X*a = \sigma(a) X.\]

This means that \(R[X, \sigma]\) is a non-commutative ring. Skew polynomials were first introduced by Ore [Ore1933].

EXAMPLES:

sage: R.<t> = ZZ[]
sage: sigma = R.hom([t+1])
sage: S.<x> = SkewPolynomialRing(R,sigma); S
Skew Polynomial Ring in x over Univariate Polynomial Ring in t over Integer Ring
 twisted by t |--> t + 1

One can also use a shorter syntax:

sage: S.<x> = R['x',sigma]; S
Skew Polynomial Ring in x over Univariate Polynomial Ring in t over Integer Ring
 twisted by t |--> t + 1

If we omit the diamond notation, the variable holding the indeterminate is not assigned:

sage: Sy = R['y',sigma]
sage: y
Traceback (most recent call last):
...
NameError: name 'y' is not defined
sage: Sy.gen()
y

Note however that contrary to usual polynomial rings, we cannot omit the variable name on the RHS, since this collides with the notation for creating polynomial rings:

sage: Sz.<z> = R[sigma]
Traceback (most recent call last):
...
ValueError: variable name 'Ring endomorphism of Univariate Polynomial Ring in t over Integer Ring\n
    Defn: t |--> t + 1' is not alphanumeric

Of course, skew polynomial rings with different twist maps are not equal either:

sage: R['x',sigma] == R['x',sigma^2]
False

Saving and loading of polynomial rings works:

sage: loads(dumps(R['x',sigma])) == R['x',sigma]
True

There is a coercion map from the base ring of the skew polynomial rings:

sage: S.has_coerce_map_from(R)
True
sage: x.parent()
Skew Polynomial Ring in x over Univariate Polynomial Ring in t over Integer Ring
 twisted by t |--> t + 1
sage: t.parent()
Univariate Polynomial Ring in t over Integer Ring
sage: y = x+t; y
x + t
sage: y.parent() is S
True
change_var(var)

Return the skew polynomial ring in variable var with the same base ring and twist map as self.

INPUT:

  • var – a string representing the name of the new variable.

OUTPUT:

self with variable name changed to var.

EXAMPLES:

sage: k.<t> = GF(5^3)
sage: Frob = k.frobenius_endomorphism()
sage: R.<x> = SkewPolynomialRing(k,Frob); R
Skew Polynomial Ring in x over Finite Field in t of size 5^3 twisted by t |--> t^5
sage: Ry = R.change_var('y'); Ry
Skew Polynomial Ring in y over Finite Field in t of size 5^3 twisted by t |--> t^5
sage: Ry is R.change_var('y')
True
characteristic()

Return the characteristic of the base ring of self.

EXAMPLES:

sage: R.<t> = QQ[]
sage: sigma = R.hom([t+1])
sage: R['x',sigma].characteristic()
0

sage: k.<u> = GF(5^3)
sage: Frob = k.frobenius_endomorphism()
sage: k['y',Frob].characteristic()
5
gen(n=0)

Return the indeterminate generator of this skew polynomial ring.

INPUT:

  • n – index of generator to return (default: 0). Exists for compatibility with other polynomial rings.

EXAMPLES:

sage: R.<t> = QQ[]
sage: sigma = R.hom([t+1])
sage: S.<x> = R['x',sigma]; S
Skew Polynomial Ring in x over Univariate Polynomial Ring in t over Rational Field twisted by t |--> t + 1
sage: y = S.gen(); y
x
sage: y == x
True
sage: y is x
True
sage: S.gen(0)
x

This is also known as the parameter:

sage: S.parameter() is S.gen()
True
gens_dict()

Return a {name: variable} dictionary of the generators of self.

EXAMPLES:

sage: R.<t> = ZZ[]
sage: sigma = R.hom([t+1])
sage: S.<x> = SkewPolynomialRing(R,sigma)
sage: S.gens_dict()
{'x': x}
is_commutative()

Return True if this skew polynomial ring is commutative, i.e. if the twist map is the identity.

EXAMPLES:

sage: k.<t> = GF(5^3)
sage: Frob = k.frobenius_endomorphism()
sage: S.<x> = k['x',Frob]
sage: S.is_commutative()
False

sage: T.<y> = k['y',Frob^3]
sage: T.is_commutative()
True
is_exact()

Return True if elements of this skew polynomial ring are exact. This happens if and only if elements of the base ring are exact.

EXAMPLES:

sage: k.<t> = GF(5^3)
sage: Frob = k.frobenius_endomorphism()
sage: S.<x> = k['x',Frob]
sage: S.is_exact()
True
sage: S.base_ring().is_exact()
True

sage: R.<u> = k[[]]
sage: sigma = R.hom([u+u^2])
sage: T.<y> = R['y',sigma]
sage: T.is_exact()
False
sage: T.base_ring().is_exact()
False
is_finite()

Return False since skew polynomial rings are not finite (unless the base ring is \(0\).)

EXAMPLES:

sage: k.<t> = GF(5^3)
sage: k.is_finite()
True
sage: Frob = k.frobenius_endomorphism()
sage: S.<x> = k['x',Frob]
sage: S.is_finite()
False
is_sparse()

Return True if the elements of this polynomial ring are sparsely represented.

Warning

Since sparse skew polynomials are not yet implemented, this function always returns False.

EXAMPLES:

sage: R.<t> = RR[]
sage: sigma = R.hom([t+1])
sage: S.<x> = R['x',sigma]
sage: S.is_sparse()
False
lagrange_polynomial(points)

Return the minimal-degree polynomial which interpolates the given points.

More precisely, given \(n\) pairs \((x_1, y_1), ..., (x_n, y_n) \in R^2\), where \(R\) is self.base_ring(), compute a skew polynomial \(p(x)\) such that \(p(x_i) = y_i\) for each \(i\), under the condition that the \(x_i\) are linearly independent over the fixed field of self.twist_map().

If the \(x_i\) are linearly independent over the fixed field of self.twist_map() then such a polynomial is guaranteed to exist. Otherwise, it might exist depending on the \(y_i\), but the algorithm used in this implementation does not support that, and so an error is always raised.

INPUT:

  • points – a list of pairs (x_1, y_1),..., (x_n, y_n) of elements of the base ring of self. The \(x_i\) should be linearly independent over the fixed field of self.twist_map().

OUTPUT:

The Lagrange polynomial.

EXAMPLES:

sage: k.<t> = GF(5^3)
sage: Frob = k.frobenius_endomorphism()
sage: S.<x> = k['x',Frob]
sage: points = [(t, 3*t^2 + 4*t + 4), (t^2, 4*t)]
sage: d = S.lagrange_polynomial(points); d
x + t

sage: R.<t> = ZZ[]
sage: sigma = R.hom([t+1])
sage: T.<x> = R['x', sigma]
sage: points = [ (1, t^2 + 3*t + 4), (t, 2*t^2 + 3*t + 1), (t^2, t^2 + 3*t + 4) ]
sage: p = T.lagrange_polynomial(points); p
((-t^4 - 2*t - 3)/-2)*x^2 + (-t^4 - t^3 - t^2 - 3*t - 2)*x + (-t^4 - 2*t^3 - 4*t^2 - 10*t - 9)/-2
sage: p.multi_point_evaluation([1, t, t^2]) == [ t^2 + 3*t + 4, 2*t^2 + 3*t + 1, t^2 + 3*t + 4 ]
True

If the \(x_i\) are linearly dependent over the fixed field of self.twist_map(), then an error is raised:

sage: T.lagrange_polynomial([ (t, 1), (2*t, 3) ])
Traceback (most recent call last):
...
ValueError: the given evaluation points are linearly dependent over the fixed field of the twist map, so a Lagrange polynomial could not be determined (and might not exist).
minimal_vanishing_polynomial(eval_pts)

Return the minimal-degree, monic skew polynomial which vanishes at all the given evaluation points.

The degree of the vanishing polynomial is at most the length of eval_pts. Equality holds if and only if the elements of eval_pts are linearly independent over the fixed field of self.twist_map().

INPUT:

  • eval_pts – list of evaluation points which are linearly independent over the fixed field of the twist map of the associated skew polynomial ring

OUTPUT:

The minimal vanishing polynomial.

EXAMPLES:

sage: k.<t> = GF(5^3)
sage: Frob = k.frobenius_endomorphism()
sage: S.<x> = k['x',Frob]
sage: eval_pts = [1, t, t^2]
sage: b = S.minimal_vanishing_polynomial(eval_pts); b
x^3 + 4

The minimal vanishing polynomial evaluates to 0 at each of the evaluation points:

sage: eval = b.multi_point_evaluation(eval_pts); eval
[0, 0, 0]

If the evaluation points are linearly dependent over the fixed field of the twist map, then the returned polynomial has lower degree than the number of evaluation points:

sage: S.minimal_vanishing_polynomial([t])
x + 3*t^2 + 3*t
sage: S.minimal_vanishing_polynomial([t, 3*t])
x + 3*t^2 + 3*t
ngens()

Return the number of generators of this skew polynomial ring, which is 1.

EXAMPLES:

sage: R.<t> = RR[]
sage: sigma = R.hom([t+1])
sage: S.<x> = R['x',sigma]
sage: S.ngens()
1
parameter(n=0)

Return the indeterminate generator of this skew polynomial ring.

INPUT:

  • n – index of generator to return (default: 0). Exists for compatibility with other polynomial rings.

EXAMPLES:

sage: R.<t> = QQ[]
sage: sigma = R.hom([t+1])
sage: S.<x> = R['x',sigma]; S
Skew Polynomial Ring in x over Univariate Polynomial Ring in t over Rational Field twisted by t |--> t + 1
sage: y = S.gen(); y
x
sage: y == x
True
sage: y is x
True
sage: S.gen(0)
x

This is also known as the parameter:

sage: S.parameter() is S.gen()
True
random_element(degree=2, monic=False, *args, **kwds)

Return a random skew polynomial in self.

INPUT:

  • degree – (default: 2) integer with degree or a tuple of integers with minimum and maximum degrees
  • monic – (default: False) if True, return a monic skew polynomial
  • *args, **kwds – passed on to the random_element method for the base ring

OUTPUT:

Skew polynomial such that the coefficients of \(x^i\), for \(i\) up to degree, are random elements from the base ring, randomized subject to the arguments *args and **kwds.

EXAMPLES:

sage: k.<t> = GF(5^3)
sage: Frob = k.frobenius_endomorphism()
sage: S.<x> = k['x', Frob]
sage: S.random_element()  # random
(2*t^2 + 3)*x^2 + (4*t^2 + t + 4)*x + 2*t^2 + 2
sage: S.random_element(monic=True)  # random
x^2 + (2*t^2 + t + 1)*x + 3*t^2 + 3*t + 2

Use degree to obtain polynomials of higher degree

sage: p = S.random_element(degree=5) # random (t^2 + 3*t)*x^4 + (4*t + 4)*x^3 + (4*t^2 + 4*t)*x^2 + (2*t^2 + 1)*x + 3

When monic is False, the returned skew polynomial may have a degree less than degree (it happens when the random leading coefficient is zero). However, if monic is True, this can’t happen:

sage: p = S.random_element(degree=4, monic=True)
sage: p.leading_coefficient() == S.base_ring().one()
True
sage: p.degree() == 4
True

If a tuple of two integers is given for the degree argument, a random integer will be chosen between the first and second element of the tuple as the degree, both inclusive:

sage: S.random_element(degree=(2,7))  # random
(3*t^2 + 1)*x^4 + (4*t + 2)*x^3 + (4*t + 1)*x^2
 + (t^2 + 3*t + 3)*x + 3*t^2 + 2*t + 2

If the first tuple element is greater than the second, a a ValueError is raised:

sage: S.random_element(degree=(5,4))
Traceback (most recent call last):
...
ValueError: first degree argument must be less or equal to the second
twist_map(n=1)

Return the twist map, the automorphism of the base ring of self, iterated n times.

INPUT:

  • n - an integer (default: 1)

OUTPUT:

n-th iterative of the twist map of this skew polynomial ring.

EXAMPLES:

sage: R.<t> = QQ[]
sage: sigma = R.hom([t+1])
sage: S.<x> = R['x',sigma]
sage: S.twist_map()
Ring endomorphism of Univariate Polynomial Ring in t over Rational Field
  Defn: t |--> t + 1
sage: S.twist_map() == sigma
True
sage: S.twist_map(10)
Ring endomorphism of Univariate Polynomial Ring in t over Rational Field
  Defn: t |--> t + 10

If n in negative, Sage tries to compute the inverse of the twist map:

sage: k.<t> = GF(5^3)
sage: Frob = k.frobenius_endomorphism()
sage: T.<y> = k['y',Frob]
sage: T.twist_map(-1)
Frobenius endomorphism t |--> t^(5^2) on Finite Field in t of size 5^3

Sometimes it fails, even if the twist map is actually invertible:

sage: S.twist_map(-1)
Traceback (most recent call last):
...
NotImplementedError: inversion of the twist map Ring endomorphism of Univariate Polynomial Ring in t over Rational Field
      Defn: t |--> t + 1