RISC JKU

ore_operator_1_1

Special classes for operators living in algebras with one generator with base rings that also have one generator.

class ore_operator_1_1.UnivariateDifferenceOperatorOverUnivariateRing(parent, *data, **kwargs)

Element of an Ore algebra K(x)[F], where F is the forward difference operator F f(x) = f(x+1) - f(x)

indicial_polynomial(*args, **kwargs)

Computes the indicial polynomial of self at (a root of) p.

The indicial polynomial is a polynomial in the given variable var with coefficients in the fraction field of the base ring’s base ring.

The precise meaning of this polynomial may depend on the parent of self. A minimum requirement is that if self has a rational solution whose denominator contains sigma.factorial(p, e) but neither sigma(p, -1)*sigma.factorial(p, e) nor sigma.factorial(p, e + 1), then -e is a root of this polynomial.

Applied to \(p=1/x\), the maximum integer root of the output should serve as a degree bound for the polynomial solutions of self.

This method is a stub. Depending on the particular subclass, restrictions on p may apply.

spread(p=0)

Returns the spread of this operator.

This is the set of integers \(i\) such that sigma(self[0], i) and sigma(self[r], -r) have a nontrivial common factor, where sigma is the shift of the parent’s algebra and \(r\) is the order of self.

If the optional argument \(p\) is given, the method is applied to gcd(self[0], p) instead of self[0].

The output set contains \(\infty\) if the constant coefficient of self is zero.

EXAMPLES:

sage: R.<x> = ZZ['x']; A.<Sx> = OreAlgebra(R, 'Sx');
sage: ((x+5)*Sx - x).spread()
[4]
sage: ((x+5)*Sx - x).lclm((x+19)*Sx - x).spread()
[3, 4, 17, 18]
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
to_D(alg)

Returns a differential operator which annihilates every power series (about the origin) whose coefficient sequence is annihilated by self. The output operator may not be minimal.

INPUT:

  • alg – the Ore algebra in which the output should be expressed. The algebra must satisfy alg.base_ring().base_ring() == self.base_ring().base_ring() and alg.is_D() is not False. Instead of an algebra object, also a string can be passed as argument. This amounts to specifying an Ore algebra over self.base_ring() with the standard derivation with respect to self.base_ring().gen().

EXAMPLES:

sage: Rn.<n> = ZZ['n']; Rx.<x> = ZZ['x']
sage: A.<Fn> = OreAlgebra(Rn, 'Fn')
sage: B.<Dx> = OreAlgebra(Rx, 'Dx')
sage: Fn.to_D(B)
(-x + 1)*Dx - 1
sage: ((n+1)*Fn - 1).to_D(B)
(-x^2 + x)*Dx^2 + (-4*x + 1)*Dx - 2
sage: (x*Dx-1).to_F(A).to_D(B)
x*Dx - 1
to_S(alg)

Returns the differential operator corresponding to self

INPUT:

  • alg – the Ore algebra in which the output should be expressed. The algebra must satisfy alg.base_ring().base_ring() == self.base_ring().base_ring() and alg.is_S() is not False. Instead of an algebra object, also a string can be passed as argument. This amounts to specifying an Ore algebra over self.base_ring() with a standard shift with respect to self.base_ring().gen().

EXAMPLES:

sage: R.<x> = ZZ['x']
sage: A.<Fx> = OreAlgebra(R, 'Fx')
sage: (Fx^4).to_S(OreAlgebra(R, 'Sx'))
Sx^4 - 4*Sx^3 + 6*Sx^2 - 4*Sx + 1
sage: (Fx^4).to_S('Sx')
Sx^4 - 4*Sx^3 + 6*Sx^2 - 4*Sx + 1
to_T(alg)

Returns a differential operator, expressed in terms of the Euler derivation, which annihilates every power series (about the origin) whose coefficient sequence is annihilated by self. The output operator may not be minimal.

INPUT:

  • alg – the Ore algebra in which the output should be expressed. The algebra must satisfy alg.base_ring().base_ring() == self.base_ring().base_ring() and alg.is_T() is not False. Instead of an algebra object, also a string can be passed as argument. This amounts to specifying an Ore algebra over self.base_ring() with the Euler derivation with respect to self.base_ring().gen().

EXAMPLES:

sage: Rn.<n> = ZZ['n']; Rx.<x> = ZZ['x']
sage: A.<Fn> = OreAlgebra(Rn, 'Fn')
sage: B.<Tx> = OreAlgebra(Rx, 'Tx')
sage: Fn.to_T(B)
(-x + 1)*Tx - x
sage: ((n+1)*Fn - 1).to_T(B)
(-x + 1)*Tx^2 - 3*x*Tx - 2*x
sage: (x*Tx-1).to_F(A).to_T(B)
x*Tx^2 + (x - 1)*Tx
to_list(*args, **kwargs)

Computes the terms of some sequence annihilated by self.

INPUT:

  • init – a vector (or list or tuple) of initial values. The components must be elements of self.base_ring().base_ring().fraction_field(). If the length is more than self.order(), we do not check whether the given terms are consistent with self.
  • n – desired number of terms.
  • start (optional) – index of the sequence term which is represented by the first entry of init. Defaults to zero.
  • append (optional) – if True, the computed terms are appended to init list. Otherwise (default), a new list is created.
  • padd (optional) – if True, the vector of initial values is implicitely prolonged to the left (!) by zeros if it is too short. Otherwise (default), the method raises a ValueError if init is too short.

OUTPUT:

A list of n terms whose \(k\) th component carries the sequence term with index start+k. Terms whose calculation causes an error are represented by None.

EXAMPLES:

sage: R = ZZ['x']['n']; x = R('x'); n = R('n')
sage: A.<Sn> = OreAlgebra(R, 'Sn')
sage: L = ((n+2)*Sn^2 - x*(2*n+3)*Sn + (n+1))
sage: L.to_list([1, x], 5)
[1, x, (3*x^2 - 1)/2, (5*x^3 - 3*x)/2, (35*x^4 - 30*x^2 + 3)/8]
sage: polys = L.to_list([1], 5, padd=True)
[1, x, (3*x^2 - 1)/2, (5*x^3 - 3*x)/2, (35*x^4 - 30*x^2 + 3)/8]
sage: L.to_list([polys[3], polys[4]], 8, start=3)
[(5*x^3 - 3*x)/2, (35*x^4 - 30*x^2 + 3)/8, (63*x^5 - 70*x^3 + 15*x)/8, (231*x^6 - 315*x^4 + 105*x^2 - 5)/16, (429*x^7 - 693*x^5 + 315*x^3 - 35*x)/16]

sage: ((n-5)*Sn - 1).to_list([1], 10)
[1, 1/-5, 1/20, 1/-60, 1/120, -1/120, None, None, None, None]
class ore_operator_1_1.UnivariateDifferentialOperatorOverUnivariateRing(parent, *data, **kwargs)

Element of an Ore algebra K(x)[D], where D acts as derivation d/dx on K(x).

annihilator_of_composition(a, solver=None)

Returns an operator \(L\) which annihilates all the functions \(f(a(x))\) where \(f\) runs through the functions annihilated by self. The output operator is not necessarily of smallest possible order.

INPUT:

  • a – either an element of the base ring of the parent of self, or an element of an algebraic extension of this ring.
  • solver (optional) – a callable object which applied to a matrix with polynomial entries returns its kernel.

EXAMPLES:

sage: R.<x> = ZZ['x']
sage: K.<y> = R.fraction_field()['y']
sage: K.<y> = R.fraction_field().extension(y^3 - x^2*(x+1))
sage: A.<Dx> = OreAlgebra(R, 'Dx')
sage: (x*Dx-1).annihilator_of_composition(y) # ann for x^(2/3)*(x+1)^(1/3)
(3*x^2 + 3*x)*Dx - 3*x - 2
sage: (x*Dx-1).annihilator_of_composition(y + 2*x) # ann for 2*x + x^(2/3)*(x+1)^(1/3)
(-3*x^3 - 3*x^2)*Dx^2 + 2*x*Dx - 2
sage: (Dx - 1).annihilator_of_composition(y) # ann for exp(x^(2/3)*(x+1)^(1/3))
(243*x^6 + 810*x^5 + 999*x^4 + 540*x^3 + 108*x^2)*Dx^3 + (162*x^3 + 270*x^2 + 108*x)*Dx^2 + (-162*x^2 - 180*x - 12)*Dx - 243*x^6 - 810*x^5 - 1080*x^4 - 720*x^3 - 240*x^2 - 32*x
annihilator_of_integral()

