RISC JKU

ore_operator

class ore_operator.OreOperator(parent, is_gen=False, construct=False)

An Ore operator. This is an abstract class whose instances represent elements of OreAlgebra.

In addition to usual RingElement features, Ore operators provide coefficient extraction functionality and the possibility of letting an operator act on another object. The latter is provided through __call__.

base_extend(R)

Return a copy of this operator but with coefficients in R, if there is a natural map from coefficient ring of self to R.

EXAMPLES:

sage: L = OreAlgebra(QQ['x'], 'Dx').random_element()
sage: L = L.base_extend(QQ['x'].fraction_field())
sage: L.parent()
Univariate Ore algebra in Dx over Fraction Field of Univariate Polynomial Ring in x over Rational Field
base_ring()

Return the base ring of the parent of self.

EXAMPLES:

sage: OreAlgebra(QQ['x'], 'Dx').random_element().base_ring()
Univariate Polynomial Ring in x over Rational Field
change_ring(R)

Return a copy of this operator but with coefficients in R, if at all possible.

EXAMPLES:

sage: R.<x> = QQ['x']
sage: A.<Dx> = OreAlgebra(R, 'Dx')
sage: op = Dx^2 + 5*x*Dx + 1
sage: op.parent()
Univariate Ore algebra in Dx over Univariate Polynomial Ring in x over Rational Field
sage: op = op.change_ring(R.fraction_field())
sage: op.parent()
Univariate Ore algebra in Dx over Fraction Field of Univariate Polynomial Ring in x over Rational Field
coefficients()

Return the coefficients of the monomials appearing in self.

constant_coefficient()

Return the coefficient of \(\partial^0\) of this operator.

content(proof=True)

Returns the content of self.

If the base ring of self‘s parent is a field, the method returns the leading coefficient.

If the base ring is not a field, then it is a polynomial ring. In this case, the method returns the greatest common divisor of the nonzero coefficients of self. If the base ring does not know how to compute gcds, the method returns \(1\).

If proof is set to False, the gcd of two random linear combinations of the coefficients is taken instead of the gcd of all the coefficients.

EXAMPLES:

sage: R.<x> = ZZ['x']
sage: A.<Dx> = OreAlgebra(R, 'Dx')
sage: (5*x^2*Dx + 10*x).content()
5*x
sage: R.<x> = QQ['x']
sage: A.<Dx> = OreAlgebra(R, 'Dx')
sage: (5*x^2*Dx + 10*x).content()
x
sage: (5*x^2*Dx + 10*x).content(proof=False)
x
sage: R.<x> = QQ['x']
sage: A.<Dx> = OreAlgebra(R.fraction_field(), 'Dx')
sage: (5*x^2*Dx + 10*x).content()
5*x^2
denominator()

Return a denominator of self.

If the base ring of the algebra of self is not a field, this returns the one element of the base ring.

If the base ring is a field, then it is the fraction field of a polynomial ring. In this case, the method returns the least common multiple of the denominators of all the coefficients of self. It is an element of the polynomial ring.

EXAMPLES:

sage: R.<x> = ZZ['x']
sage: A.<Dx> = OreAlgebra(R.fraction_field(), 'Dx')
sage: op = (5*x+3)/(3*x+5)*Dx + (7*x+1)/(2*x+5)
sage: op.denominator()
6*x^2 + 25*x + 25
sage: R.<x> = QQ['x']
sage: A.<Dx> = OreAlgebra(R.fraction_field(), 'Dx')
sage: op = (5*x+3)/(3*x+5)*Dx + (7*x+1)/(2*x+5)
sage: op.denominator()
x^2 + 25/6*x + 25/6
dict()

Return a sparse dictionary representation of this operator.

exponents()

Return the exponents of the monomials appearing in self.

is_gen()

Return True if this operator is one of the generators of the parent Ore algebra.

Important - this function doesn’t return True if self equals the generator; it returns True if self is the generator.

is_monic()

Returns True if this polynomial is monic. The zero operator is by definition not monic.

EXAMPLES:

sage: R.<x> = QQ['x']
sage: A.<Dx> = OreAlgebra(R, 'Dx')
sage: (Dx^3 + (5*x+3)*Dx + (71*x+1)).is_monic()
True
sage: ((5*x+3)*Dx^2 + (71*x+1)).is_monic()
False 
is_monomial()

Returns True if self is a monomial, i.e., a power product of the generators.

is_primitive(n=None, n_prime_divs=None)