Returns an operator \(L\) which annihilates all the indefinite integrals \(\int f\) where \(f\) runs through the functions annihilated by self. The output operator is not necessarily of smallest possible order.

EXAMPLES:

sage: R.<x> = ZZ['x']
sage: A.<Dx> = OreAlgebra(R, 'Dx')
sage: ((x-1)*Dx - 2*x).annihilator_of_integral()
(x-1)*Dx^2 - 2*x*Dx
sage: _.annihilator_of_associate(Dx)
(x-1)*Dx - 2*x
generalized_series_solutions(n=5, base_extend=True, ramification=True, exp=True)

Returns the generalized series solutions of this operator.

These are solutions of the form

\(\exp(\int_0^x \frac{p(t^{-1/s})}t dt)*q(x^{1/s},\log(x))\)

where

  • \(s\) is a positive integer (the object’s “ramification”)
  • \(p\) is in \(K[x]\) (the object’s “exponential part”)
  • \(q\) is in \(K[[x]][y]\) with \(x\nmid q\) unless \(q\) is zero (the object’s “tail”)
  • \(K\) is some algebraic extension of the base ring’s base ring.

An operator of order \(r\) has exactly \(r\) linearly independent solutions of this form. This method computes them all, unless the flags specified in the arguments rule out some of them.

At present, the method only works for operators where the base ring’s base ring is either QQ or a number field (i.e., no finite fields, no formal parameters).

INPUT:

  • n (default: 5) – minimum number of terms in the series expansions to be computed in addition to those needed to separate all solutions from each other.
  • base_extend (default: True) – whether or not the coefficients of the solutions may belong to an algebraic extension of the base ring’s base ring.
  • ramification (default: True) – whether or not the exponential parts of the solutions may involve fractional exponents.
  • exp (default: True) – set this to False if you only want solutions that have no exponential part (viz \(\deg(p)\leq0\)). If set to a positive rational number \(\alpha\), the method returns all those solutions whose exponential part involves only terms \(x^{-i/r}\) with \(i/r<\alpha\).

OUTPUT:

  • a list of ContinuousGeneralizedSeries objects forming a fundamental system for this operator.

Note

  • Different solutions may require different algebraic extensions. Thus in the list returned by this method, the coefficient fields of different series typically do not coincide.
  • If a solution involves an algebraic extension of the coefficient field, then all its conjugates are solutions, too. But only one representative is listed in the output.

ALGORITHM:

  • Ince, Ordinary Differential Equations, Chapters 16 and 17
  • Kauers/Paule, The Concrete Tetrahedron, Section 7.3

EXAMPLES:

sage: R.<x> = QQ['x']; A.<Dx> = OreAlgebra(R, 'Dx')
sage: L = (6+6*x-3*x^2) - (10*x-3*x^2-3*x^3)*Dx + (4*x^2-6*x^3+2*x^4)*Dx^2
sage: L.generalized_series_solutions()
[x^3*(1 + 3/2*x + 7/4*x^2 + 15/8*x^3 + 31/16*x^4 + O(x^5)), x^(1/2)*(1 + 3/2*x + 7/4*x^2 + 15/8*x^3 + 31/16*x^4 + O(x^5))]
sage: map(L, _)
[0, 0]

sage: L = (1-24*x+96*x^2) + (15*x-117*x^2+306*x^3)*Dx + (9*x^2-54*x^3)*Dx^2
sage: L.generalized_series_solutions(2)
sage: [x^(-1/3)*(1 + x + 8/3*x^2 + O(x^3)), x^(-1/3)*((1 + x + 8/3*x^2 + O(x^3))*log(x) + x - 59/12*x^2 + O(x^3))]
sage: map(L, _)
[0, 0]

sage: L = 216*(1+x+x^3) + x^3*(36-48*x^2+41*x^4)*Dx - x^7*(6+6*x-x^2+4*x^3)*Dx^2
sage: L.generalized_series_solutions(3)
[exp(3*x^(-2))*x^(-2)*(1 + 91/12*x^2 + O(x^3)), exp(-2*x^(-3) + x^(-1))*x^2*(1 + 41/3*x + 2849/36*x^2 + O(x^3))]
sage: map(L, _)
[0, 0]

sage: L = 9 - 49*x - 2*x^2 + 6*x^2*(7 + 5*x)*Dx + 36*(-1 + x)*x^3*Dx^2
sage: L.generalized_series_solutions()
[exp(x^(-1/2))*x^(4/3)*(1 + x^(2/2) + x^(4/2)), exp(-x^(-1/2))*x^(4/3)*(1 + x^(2/2) + x^(4/2))]
sage: L.generalized_series_solutions(ramification=False)
[]

sage: L = 2*x^3*Dx^2 + 3*x^2*Dx-1
sage: L.generalized_series_solutions()
[exp(a_0*x^(-1/2))]
sage: _[0].base_ring()
Number Field in a_0 with defining polynomial c^2 - 2
indicial_polynomial(p, var='alpha')

Computes the indicial polynomial of this operator at (a root of) \(p\).

If \(x\) is the generator of the base ring, the input may be either irreducible polynomial in \(x\) or the rational function \(1/x\).

The output is a univariate polynomial in the given variable var with coefficients in the base ring’s base ring. It has the following property: for every nonzero series solution of self in rising powers of \(p\), i.e. \(p_0 p^lpha + p_1 p^{lpha+1} + ...\), the minimal exponent \(lpha\) is a root of the indicial polynomial. The converse may not hold.

INPUT:

  • p – an irreducible polynomial in the base ring of the operator algebra, or \(1/x\).
  • var (optional) – the variable name to use for the indicial polynomial.

EXAMPLES:

sage: R.<x> = ZZ['x']; A.<Dx> = OreAlgebra(R, 'Dx');
sage: L = (x*Dx-5).lclm((x^2+1)*Dx - 7*x).lclm(Dx - 1)
sage: L.indicial_polynomial(x).factor()
(-1) * 5 * 2^2 * (alpha - 5) * (alpha - 1) * alpha
sage: L.indicial_polynomial(1/x).factor()
2 * (alpha - 7) * (alpha - 5)
sage: L.indicial_polynomial(x^2+1).factor()
5 * 7 * 2^2 * (alpha - 1) * alpha * (2*alpha - 7)
power_series_solutions(n=5)

Computes the first few terms of the power series solutions of this operator.

The method raises an error if Sage does not know how to factor univariate polynomials over the base ring’s base ring.

The base ring has to have characteristic zero.

INPUT:

  • n – minimum number of terms to be computed

OUTPUT:

A list of power series of the form \(x^lpha + ...\) with pairwise distinct exponents \(lpha\) and coefficients in the base ring’s base ring’s fraction field. All expansions are computed up to order \(k\) where \(k\) is obtained by adding the maximal \(lpha\) to the maximum of \(n\) and the order of self.

EXAMPLES:

sage: R.<x> = ZZ['x']
sage: A.<Dx> = OreAlgebra(R, 'Dx')
sage: ((1-x)*Dx - 1).power_series_solutions(10) # geometric series
[1 + x + x^2 + x^3 + x^4 + x^5 + x^6 + x^7 + x^8 + O(x^9)]
sage: (Dx - 1).power_series_solutions(5) # exp(x)
[1 + x + 1/2*x^2 + 1/6*x^3 + O(x^4)]
sage: (Dx^2 - Dx + x).power_series_solutions(5) # a 2nd order equation
[x + 1/2*x^2 + 1/6*x^3 - 1/24*x^4 + O(x^5), 1 - 1/6*x^3 - 1/24*x^4 + O(x^5)]
sage: (2*x*Dx - 1).power_series_solutions(5) # sqrt(x) is not a power series
[]
spread(p=0)

Returns the spread of this operator.

This is the set of integers \(i\) such that sigma(self[0], i) and sigma(self[r], -r) have a nontrivial common factor, where sigma is the shift of the parent’s algebra and \(r\) is the order of self.

If the optional argument \(p\) is given, the method is applied to gcd(self[0], p) instead of self[0].