Returns True if this operator’s content is a unit of the base ring.

is_unit()

Return True if this operator is a unit.

EXAMPLES:

sage: R.<x> = QQ['x']
sage: A.<Dx> = OreAlgebra(R, 'Dx')
sage: A(x).is_unit()
False
sage: A.<Dx> = OreAlgebra(R.fraction_field(), 'Dx')
sage: A(x).is_unit()
True
leading_coefficient()

Return the leading coefficient of this operator.

list()

Return a new copy of the list of the underlying elements of self.

map_coefficients(f, new_base_ring=None)

Returns the operator obtained by applying f to the non-zero coefficients of self.

monic()

Return this operator divided from the left by its leading coefficient. Does not change this operator. If the leading coefficient does not have a multiplicative inverse in the base ring of self‘s parent, the the method returns an element of a suitably extended algebra.

EXAMPLES:

sage: R.<x> = QQ['x']
sage: A.<Dx> = OreAlgebra(R, 'Dx')
sage: (x*Dx + 1).monic()
Dx + 1/x
sage: _.parent()
Univariate Ore algebra in Dx over Fraction Field of Univariate Polynomial Ring in x over Rational Field
normalize(proof=True)

Returns a normal form of self.

Call two operators \(A,B\) equivalent iff there exist nonzero elements \(p,q\) of the base ring such that \(p*A=q*B\). Then \(A\) and \(B\) are equivalent iff their normal forms as computed by this method agree.

The normal form is a left multiple of self by an element of (the fraction field of) the base ring. An attempt is made in choosing a “simple” representative of the equivalence class.

EXAMPLES:

sage: R.<x> = QQ['x']
sage: A.<Dx> = OreAlgebra(R, 'Dx')
sage: (10*(x+1)*Dx - 5*x).normalize()
(x + 1)*Dx - 1/2*x
numerator()

Return a numerator of self.

If the base ring of self‘s parent is not a field, this returns self.

If the base ring is a field, then it is the fraction field of a polynomial ring. In this case, the method returns self.denominator()*self and tries to cast the result into the Ore algebra whose base ring is just the polynomial ring. If this fails (for example, because some \(\sigma\) maps a polynomial to a rational function), the result will be returned as element of the original algebra.

EXAMPLES:

sage: R.<x> = ZZ['x']
sage: A.<Dx> = OreAlgebra(R.fraction_field(), 'Dx')
sage: op = (5*x+3)/(3*x+5)*Dx + (7*x+1)/(2*x+5)
sage: op.numerator()
(10*x^2 + 31*x + 15)*Dx + 21*x^2 + 38*x + 5
sage: R.<x> = QQ['x']
sage: A.<Dx> = OreAlgebra(R.fraction_field(), 'Dx')
sage: op = (5*x+3)/(3*x+5)*Dx + (7*x+1)/(2*x+5)
sage: op.numerator()
(5/3*x^2 + 31/6*x + 5/2)*Dx + 7/2*x^2 + 19/3*x + 5/6          
prec()

Return the precision of this operator. This is always infinity, since operators are of infinite precision by definition (there is no big-oh).

primitive_part(proof=True)

Returns the primitive part of self.

It is obtained by dividing self from the left by self.content().

The proof option is passed on to the content computation.

EXAMPLES:

sage: R.<x> = ZZ['x']
sage: A.<Dx> = OreAlgebra(R, 'Dx')
sage: (5*x^2*Dx + 10*x).primitive_part()
x*Dx + 2
sage: A.<Dx> = OreAlgebra(R.fraction_field(), 'Dx')
sage: (5*x^2*Dx + 10*x).primitive_part()
Dx + 2/x
quo_rem(other)

Quotient and remainder.

EXAMPLES:

sage: R.<x> = QQ['x']
sage: A.<Dx> = OreAlgebra(R.fraction_field(), 'Dx')
sage: U = (15*x^2 + 29*x + 5)*Dx^2 + (5*x^2 - 50*x - 41)*Dx - 2*x + 64
sage: V = (3*x+5)*Dx + (x-9)
sage: Q, R = U.quo_rem(V)
sage: Q*V + R == U
True 
class ore_operator.UnivariateOreOperator(parent, *data, **kwargs)

Element of an Ore algebra with a single generator and a commutative field as base ring.

annihilator_of_associate(other, solver=None)

Computes an operator \(L\) with \(L(other(f))=0\) for all \(f\) with \(self(f)=0\).

EXAMPLES:

sage: R.<x> = ZZ['x']
sage: A.<Dx> = OreAlgebra(R, 'Dx')
sage: (Dx^2 + x*Dx + 5).annihilator_of_associate(Dx + 7*x+3)
(-42*x^2 - 39*x - 7)*Dx^2 + (-42*x^3 - 39*x^2 + 77*x + 39)*Dx - 168*x^2 - 174*x - 61
sage: A.<Sx> = OreAlgebra(R, 'Sx')
(-42*x^2 - 88*x - 35)*Sx^2 + (-42*x^3 - 130*x^2 - 53*x + 65)*Sx - 210*x^2 - 860*x - 825
annihilator_of_polynomial(poly, solver=None)

Returns an annihilating operator of a polynomial expression evaluated at solutions of self.

INPUT:

  • poly – a multivariate polynomial, say in the variables \(x0,x1,x2,...\), with coefficients in the base ring of the parent of self. The number of variables in the parent of poly must be at least the order of self.

OUTPUT:

An operator \(L\) with the following property. Let \(A\) be the parent of self. For a function \(f\) write \(Df,D^2f,...\) for the functions obtained from \(f\) by letting the generator of \(A\) act on them. Then the output \(L\) is such that for all \(f\) annihilated by self we have \(L(p(f,Df,D^2f,...))=0\), where \(p\) is the input polynomial.

The method requires that a product rule is associated to \(A\). (See docstring of OreAlgebra for information about product rules.)

NOTE:

Even for small input, the output operator may be very large, and its computation may need a lot of time.

EXAMPLES:

sage: K.<x> = ZZ['x']; K = K.fraction_field(); R.<n> = K['n']
sage: A.<Sn> = OreAlgebra(R, 'Sn'); R.<y0,y1,y2> = K['n'][]
sage: L = (n+2)*Sn^2 - (2*n+3)*x*Sn + (n+1)
sage: L.annihilator_of_polynomial(y1^2-y2*y0) # Turan's determinant
(-2*n^4 - 31*n^3 - 177*n^2 - 440*n - 400)*Sn^3 + ((8*x^2 - 2)*n^4 + (100*x^2 - 21)*n^3 + (462*x^2 - 75)*n^2 + (935*x^2 - 99)*n + 700*x^2 - 28)*Sn^2 + ((-8*x^2 + 2)*n^4 + (-92*x^2 + 27)*n^3 + (-390*x^2 + 129)*n^2 + (-721*x^2 + 261)*n - 490*x^2 + 190)*Sn + 2*n^4 + 17*n^3 + 51*n^2 + 64*n + 28
sage: M = L.annihilator_of_associate(Sn).symmetric_power(2).lclm(L.annihilator_of_associate(Sn^2).symmetric_product(L)) # same by lower level methods. 
sage: M.order() # overshoots
7
sage: M % L.annihilator_of_polynomial(y1^2-y2*y0) 
0

sage: K = ZZ; R.<x> = K['x']
sage: A.<Dx> = OreAlgebra(R, 'Dx'); R.<y0,y1> = K['x'][]
sage: L = (2*x+3)*Dx^2 + (4*x+5)*Dx + (6*x+7)
sage: L.annihilator_of_polynomial((2*x+1)*y0^2-y1^2)
(-16*x^7 - 112*x^6 - 312*x^5 - 388*x^4 - 85*x^3 + 300*x^2 + 303*x + 90)*Dx^3 + (-96*x^7 - 600*x^6 - 1420*x^5 - 1218*x^4 + 747*x^3 + 2285*x^2 + 1623*x + 387)*Dx^2 + (-320*x^7 - 1920*x^6 - 4288*x^5 - 3288*x^4 + 2556*x^3 + 6470*x^2 + 4322*x + 1014)*Dx - 384*x^7 - 2144*x^6 - 4080*x^5 - 1064*x^4 + 7076*x^3 + 10872*x^2 + 6592*x + 1552
coefficients()

Return the coefficients of the monomials appearing in self.

coeffs(padd=-1)

Return the coefficient vector of this operator.

If the degree is less than the number given in the optional argument, the list is padded with zeros so as to ensure that the output has length padd + 1.

EXAMPLES:

sage: A.<Sx> = OreAlgebra(ZZ['x'], 'Sx')
sage: (5*Sx^3-4).coeffs()
[-4, 0, 0, 5]
sage: (5*Sx^3-4).coeffs(padd=5)
[-4, 0, 0, 5, 0, 0]
sage: (5*Sx^3-4).coeffs(padd=1)
[-4, 0, 0, 5]
companion_matrix()