The output set contains \(\infty\) if the constant coefficient of self is zero.

This method is a stub and may not be implemented for every algebra.

EXAMPLES:

sage: R.<x> = ZZ['x']; A.<Sx> = OreAlgebra(R, 'Sx');
sage: ((x+5)*Sx - x).spread()
[4]
sage: ((x+5)*Sx - x).lclm((x+19)*Sx - x).spread()
[3, 4, 17, 18]
to_F(alg)

Returns a difference operator annihilating the coefficient sequence of every power series (about the origin) annihilated by self.

INPUT:

  • alg – the Ore algebra in which the output should be expressed. The algebra must satisfy alg.base_ring().base_ring() == self.base_ring().base_ring() and alg.is_F() is not False. Instead of an algebra object, also a string can be passed as argument. This amounts to specifying an Ore algebra over self.base_ring() with the forward difference with respect to self.base_ring().gen().

EXAMPLES:

sage: R.<x> = ZZ['x']
sage: A.<Dx> = OreAlgebra(R, 'Dx')
sage: R2.<n> = ZZ['n']
sage: A2.<Sn> = OreAlgebra(R2, 'Fn')
sage: (Dx - 1).to_F(A2)
(n + 1)*Fn + n
sage: ((1+x)*Dx^2 + Dx).to_F(A2)
(n^2 + n)*Fn + 2*n^2 + n
sage: ((x^3+x^2-x)*Dx + (x^2+1)).to_F(A2)
(-n - 1)*Fn^2 + (-n - 1)*Fn + n + 1
to_S(alg)

Returns a recurrence operator annihilating the coefficient sequence of every power series (about the origin) annihilated by self.

INPUT:

  • alg – the Ore algebra in which the output should be expressed. The algebra must satisfy alg.base_ring().base_ring() == self.base_ring().base_ring() and alg.is_S() is not False. Instead of an algebra object, also a string can be passed as argument. This amounts to specifying an Ore algebra over self.base_ring() with the standard shift with respect to self.base_ring().gen().

EXAMPLES:

sage: R.<x> = ZZ['x']
sage: A.<Dx> = OreAlgebra(R, 'Dx')
sage: R2.<n> = ZZ['n']
sage: A2.<Sn> = OreAlgebra(R2, 'Sn')
sage: (Dx - 1).to_S(A2)
(n + 1)*Sn - 1
sage: ((1+x)*Dx^2 + Dx).to_S(A2)
(n^2 + n)*Sn + n^2
sage: ((x^3+x^2-x)*Dx + (x^2+1)).to_S(A2)
(-n - 1)*Sn^2 + (n + 1)*Sn + n + 1
to_T(alg)

Rewrites self in terms of the eulerian derivation \(x*d/dx\).

If the base ring of the target algebra is not a field, the operator returned by the method may not correspond exactly to self, but only to a suitable left-multiple by a term \(x^k\).

INPUT:

  • alg – the Ore algebra in which the output should be expressed. The algebra must satisfy alg.base_ring().base_ring() == self.base_ring().base_ring() and alg.is_T() is not False. Instead of an algebra object, also a string can be passed as argument. This amounts to specifying an Ore algebra over self.base_ring() with an euler derivation with respect to self.base_ring().gen().

EXAMPLES:

sage: R.<x> = ZZ['x']
sage: R2.<y> = ZZ['y']
sage: A.<Dx> = OreAlgebra(R, 'Dx')
sage: (Dx^4).to_T(OreAlgebra(R2, 'Ty'))
Ty^4 - 6*Ty^3 + 11*Ty^2 - 6*Ty
sage: (Dx^4).to_T('Tx').to_D(A)
x^4*Dx^4
sage: _.to_T('Tx')
Tx^4 - 6*Tx^3 + 11*Tx^2 - 6*Tx
class ore_operator_1_1.UnivariateEulerDifferentialOperatorOverUnivariateRing(parent, *data, **kwargs)

Element of an Ore algebra K(x)[T], where T is the Euler differential operator T = x*d/dx

power_series_solutions(*args, **kwargs)

Computes the first few terms of the power series solutions of this operator.

The method raises an error if Sage does not know how to factor univariate polynomials over the base ring’s base ring.

The base ring has to have characteristic zero.

INPUT:

  • n – minimum number of terms to be computed

OUTPUT:

A list of power series of the form \(x^lpha + ...\) with pairwise distinct exponents \(lpha\) and coefficients in the base ring’s base ring’s fraction field. All expansions are computed up to order \(k\) where \(k\) is obtained by adding the maximal \(lpha\) to the maximum of \(n\) and the order of self.

EXAMPLES:

sage: R.<x> = ZZ['x']
sage: A.<Dx> = OreAlgebra(R, 'Dx')
sage: ((1-x)*Dx - 1).power_series_solutions(10) # geometric series
[1 + x + x^2 + x^3 + x^4 + x^5 + x^6 + x^7 + x^8 + O(x^9)]
sage: (Dx - 1).power_series_solutions(5) # exp(x)
[1 + x + 1/2*x^2 + 1/6*x^3 + O(x^4)]
sage: (Dx^2 - Dx + x).power_series_solutions(5) # a 2nd order equation
[x + 1/2*x^2 + 1/6*x^3 - 1/24*x^4 + O(x^5), 1 - 1/6*x^3 - 1/24*x^4 + O(x^5)]
sage: (2*x*Dx - 1).power_series_solutions(5) # sqrt(x) is not a power series
[]
spread(p=0)

Returns the spread of this operator.

This is the set of integers \(i\) such that sigma(self[0], i) and sigma(self[r], -r) have a nontrivial common factor, where sigma is the shift of the parent’s algebra and \(r\) is the order of self.

If the optional argument \(p\) is given, the method is applied to gcd(self[0], p) instead of self[0].

The output set contains \(\infty\) if the constant coefficient of self is zero.

This method is a stub and may not be implemented for every algebra.

EXAMPLES:

sage: R.<x> = ZZ['x']; A.<Sx> = OreAlgebra(R, 'Sx');
sage: ((x+5)*Sx - x).spread()
[4]
sage: ((x+5)*Sx - x).lclm((x+19)*Sx - x).spread()
[3, 4, 17, 18]
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
to_D(alg)

Returns the differential operator corresponding to self

INPUT:

  • alg – the Ore algebra in which the output should be expressed. The algebra must satisfy alg.base_ring().base_ring() == self.base_ring().base_ring() and alg.is_D() is not False. Instead of an algebra object, also a string can be passed as argument. This amounts to specifying an Ore algebra over self.base_ring() with the standard derivation with respect to self.base_ring().gen().

EXAMPLES:

sage: R.<x> = ZZ['x']
sage: A.<Tx> = OreAlgebra(R, 'Tx')
sage: (Tx^4).to_D(OreAlgebra(R, 'Dx'))
x^4*Dx^4 + 6*x^3*Dx^3 + 7*x^2*Dx^2 + x*Dx
sage: (Tx^4).to_D('Dx').to_T(A)
Tx^4
to_F(alg)

Returns a difference operator annihilating the coefficient sequence of every power series (about the origin) annihilated by self.

INPUT:

  • alg – the Ore algebra in which the output should be expressed. The algebra must satisfy alg.base_ring().base_ring() == self.base_ring().base_ring() and alg.is_F() is not False. Instead of an algebra object, also a string can be passed as argument. This amounts to specifying an Ore algebra over self.base_ring() with the forward difference with respect to self.base_ring().gen().

EXAMPLES:

sage: R.<x> = ZZ['x']
sage: A.<Tx> = OreAlgebra(R, 'Tx')
sage: R2.<n> = ZZ['n']
sage: A2.<Fn> = OreAlgebra(R2, 'Fn')
sage: (Tx - 1).to_F(A2)
n - 1
sage: ((1+x)*Tx^2 + Tx).to_F(A2)
(n^2 + 3*n + 2)*Fn + 2*n^2 + 3*n + 2
sage: ((x^3+x^2-x)*Tx + (x^2+1)).to_F(A2)
Fn^3 + (-n + 1)*Fn^2 + (-n + 1)*Fn + n + 1
to_S(alg)

Returns a recurrence operator annihilating the coefficient sequence of every power series (at the origin) annihilated by self.