Returns the companion matrix of self.

If \(r\) is the order of self and \(y\) is annihilated by self, then the companion matrix as computed by this method is an \(r\times r\) matrix \(M\) such that \([\partial y,\partial^2 y,\dots,\partial^r y] = M [y,\partial y,\dots,\partial^{r-1}y]^T\).

In the shift case, if \(c_i\) is a sequence annihilated by self, then also \([c_{i+1}, c_{i+2}, \ldots, c_{i+r}]^T = M(i) [c_i, c_{i+1}, \ldots, c_{i+r-1}]^T\)

EXAMPLES:

sage: R.<n> = QQ['n']
sage: A.<Sn> = OreAlgebra(R, 'Sn')
sage: M = ((-n-4)*Sn**2 + (5+2*n)*Sn + (3+3*n)).companion_matrix()
sage: M
[                0                 1]
[(3*n + 3)/(n + 4) (2*n + 5)/(n + 4)]
sage: initial = Matrix([[1],[1]])
sage: [prod(M(k) for k in range(n, -1, -1)) * initial for n in range(10)]
[
[1]  [2]  [4]  [ 9]  [21]  [ 51]  [127]  [323]  [ 835]  [2188]
[2], [4], [9], [21], [51], [127], [323], [835], [2188], [5798]
]
constant_coefficient()

Return the coefficient of \(\partial^0\) of this operator.

exponents()

Return the exponents of the monomials appearing in self.

gcrd(*other, **kwargs)

Returns the GCRD of self and other. It is possible to specify which remainder sequence should be used.

is_gen()

Return True if this operator is one of the generators of the parent Ore algebra.

Important - this function doesn’t return True if self equals the generator; it returns True if self is the generator.

is_monic()

Returns True if this polynomial is monic. The zero operator is by definition not monic.

EXAMPLES:

sage: R.<x> = QQ['x']
sage: A.<Dx> = OreAlgebra(R, 'Dx')
sage: (Dx^3 + (5*x+3)*Dx + (71*x+1)).is_monic()
True
sage: ((5*x+3)*Dx^2 + (71*x+1)).is_monic()
False 
is_unit()

Return True if this operator is a unit.

EXAMPLES:

sage: R.<x> = QQ['x']
sage: A.<Dx> = OreAlgebra(R, 'Dx')
sage: A(x).is_unit()
False
sage: A.<Dx> = OreAlgebra(R.fraction_field(), 'Dx')
sage: A(x).is_unit()
True
lclm(*other, **kwargs)

Computes the least common left multiple of self and other.

That is, it returns an operator \(L\) of minimal order such that there exist \(U\) and \(V\) with \(L=U*self=V*other\). The base ring of the parent of \(U\) and \(V\) is the fraction field of the base ring of the parent of self and other. The parent of \(L\) is the same as the parent of the input operators.

If more than one operator is given, the function computes the lclm of all the operators.

The optional argument algorithm allows to select between the following methods.

  • linalg (default) – makes an ansatz for cofactors and solves a linear system over the base ring. Through the optional argument solver, a callable object can be provided which the function should use for computing the kernel of matrices with entries in the Ore algebra’s base ring.
  • euclid – uses the extended Euclidean algorithm to compute a minimal syzygy between the operators in the input. Further optional arguments can be passed as explained in the docstring of xgcrd.
  • guess – computes the first terms of a solution of self and other and guesses from these a minimal operator annihilating a generic linear combination. Unless the input are recurrence operators, an keyword argument to_list has to be present which specifies a function for computing the terms (input: an operator, a list of initial values, and the desired number of terms). This method is heuristic. It may be much faster than the others, but with low probability its output is incorrect or it aborts with an error.

EXAMPLES:

sage: R.<x> = ZZ['x']
sage: Alg.<Dx> = OreAlgebra(R, 'Dx')
sage: A = 5*(x+1)*Dx + (x - 7); B = (3*x+5)*Dx - (8*x+1)
sage: L = A.lclm(B)
(-645*x^4 - 2155*x^3 - 1785*x^2 + 475*x + 750)*Dx^2 + (1591*x^4 + 3696*x^3 + 3664*x^2 + 2380*x + 725)*Dx + 344*x^4 - 2133*x^3 - 2911*x^2 - 1383*x - 1285
sage: A*B
(15*x^2 + 40*x + 25)*Dx^2 + (-37*x^2 - 46*x - 25)*Dx - 8*x^2 + 15*x - 33
sage: B.lclm(A*B)
(-15*x^2 - 40*x - 25)*Dx^2 + (37*x^2 + 46*x + 25)*Dx + 8*x^2 - 15*x + 33
sage: B.lclm(L, A*B)
(15*x^2 + 40*x + 25)*Dx^2 + (-37*x^2 - 46*x - 25)*Dx - 8*x^2 + 15*x - 33
leading_coefficient()