INPUT:

  • alg – the Ore algebra in which the output should be expressed. The algebra must satisfy alg.base_ring().base_ring() == self.base_ring().base_ring() and alg.is_S() is not False. Instead of an algebra object, also a string can be passed as argument. This amounts to specifying an Ore algebra over self.base_ring() with the standard shift with respect to self.base_ring().gen().

EXAMPLES:

sage: R.<x> = ZZ['x']
sage: A.<Tx> = OreAlgebra(R, 'Tx')
sage: R2.<n> = ZZ['n']
sage: A2.<Sn> = OreAlgebra(R2, 'Sn')
sage: (Tx - 1).to_S(A2)
n - 1
sage: ((1+x)*Tx^2 + Tx).to_S(A2)
(n^2 + 3*n + 2)*Sn + n^2
sage: ((x^3+x^2-x)*Tx + (x^2+1)).to_S(A2)
Sn^3 + (-n - 2)*Sn^2 + (n + 2)*Sn + n
class ore_operator_1_1.UnivariateOreOperatorOverUnivariateRing(parent, *data, **kwargs)

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

associate_solutions(D, p)

If self is \(P\), this returns a list of pairs \((M, m)\) such that \(D*M = p + m*P\)

INPUT:

  • \(D\) – a first order operator with the same parent as self. Depending on the algebra, this operator may be constrained to certain choices. For example, for differential operators, it can only be \(D\) (corresponding to integration), and for recurrence operators, it can only be \(S - 1\) (corresponding to summation).
  • \(p\) – a nonzero base ring element

OUTPUT:

  • \(M\) – an operator of order self.order() - 1 with rational function coefficients.
  • \(m\) – a nonzero rational function.

Intended application: Express indefinite sums or integrals of holonomic functions in terms of the summand/integrand. For example, with \(D=S-1\) and \(P=S^2-S-1\) and \(p\) some polynomial, the output \(M\) is such that

\(\sum_{k=0}^n p(k) F_k = const + M(F_n)\)

where \(F_k\) denotes the Fibonacci sequence. The rational function \(m\) does not appear in the closed form, it can be regarded as a certificate.

The method returns the empty list if and only if no nontrivial solutions exist.

This function may not be implemented for every algebra.

EXAMPLES:

sage: R.<x> = QQ['x']; A.<Dx> = OreAlgebra(R, 'Dx');
sage: L = x*Dx^2 + Dx; p = 1  ## L(log(x)) == 0
sage: L.associate_solutions(Dx, p)
[(-x^2*Dx + x, -x)]
sage: (M, m) = _[0]
sage: Dx*M == p + m*L  ## this implies int(log(x)) == M(log(x)) = x*log(x) - x
True

sage: R.<x> = QQ['x']; A.<Dx> = OreAlgebra(R, 'Dx');
sage: L = x^2*Dx^2 + x*Dx + (x^2 - 1); p = 1  ## L(bessel(x)) == 0
sage: L.associate_solutions(Dx, p)
[(-Dx + 1/-x, 1/-x^2)]
sage: (M, m) = _[0]
sage: Dx*M == p + m*L  ## this implies int(bessel(x)) == -bessel'(x) -1/x*bessel(x)
True

sage: R.<n> = QQ['n']; A.<Sn> = OreAlgebra(R, 'Sn');
sage: L = Sn^2 - Sn - 1; p = 1  ## L(fib(n)) == 0
sage: L.associate_solutions(Sn - 1, p)
[(Sn, 1)]
sage: (M, m) = _[0]
sage: (Sn-1)*M == p + m*L  ## this implies sum(fib(n)) == fib(n+1)
True

sage: R.<n> = QQ['n']; A.<Sn> = OreAlgebra(R, 'Sn');
sage: L = Sn^3 - 2*Sn^2 - 2*Sn + 1; p = 1  ## L(fib(n)^2) == 0
sage: L.associate_solutions(Sn - 1, p)
[(1/2*Sn^2 - 1/2*Sn - 3/2, 1/2)]
sage: (Sn-1)*M == p + m*L  ## this implies sum(fib(n)^2) == 1/2*fib(n+2)^2 - 1/2*fib(n+1)^2 - 3/2*fib(n)^2
True
degree()

Returns the maximum degree among the coefficients of self

The degree of the zero operator is \(-1\).

If the base ring is not a polynomial ring, this causes an error.

desingularize(m=-1)

Returns a left multiple of self whose coefficients are polynomials and whose leading coefficient does not contain unnecessary factors.

INPUT:

  • \(m\) (optional) – If the order of self is \(r\), the output operator will have order \(r+m\). In order to ensure that all removable factors of the leading coefficient are removed in the output, \(m\) has to be chosen sufficiently large. If no \(m\) is given, a generic upper bound is determined. This feature may not be available for every class.

OUTPUT:

A left multiple of self whose coefficients are polynomials, whose order is \(m\) more than self, and whose leading coefficient has as low a degree as possible under these conditions.

The output is not unique. With low probability, the leading coefficient degree in the output may not be minimal.

EXAMPLES:

sage: R.<n> = ZZ['n']
sage: A.<Sn> = OreAlgebra(R, 'Sn')
sage: P = (-n^3 - 2*n^2 + 6*n + 9)*Sn^2 + (6*n^3 + 8*n^2 - 20*n - 30)*Sn - 8*n^3 - 12*n^2 + 20*n + 12
sage: Q = P.desingularize()
sage: Q.order()
4
sage: Q.leading_coefficient().degree()
1
sage: Q.leading_coefficient()
3114*n + 15570   # random
dispersion(p=0)

Returns the dispersion of this operator.

This is the maximum nonnegative integer \(i\) such that sigma(self[0], i) and sigma(self[r], -r) have a nontrivial common factor, where sigma is the shift of the parent’s algebra and \(r\) is the order of self.

An output \(-1\) indicates that there are no such integers \(i\) at all.

If the optional argument \(p\) is given, the method is applied to gcd(self[0], p) instead of self[0].

The output is \(\infty\) if the constant coefficient of self is zero.

EXAMPLES:

sage: R.<x> = ZZ['x']; A.<Sx> = OreAlgebra(R, 'Sx');
sage: ((x+5)*Sx - x).dispersion()
4
factor()

Returns a factorization of this operator into linear factors, if possible.

More precisely, the output will be a list \([L1,L2,...]\) of operators such that

  • \(L1*L2*...\) is equal to self
  • \(L2,L3,...\) are monic first order operators
  • \(L1\) has no first order right hand factor

This method requires the method right_factors() to be implemented.

EXAMPLE:

sage: R.<n> = ZZ['n']; A.<Sn> = OreAlgebra(R, 'Sn')
sage: L = (-2*n^4 - 17*n^3 - 45*n^2 - 33*n + 9)*Sn^3 + (6*n^4 + 57*n^3 + 168*n^2 + 148*n - 15)*Sn^2 + (-4*n^4 - 44*n^3 - 157*n^2 - 195*n - 38)*Sn + 4*n^3 + 34*n^2 + 80*n + 44
sage: L.factor()
[-2*n^4 - 17*n^3 - 45*n^2 - 33*n + 9,
Sn + (-4*n^5 - 44*n^4 - 171*n^3 - 295*n^2 - 230*n - 66)/(4*n^5 + 44*n^4 + 175*n^3 + 291*n^2 + 147*n - 45),
Sn + (2*n + 5)/(-2*n^2 - 7*n - 6),
Sn + (-2*n - 4)/(n + 1)]
sage: reduce(lambda p,q: p*q, _) - L
0
indicial_polynomial(p, var='alpha')

Computes the indicial polynomial of self at (a root of) p.

The indicial polynomial is a polynomial in the given variable var with coefficients in the fraction field of the base ring’s base ring.

The precise meaning of this polynomial may depend on the parent of self. A minimum requirement is that if self has a rational solution whose denominator contains sigma.factorial(p, e) but neither sigma(p, -1)*sigma.factorial(p, e) nor sigma.factorial(p, e + 1), then -e is a root of this polynomial.

Applied to \(p=1/x\), the maximum integer root of the output should serve as a degree bound for the polynomial solutions of self.

This method is a stub. Depending on the particular subclass, restrictions on p may apply.

polynomial_solutions(rhs=(), degree=None, solver=None)

Computes the polynomial solutions of this operator.