Return the leading coefficient of this operator.

map_coefficients(f, new_base_ring=None)

Returns the polynomial obtained by applying f to the non-zero coefficients of self.

order()

Returns the order of this operator, which is defined as the maximal power \(i\) of the generator which has a nonzero coefficient. The zero operator has order \(-1\).

quo_rem(other, fractionFree=False)

Quotient and remainder.

EXAMPLES:

sage: R.<x> = QQ['x']
sage: A.<Dx> = OreAlgebra(R.fraction_field(), 'Dx')
sage: U = (15*x^2 + 29*x + 5)*Dx^2 + (5*x^2 - 50*x - 41)*Dx - 2*x + 64
sage: V = (3*x+5)*Dx + (x-9)
sage: Q, R = U.quo_rem(V)
sage: Q*V + R == U
True 
resultant(other)

Returns the resultant between this operator and other.

INPUT:

  • other – some operator that lives in the same algebra as self.

OUTPUT:

The resultant between self and other, which is defined as the determinant of the \((n+m) x (n+m)\) matrix \([ A, D*A, ..., D^{m-1}*A, B, D*B, ..., D^{n-1}*B ]\) constructed from the coefficient vectors of the operators obtained from self and other by multiplying them by powers of the parent’s generator.

EXAMPLES:

sage: R.<x> = ZZ['x']
sage: A.<Dx> = OreAlgebra(R, 'Dx')
sage: L1 = (5*x+3)*Dx^3 + (7*x+4)*Dx^2 + (3*x+2)*Dx + (4*x-1)
sage: L2 = (8*x-7)*Dx^2 + (x+9)*Dx + (2*x-3)
sage: L1.resultant(L2)
2934*x^5 - 8200*x^4 + 32161*x^3 - 83588*x^2 - 67505*x + 42514
sage: L2.resultant(L1)
-2934*x^5 + 8200*x^4 - 32161*x^3 + 83588*x^2 + 67505*x - 42514

sage: R.<x> = ZZ['x']
sage: A.<Sx> = OreAlgebra(R, 'Sx')
sage: L1 = (5*x+3)*Sx^3 + (7*x+4)*Sx^2 + (3*x+2)*Sx + (4*x-1)
sage: L2 = (8*x-7)*Sx^2 + (x+9)*Sx + (2*x-3)
sage: L1.resultant(L2)
2934*x^5 + 3536*x^4 + 11142*x^3 + 16298*x^2 - 1257*x - 2260
sage: L2.resultant(L1)
-2934*x^5 - 3536*x^4 - 11142*x^3 - 16298*x^2 + 1257*x + 2260
symmetric_power(exp, solver=None)

Returns a symmetric power of this operator.

The \(n\) th symmetric power of an operator \(L\) is a minimal order operator \(Q\) such that for all “functions” \(f\) annihilated by \(L\) the operator \(Q\) annihilates the function \(f^n\).

For further information, see the docstring of symmetric_product.

EXAMPLES:

sage: R.<x> = ZZ['x']
sage: A.<Dx> = OreAlgebra(R, 'Dx')
sage: (Dx^2 + x*Dx - 2).symmetric_power(3)
Dx^4 + 6*x*Dx^3 + (11*x^2 - 16)*Dx^2 + (6*x^3 - 53*x)*Dx - 36*x^2 + 24
sage: A.<Sx> = OreAlgebra(R, 'Sx')
sage: (Sx^2 + x*Sx - 2).symmetric_power(2)
-x*Sx^3 + (x^3 + 2*x^2 + 3*x + 2)*Sx^2 + (2*x^3 + 2*x^2 + 4*x)*Sx - 8*x - 8
sage: A.random_element().symmetric_power(0)
Sx - 1
symmetric_product(other, solver=None)

Returns the symmetric product of self and other.

The symmetric product of two operators \(A\) and \(B\) is a minimal order operator \(C\) such that for all “functions” \(f\) and \(g\) with \(A.f=B.g=0\) we have \(C.(fg)=0\).