INPUT:

  • rhs (optional) – a list of base ring elements
  • degree (optional) – bound on the degree of interest.
  • solver (optional) – a callable for computing the right kernel of a matrix over the base ring’s base ring.

OUTPUT:

A list of tuples \((p, c_0,...,c_r)\) such that \(self(p) == c_0*rhs[0] + ... + c_r*rhs[r]\), where \(p\) is a polynomial and \(c_0,...,c_r\) are constants.

Note

  • Even if no rhs is given, the output will be a list of tuples [(p1,), (p2,),...] and not just a list of plain polynomials.
  • If no degree is given, a basis of all the polynomial solutions is returned. This feature may not be implemented for all algebras.

EXAMPLES:

sage: R.<n> = ZZ['n']; A.<Sn> = OreAlgebra(R, 'Sn')
sage: L = 2*Sn^2 + 3*(n-7)*Sn + 4
sage: L.polynomial_solutions((n^2+4*n-8, 4*n^2-5*n+3))
[(-70*n + 231, 242, -113)]
sage: L(-70*n + 231)
-210*n^2 + 1533*n - 2275
sage: 242*(n^2+4*n-8) - 113*(4*n^2-5*n+3)
-210*n^2 + 1533*n - 2275

sage: R.<x> = ZZ['x']; A.<Dx> = OreAlgebra(R, 'Dx')
sage: L = (x*Dx - 19).lclm( x*Dx - 4 )
sage: L.polynomial_solutions()
[(x^4,), (x^19,)]
rational_solutions(rhs=(), denominator=None, degree=None, solver=None)

Computes the rational solutions of this operator.

INPUT:

  • rhs (optional) – a list of base ring elements
  • denominator (optional) – bound on the degree of interest.
  • degree (optional) – bound on the degree of interest.
  • solver (optional) – a callable for computing the right kernel of a matrix over the base ring’s base ring.

OUTPUT:

A list of tuples \((r, c_0,...,c_r)\) such that \(self(r) == c_0*rhs[0] + ... + c_r*rhs[r]\), where \(r\) is a rational function and \(c_0,...,c_r\) are constants.

Note

  • Even if no rhs is given, the output will be a list of tuples [(p1,), (p2,),...] and not just a list of plain polynomials.
  • If no denominator is given, a basis of all the rational solutions is returned. This feature may not be implemented for all algebras.
  • If no degree is given, a basis of all the polynomial solutions is returned. This feature may not be implemented for all algebras.

EXAMPLES:

sage: R.<x> = ZZ['x']; A.<Dx> = OreAlgebra(R, 'Dx')
sage: L = ((x+3)*Dx + 2).lclm(x*Dx + 3).symmetric_product((x+4)*Dx-2)
sage: L.rational_solutions()
[((x^2 + 8*x + 16)/x^3,), ((x^2 + 8*x + 16)/(x^2 + 6*x + 9),)]
sage: L.rational_solutions((1, x))
[((7*x^5 + 21*x^4 + 73*x^2 + 168*x + 144)/(x^5 + 6*x^4 + 9*x^3), 5184, 756), ((4*x^2 + 14*x + 1)/(x^2 + 6*x + 9), 2592, 378), ((7*x^2 + 24*x)/(x^2 + 6*x + 9), 4608, 672)]
sage: L(_[0][0]) == _[0][1] + _[0][2]*x
True

sage: R.<n> = ZZ['n']; A.<Sn> = OreAlgebra(R, 'Sn');
sage: L = ((n+3)*Sn - n).lclm((2*n+5)*Sn - (2*n+1))
sage: L.rational_solutions()
[((-4*n^3 - 8*n^2 + 3)/(4*n^5 + 20*n^4 + 35*n^3 + 25*n^2 + 6*n),), (1/(4*n^2 + 8*n + 3),)]

sage: L = (2*n^2 - n - 2)*Sn^2 + (-n^2 - n - 1)*Sn + n^2 - 14
sage: y = (-n + 1)/(n^2 + 2*n - 2)
sage: L.rational_solutions((L(y),))
[((-n + 1)/(n^2 + 2*n - 2), 1)]          
spread(p=0)

Returns the spread of this operator.

This is the set of integers \(i\) such that sigma(self[0], i) and sigma(self[r], -r) have a nontrivial common factor, where sigma is the shift of the parent’s algebra and \(r\) is the order of self.

If the optional argument \(p\) is given, the method is applied to gcd(self[0], p) instead of self[0].

The output set contains \(\infty\) if the constant coefficient of self is zero.

This method is a stub and may not be implemented for every algebra.

EXAMPLES:

sage: R.<x> = ZZ['x']; A.<Sx> = OreAlgebra(R, 'Sx');
sage: ((x+5)*Sx - x).spread()
[4]
sage: ((x+5)*Sx - x).lclm((x+19)*Sx - x).spread()
[3, 4, 17, 18]
class ore_operator_1_1.UnivariateQDifferentialOperatorOverUnivariateRing(parent, *data, **kwargs)

Element of an Ore algebra K(x)[J], where J is the Jackson q-differentiation J f(x) = (f(q*x) - f(x))/(q*(x-1))

annihilator_of_integral()

Returns an operator \(L\) which annihilates all the indefinite \(q\)-integrals \(\int_q f\) where \(f\) runs through the functions annihilated by self. The output operator is not necessarily of smallest possible order.

EXAMPLES:

sage: R.<x> = ZZ['q'].fraction_field()['x']
sage: A.<Jx> = OreAlgebra(R, 'Jx')
sage: ((x-1)*Jx - 2*x).annihilator_of_integral()
(x - 1)*Jx^2 - 2*x*Jx
sage: _.annihilator_of_associate(Jx)
(x - 1)*Jx - 2*x
power_series_solutions(n=5)

Computes the first few terms of the power series solutions of this operator.

The method raises an error if Sage does not know how to factor univariate polynomials over the base ring’s base ring.

The base ring has to have characteristic zero.

INPUT:

  • n – minimum number of terms to be computed

OUTPUT:

A list of power series of the form \(x^lpha + ...\) with pairwise distinct exponents \(lpha\) and coefficients in the base ring’s base ring’s fraction field. All expansions are computed up to order \(k\) where \(k\) is obtained by adding the maximal \(lpha\) to the maximum of \(n\) and the order of self.

EXAMPLES:

sage: R.<x> = QQ['x']
sage: A.<Jx> = OreAlgebra(R, 'Jx', q=2)
sage: (Jx-1).lclm((1-x)*Jx-1).power_series_solutions()
[x^2 + x^3 + 3/5*x^4 + 11/35*x^5 + O(x^6), 1 + x - 2/7*x^3 - 62/315*x^4 - 146/1395*x^5 + O(x^6)]
spread(p=0)

Returns the spread of this operator.

This is the set of integers \(i\) such that sigma(self[0], i) and sigma(self[r], -r) have a nontrivial common factor, where sigma is the shift of the parent’s algebra and \(r\) is the order of self.

If the optional argument \(p\) is given, the method is applied to gcd(self[0], p) instead of self[0].

The output set contains \(\infty\) if the constant coefficient of self is zero.

This method is a stub and may not be implemented for every algebra.

EXAMPLES:

sage: R.<x> = ZZ['x']; A.<Sx> = OreAlgebra(R, 'Sx');
sage: ((x+5)*Sx - x).spread()
[4]
sage: ((x+5)*Sx - x).lclm((x+19)*Sx - x).spread()
[3, 4, 17, 18]
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
to_Q(alg)

Returns a q-recurrence operator which annihilates the coefficient sequence of every power series (about the origin) annihilated by self. The output operator may not be minimal.

INPUT:

  • alg – the Ore algebra in which the output should be expressed. The algebra must satisfy alg.base_ring().base_ring() == self.base_ring().base_ring() and alg.is_Q() == self.parent().is_J(). Instead of an algebra object, also a string can be passed as argument. This amounts to specifying an Ore algebra over self.base_ring() with the q-shift with respect to self.base_ring().gen().

EXAMPLES:

sage: Rn.<n> = ZZ['n']; Rx.<x> = ZZ['x']
sage: A.<Jx> = OreAlgebra(Rx, 'Jx', q=2)
sage: B.<Qn> = OreAlgebra(Rn, 'Qn', q=2)
sage: (Jx - 1).to_Q(B)
(2*n - 1)*Qn - 1
sage: ((x+1)*Jx - 1).to_Q(B)
(4*n - 1)*Qn^2 + (2*n - 2)*Qn
sage: (n*Qn-1).to_J(A).to_Q(B) % (n*Qn - 1)
0 
class ore_operator_1_1.UnivariateQRecurrenceOperatorOverUnivariateRing(parent, *data, **kwargs)

Element of an Ore algebra K(x)[S], where S is the shift x->q*x for some q in K.

annihilator_of_composition(a, solver=None)

Returns an operator \(L\) which annihilates all the sequences \(f(a(n))\) where \(f\) runs through the functions annihilated by self. The output operator is not necessarily of smallest possible order.

INPUT:

  • a – a polynomial \(u*x+v\) where \(x\) is the generator of the base ring, \(u\) and \(v\) are integers.
  • solver (optional) – a callable object which applied to a matrix with polynomial entries returns its kernel.

EXAMPLES:

sage: R.<x> = QQ['x']
sage: A.<Qx> = OreAlgebra(R, 'Qx', q=3)
sage: L = (x+3)*Qx^2 - (5*x+3)*Qx + 2*x-1
sage: data = L.to_list([1,2], 11)
sage: data
[1, 2, 15/4, 115/12, 1585/48, 19435/144, 2387975/4032, 188901875/70848, 488427432475/40336128, 1461633379710215/26500836096, 14580926901721431215/57983829378048]
sage: L2 = L.annihilator_of_composition(2*x)
sage: L2.to_list([1,15/4], 5)
[1, 15/4, 1585/48, 2387975/4032, 488427432475/40336128]
sage: Lrev = L.annihilator_of_composition(10 - x)
sage: Lrev.to_list([data[10], data[9]], 11)
[14580926901721431215/57983829378048, 1461633379710215/26500836096, 488427432475/40336128, 188901875/70848, 2387975/4032, 19435/144, 1585/48, 115/12, 15/4, 2, 1]
annihilator_of_sum()

Returns an operator \(L\) which annihilates all the indefinite sums \(\sum_{k=0}^n a_k\) where \(a_n\) runs through the sequences annihilated by self. The output operator is not necessarily of smallest possible order.

EXAMPLES:

sage: R.<x> = ZZ['q'].fraction_field()['x']
sage: A.<Qx> = OreAlgebra(R, 'Qx')
sage: ((x+1)*Qx - x).annihilator_of_sum()
(q*x + 1)*Qx - q*x
spread(p=0)

Returns the spread of this operator.

This is the set of integers \(i\) such that sigma(self[0], i) and sigma(self[r], -r) have a nontrivial common factor, where sigma is the shift of the parent’s algebra and \(r\) is the order of self.

If the optional argument \(p\) is given, the method is applied to gcd(self[0], p) instead of self[0].

The output set contains \(\infty\) if the constant coefficient of self is zero.

This method is a stub and may not be implemented for every algebra.

EXAMPLES:

sage: R.<x> = ZZ['x']; A.<Sx> = OreAlgebra(R, 'Sx');
sage: ((x+5)*Sx - x).spread()
[4]
sage: ((x+5)*Sx - x).lclm((x+19)*Sx - x).spread()
[3, 4, 17, 18]
to_J(alg)

Returns a q-differential operator which annihilates every power series (about the origin) whose coefficient sequence is annihilated by self. The output operator may not be minimal.

INPUT:

  • alg – the Ore algebra in which the output should be expressed. The algebra must satisfy alg.base_ring().base_ring() == self.base_ring().base_ring() and alg.is_J() is not False. Instead of an algebra object, also a string can be passed as argument. This amounts to specifying an Ore algebra over self.base_ring() with the q-derivation with respect to self.base_ring().gen().

EXAMPLES:

sage: Rn.<n> = ZZ['n']; Rx.<x> = ZZ['x']
sage: A.<Qn> = OreAlgebra(Rn, 'Qn', q=2)
sage: B.<Jx> = OreAlgebra(Rx, 'Jx', q=2)
sage: (Qn - 1).to_J(B)
(-2*x + 1)*Jx - 1
sage: ((n+1)*Qn - 1).to_J(B)
2*x*Jx^2 + (-4*x + 4)*Jx - 2
sage: (x*Jx-1).to_Q(A).to_J(B) % (x*Jx - 1)
0
to_list(init, n, start=0, append=False, padd=False)

Computes the terms of some sequence annihilated by self.

INPUT:

  • init – a vector (or list or tuple) of initial values. The components must be elements of self.base_ring().base_ring().fraction_field(). If the length is more than self.order(), we do not check whether the given terms are consistent with self.
  • n – desired number of terms.
  • start (optional) – index of the sequence term which is represented by the first entry of init. Defaults to zero.
  • append (optional) – if True, the computed terms are appended to init list. Otherwise (default), a new list is created.
  • padd (optional) – if True, the vector of initial values is implicitely prolonged to the left (!) by zeros if it is too short. Otherwise (default), the method raises a ValueError if init is too short.

OUTPUT:

A list of n terms whose \(k\) th component carries the sequence term with index start+k. Terms whose calculation causes an error are represented by None.

EXAMPLES:

sage: R.<x> = QQ['x']; A.<Qx> = OreAlgebra(R, 'Qx', q=3)
sage: (Qx^2-x*Qx + 1).to_list([1,1], 10)
[1, 1, 0, -1, -9, -242, -19593, -4760857, -3470645160, -7590296204063]
sage: (Qx^2-x*Qx + 1)(_)
[0, 0, 0, 0, 0, 0, 0, 0]
class ore_operator_1_1.UnivariateRecurrenceOperatorOverUnivariateRing(parent, *data, **kwargs)

Element of an Ore algebra K(x)[S], where S is the shift x->x+1.

annihilator_of_composition(a, solver=None)

Returns an operator \(L\) which annihilates all the sequences \(f(floor(a(n)))\) where \(f\) runs through the functions annihilated by self. The output operator is not necessarily of smallest possible order.

INPUT:

  • a – a polynomial \(u*x+v\) where \(x\) is the generator of the base ring, \(u\) and \(v\) are integers or rational numbers. If they are rational, the base ring of the parent of self must contain QQ.
  • solver (optional) – a callable object which applied to a matrix with polynomial entries returns its kernel.

EXAMPLES:

sage: R.<x> = QQ['x']
sage: A.<Sx> = OreAlgebra(R, 'Sx')
sage: ((2+x)*Sx^2-(2*x+3)*Sx+(x+1)).annihilator_of_composition(2*x+5) 
(16*x^3 + 188*x^2 + 730*x + 936)*Sx^2 + (-32*x^3 - 360*x^2 - 1340*x - 1650)*Sx + 16*x^3 + 172*x^2 + 610*x + 714
sage: ((2+x)*Sx^2-(2*x+3)*Sx+(x+1)).annihilator_of_composition(1/2*x)
(1/2*x^2 + 11/2*x + 15)*Sx^6 + (-3/2*x^2 - 25/2*x - 27)*Sx^4 + (3/2*x^2 + 17/2*x + 13)*Sx^2 - 1/2*x^2 - 3/2*x - 1
sage: ((2+x)*Sx^2-(2*x+3)*Sx+(x+1)).annihilator_of_composition(100-x)
(-x + 99)*Sx^2 + (2*x - 199)*Sx - x + 100
annihilator_of_interlacing(*other)

Returns an operator \(L\) which annihilates any sequence which can be obtained by interlacing sequences annihilated by self and the operators given in the arguments.

More precisely, if self and the operators given in the arguments are denoted \(L_1,L_2,\dots,L_m\), and if \(f_1(n),\dots,f_m(n)\) are some sequences such that \(L_i\) annihilates \(f_i(n)\), then the output operator \(L\) annihilates sequence \(f_1(0),f_2(0),\dots,f_m(0),f_1(1),f_2(1),\dots,f_m(1),\dots\), the interlacing sequence of \(f_1(n),\dots,f_m(n)\).

The output operator is not necessarily of smallest possible order.

The other operators must be coercible to the parent of self.

EXAMPLES:

sage: R.<x> = QQ['x']
sage: A.<Sx> = OreAlgebra(R, 'Sx')
sage: (x*Sx - (x+1)).annihilator_of_interlacing(Sx - (x+1), Sx + 1)
(-x^7 - 45/2*x^6 - 363/2*x^5 - 1129/2*x^4 - 45/2*x^3 + 5823/2*x^2 + 5751/2*x - 2349)*Sx^9 + (1/3*x^8 + 61/6*x^7 + 247/2*x^6 + 4573/6*x^5 + 14801/6*x^4 + 7173/2*x^3 + 519/2*x^2 - 3051*x + 756)*Sx^6 + (-7/2*x^6 - 165/2*x^5 - 1563/2*x^4 - 7331/2*x^3 - 16143/2*x^2 - 9297/2*x + 5535)*Sx^3 - 1/3*x^8 - 67/6*x^7 - 299/2*x^6 - 6157/6*x^5 - 22877/6*x^4 - 14549/2*x^3 - 10839/2*x^2 + 1278*x + 2430
annihilator_of_sum()

Returns an operator \(L\) which annihilates all the indefinite sums \(\sum_{k=0}^n a_k\) where \(a_n\) runs through the sequences annihilated by self. The output operator is not necessarily of smallest possible order.

EXAMPLES:

sage: R.<x> = ZZ['x']
sage: A.<Sx> = OreAlgebra(R, 'Sx')
sage: ((x+1)*Sx - x).annihilator_of_sum() # constructs L such that L(H_n) == 0
(x + 2)*Sx^2 + (-2*x - 3)*Sx + x + 1
forward_matrix_bsplit(n, start=0)

Uses division-free binary splitting to compute a product of n consecutive companion matrices of self.

If self annihilates some sequence \(c\) of order \(r\), this allows rapidly computing \(c_n, \ldots, c_{n+r-1}\) (or just \(c_n\)) without generating all the intermediate values.

INPUT:

  • n – desired number of terms to move forward
  • start (optional) – starting index. Defaults to zero.

OUTPUT:

A pair \((M, Q)\) where \(M\) is an \(r\) by \(r\) matrix and \(Q\) is a scalar, such that \(M / Q\) is the product of the companion matrix at \(n\) consecutive indices.

We have \(Q [c_{s+n}, \ldots, c_{s+r-1+n}]^T = M [c_s, c_{s+1}, \ldots, c_{s+r-1}]^T\), where \(s\) is the initial position given by start.

EXAMPLES:

sage: R = ZZ
sage: Rx.<x> = R[]
sage: Rxk.<k> = Rx[]
sage: Rxks = OreAlgebra(Rxk, 'Sk')
sage: ann = Rxks([1+k, -3*x - 2*k*x, 2+k])
sage: initial = Matrix([[1], [x]])
sage: M, Q = ann.forward_matrix_bsplit(5)
sage: (M * initial).change_ring(QQ['x']) / Q
[               63/8*x^5 - 35/4*x^3 + 15/8*x]
[231/16*x^6 - 315/16*x^4 + 105/16*x^2 - 5/16]

sage: Matrix([[legendre_P(5, x)], [legendre_P(6, x)]])
[               63/8*x^5 - 35/4*x^3 + 15/8*x]
[231/16*x^6 - 315/16*x^4 + 105/16*x^2 - 5/16]

TODO: this should detect if the base coefficient ring is QQ (etc.) and then switch to ZZ (etc.) internally.

forward_matrix_param_rectangular(value, n, start=0, m=None)

Assuming the coefficients of self are in \(R[x][k]\), computes the nth forward matrix with the parameter \(x\) evaluated at value, using rectangular splitting with a step size of \(m\).

TESTS:

sage: from sage.all import Matrix, randrange sage: R = ZZ sage: Rx = R[‘x’]; x = Rx.gen() sage: Rxk = Rx[‘k’]; k = Rxk.gen() sage: Rxks = OreAlgebra(Rxk, ‘Sk’) sage: V = QQ sage: Vks = OreAlgebra(V[‘k’], ‘Sk’) sage: for i in range(1000): sage: A = Rxks.random_element(randrange(1,4)) sage: r = A.order() sage: v = V.random_element() sage: initial = [V.random_element() for i in range(r)] sage: start = randrange(0,5) sage: n = randrange(0,30) sage: m = randrange(0,10) sage: B = Vks(list(A.polynomial()(x=v))) sage: singular = any([B[r](i) == 0 for i in range(n+r)]) sage: M, Q = A.forward_matrix_param_rectangular(v, n, m=m, start=start) sage: if Q == 0: sage: assert singular sage: else: sage: V1 = M * Matrix(initial).transpose() / Q sage: values = B.to_list(initial, n + r, start) sage: V2 = Matrix(values[-r:]).transpose() sage: if V1 != V2: sage: raise ValueError
generalized_series_solutions(n=5, dominant_only=False)

Returns the generalized series solutions of this operator.

These are solutions of the form

\((x/e)^{x u/v}\rho^x\exp\bigl(c_1 x^{1/m} +...+ c_{v-1} x^{1-1/m}\bigr)x^\alpha p(x^{-1/m},\log(x))\)

where

  • \(e\) is Euler’s constant (2.71...)
  • \(v\) is a positive integer
  • \(u\) is an integer; the term \((x/e)^(v/u)\) is called the “superexponential part” of the solution
  • \(\rho\) is an element of an algebraic extension of the coefficient field \(K\) (the algebra’s base ring’s base ring); the term \(\rho^x\) is called the “exponential part” of the solution
  • \(c_1,...,c_{v-1}\) are elements of \(K(\rho)\); the term \(\exp(...)\) is called the “subexponential part” of the solution
  • \(m\) is a positive integer multiple of \(v\), it is called the object’s “ramification”
  • \(\alpha\) is an element of some algebraic extension of \(K(\rho)\); the term \(n^\alpha\) is called the “polynomial part” of the solution (even if \(\alpha\) is not an integer)
  • \(p\) is an element of \(K(\rho)(\alpha)[[x]][y]\). It is called the “expansion part” of the solution.

An operator of order \(r\) has exactly \(r\) linearly independent solutions of this form. This method computes them all, unless the flags specified in the arguments rule out some of them.

Generalized series solutions are asymptotic expansions of sequences annihilated by the operator.

At present, the method only works for operators where \(K\) is some field which supports coercion to QQbar.

INPUT:

  • n (default: 5) – minimum number of terms in the expansions parts to be computed.
  • dominant_only (default: False) – if set to True, only compute solution(s) with maximal growth.

OUTPUT:

  • a list of DiscreteGeneralizedSeries objects forming a fundamental system for this operator.

EXAMPLES:

sage: R.<n> = QQ['n']; A.<Sn> = OreAlgebra(R, 'Sn')
sage: (Sn - (n+1)).generalized_series_solutions()
[(n/e)^n*n^(1/2)*(1 + 1/12*n^(-1) + 1/288*n^(-2) - 139/51840*n^(-3) - 571/2488320*n^(-4) + O(n^(-5)))]
sage: map(Sn - (n+1), _)
[0]

sage: L = ((n+1)*Sn - n).annihilator_of_sum().symmetric_power(2)
sage: L.generalized_series_solutions()
[1 + O(n^(-5)), (1 + O(n^(-5)))*log(n) + 1/2*n^(-1) - 1/12*n^(-2) + 1/120*n^(-4) + O(n^(-5)), (1 + O(n^(-5)))*log(n)^2 + (1/2*n^(-1) - 1/12*n^(-2) + 1/120*n^(-4) + O(n^(-5)))*log(n) - 3/2*n^(-1) + 5/8*n^(-2) - 7/54*n^(-3) - 59/3840*n^(-4) + O(n^(-5))]
sage: map(L, _)
[0, 0, 0]

sage: L = n^2*(1-2*Sn+Sn^2) + (n+1)*(1+Sn+Sn^2)
sage: L.generalized_series_solutions()
[exp(3.464101615137755?*I*n^(1/2))*n^(1/4)*(1 - 2.056810333988042?*I*n^(-1/2) - 1107/512*n^(-2/2) + (0.?e-19 + 1.489453749877895?*I)*n^(-3/2) + 2960239/2621440*n^(-4/2) + (0.?e-19 - 0.926161373412572?*I)*n^(-5/2) - 16615014713/46976204800*n^(-6/2) + (0.?e-20 + 0.03266142931818572?*I)*n^(-7/2) + 16652086533741/96207267430400*n^(-8/2) + (0.?e-20 - 0.1615093987591473?*I)*n^(-9/2) + O(n^(-10/2))), exp(-3.464101615137755?*I*n^(1/2))*n^(1/4)*(1 + 2.056810333988042?*I*n^(-1/2) - 1107/512*n^(-2/2) + (0.?e-19 - 1.489453749877895?*I)*n^(-3/2) + 2960239/2621440*n^(-4/2) + (0.?e-19 + 0.926161373412572?*I)*n^(-5/2) - 16615014713/46976204800*n^(-6/2) + (0.?e-20 - 0.03266142931818572?*I)*n^(-7/2) + 16652086533741/96207267430400*n^(-8/2) + (0.?e-20 + 0.1615093987591473?*I)*n^(-9/2) + O(n^(-10/2)))]