The method requires that a product rule is associated to the Ore algebra where self and other live. (See docstring of OreAlgebra for information about product rules.)

If no solver is specified, the the Ore algebra’s solver is used.

EXAMPLES:

sage: R.<x> = ZZ['x']
sage: A.<Dx> = OreAlgebra(R, 'Dx')
sage: (Dx - 1).symmetric_product(x*Dx - 1)
x*Dx - x - 1
sage: (x*Dx - 1).symmetric_product(Dx - 1)
x*Dx - x - 1
sage: ((x+1)*Dx^2 + (x-1)*Dx + 8).symmetric_product((x-1)*Dx^2 + (2*x+3)*Dx + (8*x+5))
(-29*x^8 + 4*x^7 + 55*x^6 + 34*x^5 + 23*x^4 - 80*x^3 - 95*x^2 + 42*x + 46)*Dx^4 + (-174*x^8 - 150*x^7 - 48*x^6 + 294*x^5 + 864*x^4 + 646*x^3 - 232*x^2 - 790*x - 410)*Dx^3 + (-783*x^8 - 1661*x^7 + 181*x^6 + 1783*x^5 + 3161*x^4 + 3713*x^3 - 213*x^2 - 107*x + 1126)*Dx^2 + (-1566*x^8 - 5091*x^7 - 2394*x^6 - 2911*x^5 + 10586*x^4 + 23587*x^3 + 18334*x^2 + 2047*x - 5152)*Dx - 2552*x^8 - 3795*x^7 - 8341*x^6 - 295*x^5 + 6394*x^4 + 24831*x^3 + 35327*x^2 + 23667*x + 13708

sage: A.<Sx> = OreAlgebra(R, 'Sx')
sage: (Sx - 2).symmetric_product(x*Sx - (x+1))
x*Sx - 2*x - 2
sage: (x*Sx - (x+1)).symmetric_product(Sx - 2)
x*Sx - 2*x - 2
sage: ((x+1)*Sx^2 + (x-1)*Sx + 8).symmetric_product((x-1)*Sx^2 + (2*x+3)*Sx + (8*x+5))
(8*x^8 + 13*x^7 - 300*x^6 - 1640*x^5 - 3698*x^4 - 4373*x^3 - 2730*x^2 - 720*x)*Sx^4 + (-16*x^8 - 34*x^7 + 483*x^6 + 1947*x^5 + 2299*x^4 + 2055*x^3 + 4994*x^2 + 4592*x)*Sx^3 + (64*x^8 - 816*x^7 - 1855*x^6 + 21135*x^5 + 76919*x^4 + 35377*x^3 - 179208*x^2 - 283136*x - 125440)*Sx^2 + (-1024*x^7 - 1792*x^6 + 39792*x^5 + 250472*x^4 + 578320*x^3 + 446424*x^2 - 206528*x - 326144)*Sx + 32768*x^6 + 61440*x^5 - 956928*x^4 - 4897984*x^3 - 9390784*x^2 - 7923200*x - 2329600
valuation()

Returns the valuation of this operator, which is defined as the minimal power \(i\) of the generator which has a nonzero coefficient. The zero operator has order \(\infty\).

xgcrd(other, prs=None)

When called for two operators p,q, this will return their GCRD g together with two operators s and t such that sp+tq=g. It is possible to specify which remainder sequence should be used.

xlclm(other)

Computes the least common left multiple of self and other along with the appropriate cofactors.

That is, it returns a triple \((L,U,V)\) such that \(L=U*self=V*other\) and \(L\) has minimal possible order. The base ring of the parent of \(U\) and \(V\) is the fraction field of the base ring of the parent of self and other. The parent of \(L\) is the same as the parent of the input operators.

EXAMPLES:

sage: R.<x> = QQ['x']
sage: Alg.<Dx> = OreAlgebra(R, 'Dx')
sage: A = 5*(x+1)*Dx + (x - 7); B = (3*x+5)*Dx - (8*x+1)
sage: L, U, V = A.xlclm(B)
sage: L == U*A
True
sage: L == V*B
True
sage: L.parent()
Univariate Ore algebra in Dx over Univariate Polynomial Ring in x over Rational Field
sage: U.parent()
Univariate Ore algebra in Dx over Fraction Field of Univariate Polynomial Ring in x over Rational Field

Previous topic

ore_algebra

Next topic

ore_operator_1_1