sage: L = guess([(-3)^n*(n+1)/(2*n+4) - 2^n*n^3/(n+3) for n in xrange(500)], A)
sage: L.generalized_series_solutions()
[2^n*n^2*(1 - 3*n^(-1) + 9*n^(-2) - 27*n^(-3) + 81*n^(-4) + O(n^(-5))), (-3)^n*(1 - n^(-1) + 2*n^(-2) - 4*n^(-3) + 8*n^(-4) + O(n^(-5)))]
sage: L.generalized_series_solutions(dominant_only=True)
[(-3)^n*(1 - n^(-1) + 2*n^(-2) - 4*n^(-3) + 8*n^(-4) + O(n^(-5)))]
right_factors(early_termination=False)

Returns a list of first-order right hand factors of this operator.

INPUT:

  • early_termination (optional) – if set to True, the search for factors will be aborted as soon as one factor has been found. A list containing this single factor will be returned (or the empty list if there are no first order factors). If set to False (default), a complete list will be computed.

OUTPUT:

A list of first order operators living in the parent of self of which self is a left multiple.

Note that this implementation does not construct factors that involve algebraic extensions of the constant field.

EXAMPLES:

sage: R.<n> = ZZ['n']; A.<Sn> = OreAlgebra(R, 'Sn');
sage: L = (-25*n^6 - 180*n^5 - 584*n^4 - 1136*n^3 - 1351*n^2 - 860*n - 220)*Sn^2 + (50*n^6 + 560*n^5 + 2348*n^4 + 5368*n^3 + 7012*n^2 + 4772*n + 1298)*Sn - 200*n^5 - 1540*n^4 - 5152*n^3 - 8840*n^2 - 7184*n - 1936
sage: L.right_factors()
[(n^2 + 2*n + 1)*Sn - 4*n - 2, (-25*n^2 - 30*n - 9)*Sn + 50*n^2 + 160*n + 128]
spread(p=0)

Returns the spread of this operator.

This is the set of integers \(i\) such that sigma(self[0], i) and sigma(self[r], -r) have a nontrivial common factor, where sigma is the shift of the parent’s algebra and \(r\) is the order of self.

If the optional argument \(p\) is given, the method is applied to gcd(self[0], p) instead of self[0].

The output set contains \(\infty\) if the constant coefficient of self is zero.

EXAMPLES:

sage: R.<x> = ZZ['x']; A.<Sx> = OreAlgebra(R, 'Sx');
sage: ((x+5)*Sx - x).spread()
[4]
sage: ((x+5)*Sx - x).lclm((x+19)*Sx - x).spread()
[3, 4, 17, 18]
to_D(alg)

Returns a differential operator which annihilates every power series whose coefficient sequence is annihilated by self. The output operator may not be minimal.

INPUT:

  • alg – the Ore algebra in which the output should be expressed. The algebra must satisfy alg.base_ring().base_ring() == self.base_ring().base_ring() and alg.is_D() is not False. Instead of an algebra object, also a string can be passed as argument. This amounts to specifying an Ore algebra over self.base_ring() with the standard derivation with respect to self.base_ring().gen().

EXAMPLES:

sage: Rn.<n> = ZZ['n']; Rx.<x> = ZZ['x']
sage: A.<Sn> = OreAlgebra(Rn, 'Sn')
sage: B.<Dx> = OreAlgebra(Rx, 'Dx')
sage: (Sn - 1).to_D(B)
(-x + 1)*Dx - 1
sage: ((n+1)*Sn - 1).to_D(B)
x*Dx^2 + (-x + 1)*Dx - 1
sage: (x*Dx-1).to_S(A).to_D(B)
x*Dx - 1
to_F(alg)

Returns the difference operator corresponding to self

INPUT:

  • alg – the Ore algebra in which the output should be expressed. The algebra must satisfy alg.base_ring().base_ring() == self.base_ring().base_ring() and alg.is_F() is not False. Instead of an algebra object, also a string can be passed as argument. This amounts to specifying an Ore algebra over self.base_ring() with the forward difference with respect to self.base_ring().gen().

EXAMPLES:

sage: R.<x> = ZZ['x']
sage: A.<Sx> = OreAlgebra(R, 'Sx')
sage: (Sx^4).to_F(OreAlgebra(R, 'Fx'))
Fx^4 + 4*Fx^3 + 6*Fx^2 + 4*Fx + 1
sage: (Sx^4).to_F('Fx').to_S(A)
Sx^4
to_T(alg)

Returns a differential operator, expressed in terms of the Euler derivation, which annihilates every power series (about the origin) whose coefficient sequence is annihilated by self. The output operator may not be minimal.

INPUT:

  • alg – the Ore algebra in which the output should be expressed. The algebra must satisfy alg.base_ring().base_ring() == self.base_ring().base_ring() and alg.is_T() is not False. Instead of an algebra object, also a string can be passed as argument. This amounts to specifying an Ore algebra over self.base_ring() with the Euler derivation with respect to self.base_ring().gen().

EXAMPLES:

sage: Rn.<n> = ZZ['n']; Rx.<x> = ZZ['x']
sage: A.<Sn> = OreAlgebra(Rn, 'Sn')
sage: B.<Tx> = OreAlgebra(Rx, 'Tx')
sage: (Sn - 1).to_T(B)
(-x + 1)*Tx - x
sage: ((n+1)*Sn - 1).to_T(B)
Tx^2 - x*Tx - x
sage: (x*Tx-1).to_S(A).to_T(B)
x*Tx^2 + (x - 1)*Tx
to_list(init, n, start=0, append=False, padd=False)

Computes the terms of some sequence annihilated by self.

INPUT:

  • init – a vector (or list or tuple) of initial values. The components must be elements of self.base_ring().base_ring().fraction_field(). If the length is more than self.order(), we do not check whether the given terms are consistent with self.
  • n – desired number of terms.
  • start (optional) – index of the sequence term which is represented by the first entry of init. Defaults to zero.
  • append (optional) – if True, the computed terms are appended to init list. Otherwise (default), a new list is created.
  • padd (optional) – if True, the vector of initial values is implicitely prolonged to the left (!) by zeros if it is too short. Otherwise (default), the method raises a ValueError if init is too short.

OUTPUT:

A list of n terms whose \(k\) th component carries the sequence term with index start+k. Terms whose calculation causes an error are represented by None.

EXAMPLES:

sage: R = ZZ['x']['n']; x = R('x'); n = R('n')
sage: A.<Sn> = OreAlgebra(R, 'Sn')
sage: L = ((n+2)*Sn^2 - x*(2*n+3)*Sn + (n+1))
sage: L.to_list([1, x], 5)
[1, x, (3*x^2 - 1)/2, (5*x^3 - 3*x)/2, (35*x^4 - 30*x^2 + 3)/8]
sage: polys = L.to_list([1], 5, padd=True)
[1, x, (3*x^2 - 1)/2, (5*x^3 - 3*x)/2, (35*x^4 - 30*x^2 + 3)/8]
sage: L.to_list([polys[3], polys[4]], 8, start=3)
[(5*x^3 - 3*x)/2, (35*x^4 - 30*x^2 + 3)/8, (63*x^5 - 70*x^3 + 15*x)/8, (231*x^6 - 315*x^4 + 105*x^2 - 5)/16, (429*x^7 - 693*x^5 + 315*x^3 - 35*x)/16]

sage: ((n-5)*Sn - 1).to_list([1], 10)
[1, 1/-5, 1/20, 1/-60, 1/120, -1/120, None, None, None, None]

Previous topic

ore_operator

Next topic

guessing