DFiniteSequence and DFiniteSequenceRing
D-finite sequences
Sublasses rec_sequences.RecurrenceSequenceRing and defines
sequences satisfying a linear
recurrence equation with polynomial coefficients. Such a D-finite sequence
\(a(n)\) is defined by a recurrence
and initial values \(a(0),...,a(r-1)\).
This is just a wrapper of UnivariateDFiniteSequence from the package
ore_algebra.
The only difference is, that we make sure that the leading coefficient
of the recurrences have no zeros by shifting the recurrence appropriately
(which might increase the order of the recurrence).
EXAMPLES:
sage: from rec_sequences.DFiniteSequenceRing import *
sage: R.<n> = PolynomialRing(QQ)
sage: D = DFiniteSequenceRing(R) # create D-finite sequence ring over QQ
sage: fac = D([n+1,-1],[1]) # define factorials
sage: fac[:10]
[1, 1, 2, 6, 24, 120, 720, 5040, 40320, 362880]
sage: # define harmonic numbers using guessing
sage: harm_terms = [sum(1/i for i in range(1,k)) for k in range(1,20)]
sage: harm = D(harm_terms)
sage: harm
D-finite sequence a(n): (-n - 1)*a(n) + (2*n + 3)*a(n+1) + (-n - 2)*a(n+2)
= 0 and a(0)=0 , a(1)=1
sage: a = fac+harm
sage: a.order(), a.degree()
(3, 4)
sage: harm.is_eventually_positive() # harm. numbers positive from term 1 on
1
- class rec_sequences.DFiniteSequenceRing.DFiniteSequence(parent, coefficients, initial_values, name='a', is_gen=False, construct=False, cache=True)
Bases:
rec_sequences.RecurrenceSequenceRing.RecurrenceSequenceElementA D-finite sequence, i.e. a sequence where every term can be determined by a linear recurrence with polynomial coefficients and finitely many initial values. We assume that this recurrence holds for all values.
- __call__(expr)
Returns the sequence
self[expr(n)]if expr is a symbolic expression in one variable representing a linear polynomial.INPUT:
expr– a linear rational univariate polynomial (can be in the symbolic ring)
OUTPUT:
The sequence
self[expr(n)].EXAMPLES:
sage: from rec_sequences.DFiniteSequenceRing import * sage: R.<n> = PolynomialRing(QQ) sage: D = DFiniteSequenceRing(R) sage: fac = D([n+1,-1],[1]) sage: harm = D([-n-1,2*n+3,-n-2],[0,1]) sage: harm.subsequence(2) == harm(2*n) True
- __contains__(item)
Checks whether
itemis a term of the sequence.INPUT:
item– an element in the ground field of the sequence
OUTPUT:
Returns
Trueif there is a termitemin the sequence andFalseotherwise.EXAMPLES:
sage: from rec_sequences.DFiniteSequenceRing import * sage: R.<n> = PolynomialRing(QQ) sage: D = DFiniteSequenceRing(R) sage: fac = D([n+1,-1],[1]) sage: factorial(5) in fac # long time True sage: from rec_sequences.CFiniteSequenceRing import * sage: C = CFiniteSequenceRing(QQ) sage: 2 in C(10*[1,-1]) False sage: 4 in C.an_element() False sage: 13 in C.an_element() True
- __init__(parent, coefficients, initial_values, name='a', is_gen=False, construct=False, cache=True)
Construct a D-finite sequence \(a(n)\) with recurrence
\[p_0(n) a(n) + \dots + p_r(n) a(n+r) = 0 \text{ for all } n \geq 0\]from given list of coefficients \(p_0, ... , p_r\) and given list of initial values \(a(0), ..., a(r-1)\).
INPUT:
parent– aDFiniteSequenceRingcoefficients– the coefficients of the recurrenceinitial_values– a list of initial values, determining the sequence with at least order of the recurrence many valuesname(default “a”) – a name for the sequence
OUTPUT:
A D-finite sequence determined by the given recurrence and initial values.
EXAMPLES:
sage: from rec_sequences.DFiniteSequenceRing import * sage: R.<n> = PolynomialRing(QQ) sage: D = DFiniteSequenceRing(R) sage: fac = D([n+1,-1],[1]) sage: fac D-finite sequence a(n): (n + 1)*a(n) + (-1)*a(n+1) = 0 and a(0)=1
- __invert__()
Tries to compute the multiplicative inverse, if this is possible, by guessing the inverse using
100terms. Such an inverse exists if and only if the sequence is an interlacing of hypergeometric sequences.The method can be called by ~self or self.inverse_of_unit().
OUTPUT:
The multiplicative inverse of the sequence if it exists. Raises an
ValueErrorif the sequence is not invertible or the inverse could not be found.EXAMPLES:
sage: from rec_sequences.DFiniteSequenceRing import * sage: R.<n> = PolynomialRing(QQ) sage: D = DFiniteSequenceRing(R) sage: ~D([n+1,-1],[1]) D-finite sequence a(n): (1)*a(n) + (-n - 1)*a(n+1) = 0 and a(0)=1 sage: ~D([-n-1,2*n+3,-n-2],[0,1]) Traceback (most recent call last): ... ValueError: Sequence is not invertible
- __le__(other)
Checks whether
selfis less than or equal toothertermwise, i.e., checks whetherself[n] <= other[n]for all natural numbersncan be proven. Only returnsTrueif the inequality can be proven.Note
Since the inequality is checked termwise, this method is not equivalent to
!(self > other)or toself < other or self==other.ALGORITHM:
The first 20 initial values are used to check whether the inequality can be falsified. Then,
is_positive()is used using the amount of time specified in the initialization of the parent class (default: 2 seconds).INPUT:
other– a D-finite sequence.
OUTPUT:
Trueifself[n] <= other[n]for all natural numbersnandFalseotherwise. Raises aValueErrorif neither could be shown.EXAMPLES:
sage: from rec_sequences.DFiniteSequenceRing import * sage: from rec_sequences.CFiniteSequenceRing import * sage: R.<n> = PolynomialRing(QQ) sage: D = DFiniteSequenceRing(R) sage: C = CFiniteSequenceRing(QQ, time_limit = 10) sage: harm = D([-n-1,2*n+3,-n-2],[0,1]) sage: harm >= 1 False sage: fac = D([n+1,-1],[1]) sage: 0 <= fac True sage: fib = C([1,1,-1], [1,1]) sage: luc = C([1,1,-1], [2,1]) sage: luc >= fib True
- __lt__(other)
Checks whether
selfis less thanothertermwise, i.e., checks whetherself[n] < other[n]for all natural numbersncan be proven. Only returnsTrueif the inequality can be proven.Note
Since the inequality is checked termwise, this method is not equivalent to
!(self >= other).ALGORITHM:
The first 20 initial values are used to check whether the inequality can be falsified. Then,
is_positive()is used using the amount of time specified in the initialization of the parent class (default: 2 seconds).INPUT:
other– a D-finite sequence.
OUTPUT:
Trueifself[n] < other[n]for all natural numbersnandFalseotherwise. Raises aValueErrorif neither could be shown.EXAMPLES:
sage: from rec_sequences.DFiniteSequenceRing import * sage: from rec_sequences.CFiniteSequenceRing import * sage: R.<n> = PolynomialRing(QQ) sage: D = DFiniteSequenceRing(R) sage: C = CFiniteSequenceRing(QQ, time_limit = 10) sage: harm = D([-n-1,2*n+3,-n-2],[0,1]) sage: 0 < harm False sage: fac = D([n+1,-1],[1]) sage: fac > 0 True sage: harm < fac False sage: fib = C([1,1,-1], [0,1]) sage: c = C([2,-1], [1]) sage: c > fib True
- _add_(right)
Return the sum of
selfandright.INPUTS:
right– D-finite sequences over the same D-finite sequence ring asself.
OUTPUTS:
The addition of
selfwithright.EXAMPLES:
sage: from rec_sequences.DFiniteSequenceRing import * sage: R.<n> = PolynomialRing(QQ) sage: D = DFiniteSequenceRing(R) sage: fac = D([n+1,-1],[1]) sage: harm = D([-n-1,2*n+3,-n-2],[0,1]) sage: (fac+harm)[:8] [1, 2, 7/2, 47/6, 313/12, 7337/60, 14449/20, 705963/140]
- _latex_(name=None)
Not yet implemented.
- _mul_(right)
Return the product of
selfandright. The result is the termwise product (Hadamard product) ofselfandright. To get the cauchy product use the methodcauchy().INPUTS:
right– D-finite sequences over the same D-finite sequence ring asself
OUTPUTS:
The product of
selfwithright.EXAMPLES:
sage: from rec_sequences.DFiniteSequenceRing import * sage: R.<n> = PolynomialRing(QQ) sage: D = DFiniteSequenceRing(R) sage: fac = D([n+1,-1],[1]) sage: harm = D([-n-1,2*n+3,-n-2],[0,1]) sage: (fac*harm)[:10] [0, 1, 3, 11, 50, 274, 1764, 13068, 109584, 1026576]
- _repr_(name=None)
Produces a string representation of the sequence.
INPUT:
name(optional) – a string used as the name of the sequence; if not given,self.name()is used.
OUTPUT:
A string representation of the sequence consisting of the recurrence and enough initial values to uniquely define the sequence.
EXAMPLES:
sage: from rec_sequences.DFiniteSequenceRing import * sage: R.<n> = PolynomialRing(QQ) sage: D = DFiniteSequenceRing(R) sage: D([-n-1,2*n+3,-n-2],[0,1],name="harm") D-finite sequence harm(n): (-n - 1)*harm(n) + (2*n + 3)*harm(n+1) + (-n - 2)*harm(n+2) = 0 and harm(0)=0 , harm(1)=1
- ann()
Alias of
recurrence().
- cauchy(right)
Computes the Cauchy product of two sequences.
INPUT:
right– other D-finite sequences over the same D-finite sequence ring
OUTPUT:
The cauchy product of
selfwithright.EXAMPLES:
sage: from rec_sequences.DFiniteSequenceRing import * sage: R.<n> = PolynomialRing(QQ) sage: D = DFiniteSequenceRing(R) sage: fac = D([n+1,-1],[1]) sage: harm = D([-n-1,2*n+3,-n-2],[0,1]) sage: cauchy = harm.cauchy(fac) sage: cauchy.order(), cauchy.degree() (16, 22)
- coefficients()
Returns the list of polynomial coefficients of the recurrence of
selfin the format[p0,...,pr]representing the recurrence\[p_0(n) a(n) + \dots + p_r(n) a(n+r) = 0.\]OUTPUT:
The coefficients of the recurrence of the sequence.
EXAMPLES:
sage: from rec_sequences.DFiniteSequenceRing import * sage: R.<n> = PolynomialRing(QQ) sage: D = DFiniteSequenceRing(R) sage: harm = D([-n-1,2*n+3,-n-2],[0,1]) sage: harm.coefficients() [-n - 1, 2*n + 3, -n - 2]
- companion_matrix()
Returns the \(r \times r\) companion matrix
\[\begin{split}\begin{pmatrix} 0 & 0 & \dots & 0 & -p_0/p_r \\ 1 & 0 & \dots & 0 & -p_1/p_r \\ 0 & 1 & \dots & 0 & -p_2/p_r \\ \vdots & \vdots & \ddots & \vdots & \vdots \\ 0 & 0 & \dots & 1 & -p_{r-1}/p_r \end{pmatrix} .\end{split}\]of
selfwith entries in the fraction field of the base.OUTPUT:
The companion matrix.
EXAMPLES:
sage: from rec_sequences.DFiniteSequenceRing import * sage: R.<n> = PolynomialRing(QQ) sage: D = DFiniteSequenceRing(R) sage: harm = D([-n-1,2*n+3,-n-2],[0,1]) sage: harm.companion_matrix() [ 0 (-n - 1)/(n + 2)] [ 1 (2*n + 3)/(n + 2)]
- compress(proof=False)
Tries to compress the sequence
selfas much as possible by trying to find a smaller recurrence.INPUT:
proof(default:False) ifTrue, then the result is guaranteed to be true. Otherwise it can happen, although unlikely, that the sequences are different.
OUTPUT:
A sequence which is equal to
selfbut may consist of a smaller operator (in terms of the order). In the worst case if no compression is possibleselfis returned.
- degree()
The maximal degree of the coefficient polynomials.
OUTPUT:
Returns the degree of the sequence.
EXAMPLES:
sage: from rec_sequences.DFiniteSequenceRing import * sage: R.<n> = PolynomialRing(QQ) sage: D = DFiniteSequenceRing(R) sage: fac = D([n+1,-1],[1]) sage: fac.degree() 1
- divides(right, bound=100, divisor=False)
Checks whether
selfdividesrightin the sequence ring, i.e. checks whether there is a sequencedivin this ring such thatright/self == div.Warning
We assume that
selfdoes not contain any zero terms.INPUT:
right– a sequence in the same ring asself.bound(default:100) – maximal number of terms used to guess the divisordivisor(default:False) – ifTrue, then the divisor is returned instead ofTrueif it could be found.
OUTPUT:
Returns
Trueifselfdividesrightcould be proven. If it could not be proven,Falseis returned. IfdivisorisTrue, thendivis returned instead ofTrue.EXAMPLES:
sage: from rec_sequences.DFiniteSequenceRing import * sage: R.<n> = PolynomialRing(QQ) sage: D = DFiniteSequenceRing(R) sage: D([n+1,n^2-7,n+2], [3,1]).divides(D([n^2+3, n+1], [7])) False sage: fac = D([n+1,-1],[1]) sage: harm = D([-n-1,2*n+3,-n-2],[0,1]) sage: div = fac.divides(harm, divisor=True) sage: div*fac == harm True
- get_ore_algebra_sequence()
Returns the underlying
UnivariateDFiniteSequence.OUTPUT:
A
UnivariateDFiniteSequencedescribing the same sequence.See also
Check the documentation of
UnivariateDFiniteSequencein theore_algebrapackageEXAMPLES:
sage: from rec_sequences.DFiniteSequenceRing import * sage: R.<n> = PolynomialRing(QQ) sage: D = DFiniteSequenceRing(R) sage: fac = D([n+1,-1],[1]) sage: fac.get_ore_algebra_sequence() Univariate D-finite sequence defined by the annihilating operator -Sn + n + 1 and the initial conditions {0: 1}
- get_shift(i)
Returns a vector v, s.t.
self[n+i] = v[0]*self[n]+...+v[r-1]*self[n+r-1]for allnifselfhas orderr.INPUT: -
i– a non-negative integerOUTPUT: A vector describing the components of
self[n+i]w.r.t. the generating systemself[n],...,self[n+r-1].EXAMPLES:
sage: from rec_sequences.DFiniteSequenceRing import * sage: R.<n> = PolynomialRing(QQ) sage: D = DFiniteSequenceRing(R) sage: harm = D([-n-1,2*n+3,-n-2],[0,1]) sage: v = harm.get_shift(3) sage: v ((-2*n^2 - 7*n - 5)/(n^2 + 5*n + 6), (3*n^2 + 12*n + 11)/(n^2 + 5*n + 6)) sage: all([harm[n+3]==sum(v[i](n)*harm[n+i] for i in [0,1]) ....: for n in range(10)]) True
- has_no_zeros(bound=0, time=- 1, bound_n=5)
Tries to prove that the sequence has no zeros. This is done using different algorithms (if
timeis specified,time/2is used for each of the algorithm):Determine the sign pattern using
sign_pattern()and check whether it contains zeros.Uses
is_eventually_positive()to show that the squared sequence is positive.
INPUT:
bound(default:0) – length of induction hypothesistime(default:-1) – if positive, this is the maximal time (in seconds) after computations are abortedbound_n(default:5) – index up to which it is checked whether the sequences is positive from that term on for the algorithms usingis_eventually_positive().
OUTPUT:
Returns
Trueif every term of the sequence is not equal to zero andFalseotherwise. Raises aTimeoutErrorif neither could be proven. If the formula for CAD is too big aRuntimeErrormight be raised.EXAMPLES:
sage: from rec_sequences.DFiniteSequenceRing import * sage: R.<n> = PolynomialRing(QQ) sage: D = DFiniteSequenceRing(R) sage: fac = D([n+1,-1],[1]) sage: fac.has_no_zeros() # long time True sage: harm = D([-n-1,2*n+3,-n-2],[0,1]) sage: harm.has_no_zeros() False sage: from rec_sequences.CFiniteSequenceRing import * sage: C = CFiniteSequenceRing(QQ) sage: n = var("n") sage: C(3^n-n*2^n).has_no_zeros(time=10) # long time True sage: C(10*[1,-1]).has_no_zeros(time=10) # long time True sage: C(n-4).has_no_zeros(time=10) # long time False
- interlace(*others)
Returns the interlaced sequence of self with
others.INPUT:
others– other D-finite sequences over the same D-finite sequence ring
OUTPUT:
The interlaced sequence of self with
others.EXAMPLES:
sage: from rec_sequences.DFiniteSequenceRing import * sage: R.<n> = PolynomialRing(QQ) sage: D = DFiniteSequenceRing(R) sage: fac = D([n+1,-1],[1]) sage: harm = D([-n-1,2*n+3,-n-2],[0,1]) sage: fac.interlace(harm)[:10] [1, 0, 1, 1, 2, 3/2, 6, 11/6, 24, 25/12]
- is_eventually_positive(bound_n=5, bound=0, strict=True, time=- 1)
Uses the Gerhold-Kauers methods (Algorithm 2 in [KP10]) to check whether the sequence is eventually positive. This is done by checking whether
self.shift(n)is positive forn <= bound_n. The smallest such index is returned. For everynthe giventime(if given) andboundis used.INPUT:
bound_n(default:5) – index up to which it is checked whether the sequences is positive from that term on.bound(default:0) – length of induction hypothesisstrict(default:True) – ifFalsenon-negativity instead of positivity is checkedtime(default:-1) – if positive, this is the maximal time (in seconds) after computations are aborted
OUTPUT:
Returns
nfrom where the sequence on is positive. If no such index could be found, aValueErroris raised.EXAMPLES:
sage: from rec_sequences.DFiniteSequenceRing import * sage: R.<n> = PolynomialRing(QQ) sage: D = DFiniteSequenceRing(R) sage: harm = D([-n-1,2*n+3,-n-2],[0,1]) sage: harm.is_eventually_positive() 1 sage: fac = D([n+1,-1],[1]) sage: fac.is_eventually_positive() 0
- is_positive(bound=0, strict=True, time=- 1)
Uses the Gerhold-Kauers methods (Algorithm 2 in [KP10]) to check whether the sequence is positive.
INPUT:
bound(default:0) – length of induction hypothesisstrict(default:True) – ifFalsenon-negativity instead of positivity is checkedtime(default:-1) – if positive, this is the maximal time (in seconds) after computations are aborted
OUTPUT:
Returns
Trueif it is positive,Falseif it is not positive and raises aValueErrorexception if it could neither prove or disprove positivity. If the time runs out, aTimeoutErroris raised.EXAMPLES:
sage: from rec_sequences.DFiniteSequenceRing import * sage: R.<n> = PolynomialRing(QQ) sage: D = DFiniteSequenceRing(R) sage: harm = D([-n-1,2*n+3,-n-2],[0,1]) sage: harm.is_positive() False sage: fac = D([n+1,-1],[1]) sage: fac.is_positive() True
- log = <Logger DFin (WARNING)>
- prepend(values)
Prepends the given values to the sequence.
Input -
values– list of values in the base ringOUTPUT:
A sequence having the same terms with the additional
valuesat the beginning.EXAMPLES:
sage: from rec_sequences.DFiniteSequenceRing import * sage: R.<n> = PolynomialRing(QQ) sage: D = DFiniteSequenceRing(R) sage: fac = D([n+1,-1],[1]) sage: fac.prepend([-17,32,5])[:10] [-17, 32, 5, 1, 1, 2, 6, 24, 120, 720]
- recurrence()
The annihilating operator of
selfas anOreOperator.OUTPUT:
Annihilating
OreOperatorofself.EXAMPLES:
sage: from rec_sequences.DFiniteSequenceRing import * sage: R.<n> = PolynomialRing(QQ) sage: D = DFiniteSequenceRing(R) sage: fac = D([n+1,-1],[1]) sage: fac.recurrence() -Sn + n + 1
- shift(k=1)
Shifts
selfk-times.INPUT:
k(default:1) – an integer
OUTPUT:
The sequence \((a(n+k))_{k \in \mathbb{N}}\).
EXAMPLES:
sage: from rec_sequences.DFiniteSequenceRing import * sage: R.<n> = PolynomialRing(QQ) sage: D = DFiniteSequenceRing(R) sage: a = D([n+1,-1],[1]).shift(2) sage: a[:10] [2, 6, 24, 120, 720, 5040, 40320, 362880, 3628800, 39916800]
- sign_pattern(bound=0, time=- 1, data=100)
Suppose that
selfhas a sign pattern which is cyclic. We try to guess this pattern and then verify it using the Gerhold-Kauers [KP10] method.INPUT:
bound(default:0) – length of induction hypothesis.time(default:-1) – if positive, this is the maximal time (in seconds) used to verify the sign patterns.data(default:100) – number of terms used to guess the sign-pattern.
OUTPUT:
The sign pattern of the sequence as an object of type
rec_sequences.SignPattern. If no pattern could be guessed or this pattern could not be verified, aValueErroris raised.EXAMPLES:
sage: from rec_sequences.DFiniteSequenceRing import * sage: R.<n> = PolynomialRing(QQ) sage: D = DFiniteSequenceRing(R) sage: harm = D([-n-1,2*n+3,-n-2],[0,1]) sage: harm.sign_pattern() Sign pattern: initial values <0> cycle <+> sage: n = var("n") sage: D(3^n).interlace(D(-2^n)).prepend([1,3,-2]).sign_pattern() Sign pattern: initial values <+> cycle <+->
- subsequence(u, v=0)
Returns the sequence
self[floor(u*n+v)].INPUT:
u– a rational numberv(optional) – a rational number
OUTPUT:
The sequence
self[floor(u*n+v)].EXAMPLES:
sage: from rec_sequences.DFiniteSequenceRing import * sage: R.<n> = PolynomialRing(QQ) sage: D = DFiniteSequenceRing(R) sage: fac = D([n+1,-1],[1]) sage: harm = D([-n-1,2*n+3,-n-2],[0,1]) sage: harm[:10] [0, 1, 3/2, 11/6, 25/12, 137/60, 49/20, 363/140, 761/280, 7129/2520] sage: harm.subsequence(2)[:5] [0, 3/2, 25/12, 49/20, 761/280]
- sum()
Returns the sequence \(\sum_{i=0}^n c(i)\), the sequence describing the partial sums.
OUTPUT:
The D-finite sequence \(\sum_{i=0}^n c(i)\).
EXAMPLES:
sage: from rec_sequences.DFiniteSequenceRing import * sage: R.<n> = PolynomialRing(QQ) sage: D = DFiniteSequenceRing(R) sage: fac = D([n+1,-1],[1]) sage: harm = D([-n-1,2*n+3,-n-2],[0,1]) sage: harm[:10] [0, 1, 3/2, 11/6, 25/12, 137/60, 49/20, 363/140, 761/280, 7129/2520] sage: harm.sum()[:10] [0, 1, 5/2, 13/3, 77/12, 87/10, 223/20, 481/35, 4609/280, 4861/252]
- zeros(bound=0, time=- 1, data=100)
Computes the zeros of the sequence provided that the sequence satisfies the Skolem-Mahler-Lech-Theorem, i.e., the zeros consist of finitely many zeros together with a finite number of arithmetic progressions. The method
sign_pattern()is used to derive the sign pattern from which the zeros are extracted. All the parameters correspond to the parameters ofsign_pattern().INPUT:
bound(default:0) – length of induction hypothesis.time(default:-1) – if positive, this is the maximal time (in seconds) used to verify the sign patterns.data(default:100) – number of terms used to guess the sign-pattern.
OUTPUT:
The zero pattern of the sequence as an object of type
rec_sequences.ZeroPattern. If no pattern could be guessed or this pattern could not be verified, aValueErroris raised.EXAMPLES:
sage: from rec_sequences.CFiniteSequenceRing import * sage: from rec_sequences.DFiniteSequenceRing import * sage: R.<n> = PolynomialRing(QQ) sage: D = DFiniteSequenceRing(R) sage: harm = D([-n-1,2*n+3,-n-2],[0,1]) sage: harm.zeros() Zero pattern with finite set {0} and no arithmetic progressions sage: n = var("n") sage: D(3^n).interlace(D(-2^n)).prepend([1,3,-2]).zeros() Zero pattern with finite set {} and no arithmetic progressions sage: C = CFiniteSequenceRing(QQ) sage: C(10*[1,0,-1,0]).prepend([0,0,2]).zeros() # random Zero pattern with finite set {0, 1} and arithmetic progressions: - Arithmetic progression (4*n+6)_n - Arithmetic progression (4*n+4)_n
- class rec_sequences.DFiniteSequenceRing.DFiniteSequenceRing(ring, time_limit=2, name=None, element_class=None, category=None)
Bases:
rec_sequences.RecurrenceSequenceRing.RecurrenceSequenceRingA Ring of D-finite sequences over a field.
- Element
- __init__(ring, time_limit=2, name=None, element_class=None, category=None)
Constructor for a D-finite sequence ring.
INPUT:
ring– a polynomial ring over a field containing the coefficients of the recurrences.time_limit(default:2) – a positive number indicating the time limit in seconds used to prove inequalities.
OUTPUT:
A ring of D-finite sequences over the given ring.
EXAMPLES:
sage: from rec_sequences.DFiniteSequenceRing import * sage: R.<n> = PolynomialRing(QQ) sage: DFiniteSequenceRing(R) Ring of D-finite sequences over Rational Field sage: S.<n> = PolynomialRing(NumberField(x^2-2, "z")) sage: DFiniteSequenceRing(S) Ring of D-finite sequences over Number Field in z with defining polynomial x^2 - 2
- _element_constructor_(x, y=None, name='a', check=True, is_gen=False, construct=False, **kwds)
Tries to construct a sequence \(a(n)\).
This is possible if:
xis already a sequence in the right ring.xis aUnivariateDFiniteSequence.xis a list of field elements andyis a list of field elements. Thenxis interpreted as the coefficients of the sequence andyas the initial values of the sequence, i.e. \(a(0), ..., a(r-1)\).xcan be converted into a field element. Then it is interpreted as the constant sequence \((x)_{n \in \mathbb{N}}\)xis in the symbolic ring andya variable, then values ofx.subs(y=n)for integersnare created and a recurrence for this sequence is guessed. Ifyis not given, we try to extract it fromx.xis a list, then guessing on this list is used to determine a recurrence.xis a univariate polynomial, then the sequence represents the polynomial sequence.xis aRecurrenceSequence, the firstyterms are used to guess a D-finite sequence (ifyis not given,100terms are used).
EXAMPLES:
sage: from rec_sequences.DFiniteSequenceRing import * sage: R.<n> = PolynomialRing(QQ) sage: D = DFiniteSequenceRing(R) sage: fac = D([n+1,-1],[1]) sage: fac D-finite sequence a(n): (n + 1)*a(n) + (-1)*a(n+1) = 0 and a(0)=1 sage: fac==D(fac) True sage: D(1/7)[:5] [1/7, 1/7, 1/7, 1/7, 1/7] sage: n = var("n") sage: D(2^n)[:5] [1, 2, 4, 8, 16] sage: D(10*[1,0]) D-finite sequence a(n): (1)*a(n) + (-1)*a(n+2) = 0 and a(0)=1 , a(1)=0
- _latex_()
OUTPUT:
A latex representation of the D-finite sequence ring.
EXAMPLES:
sage: from rec_sequences.DFiniteSequenceRing import * sage: R.<n> = PolynomialRing(QQ) sage: print(latex(DFiniteSequenceRing(R))) \mathcal{D}(\Bold{Q}) sage: S.<n> = PolynomialRing(NumberField(x^2-2, "z")) sage: print(latex(DFiniteSequenceRing(S))) \mathcal{D}(\Bold{Q}[z]/(z^{2} - 2))
- _repr_()
OUTPUT:
A string representation of the D-finite sequence ring.
EXAMPLES:
sage: from rec_sequences.DFiniteSequenceRing import * sage: R.<n> = PolynomialRing(QQ) sage: DFiniteSequenceRing(R) Ring of D-finite sequences over Rational Field sage: S.<n> = PolynomialRing(NumberField(x^2-2, "z")) sage: DFiniteSequenceRing(S) Ring of D-finite sequences over Number Field in z with defining polynomial x^2 - 2
- associated_ore_algebra()
Returns the ore algebra associated to the recurrences.
OUTPUT:
The ore algebra containing recurrences of this D-finite sequence ring.
EXAMPLES:
sage: from rec_sequences.DFiniteSequenceRing import * sage: R.<n> = PolynomialRing(QQ) sage: D = DFiniteSequenceRing(R) sage: D.associated_ore_algebra() Univariate Ore algebra in Sn over Univariate Polynomial Ring in n over Rational Field
- base()
The polynomial base ring.
OUTPUT:
The polynomial base ring.
EXAMPLES:
sage: from rec_sequences.DFiniteSequenceRing import * sage: R.<n> = PolynomialRing(QQ) sage: D = DFiniteSequenceRing(R) sage: D.base() Univariate Polynomial Ring in n over Rational Field
- change_base_ring(R)
Return a copy of
selfbut with the base ringR.OUTPUT:
A D-finite sequence ring with base
R.
- construction()
Shows how the given ring can be constructed using functors.
OUTPUT:
A functor
Fand a ringRsuch thatF(R)==selfEXAMPLES:
sage: from rec_sequences.DFiniteSequenceRing import * sage: R.<n> = PolynomialRing(QQ) sage: D = DFiniteSequenceRing(R) sage: F, R = D.construction() sage: F._apply_functor(R) == D True
- get_ore_algebra_ring()
Return the associated
DFiniteFunctionRing.OUTPUT:
Associated
DFiniteFunctionRing.EXAMPLES:
sage: from rec_sequences.DFiniteSequenceRing import * sage: R.<n> = PolynomialRing(QQ) sage: D = DFiniteSequenceRing(R) sage: D.get_ore_algebra_ring() Ring of D-finite sequences over Univariate Polynomial Ring in n over Rational Field
- guess(data, name='a', *args, **kwds)
From given values guess a D-finite sequence using the
ore_algebrapackage.INPUT:
data– list of elements in the base field of the C-finite sequence ringname(default: “a”) – a name for the resulting sequencealgorithm(optional) – if “rec_sequences”, then another straightforward implementation for sequences over QQ is used.operator_only(optional) – ifTrueonly the ore-operator of the recurrence is returned.
All additional arguments are passed to the
ore_algebraguessing routine.OUTPUT:
A D-finite sequence with the specified terms as initial values. If no such sequence is found, a
ValueErroris raised.EXAMPLES:
sage: from rec_sequences.DFiniteSequenceRing import * sage: R.<n> = PolynomialRing(QQ) sage: D = DFiniteSequenceRing(R) sage: D.guess([factorial(n) for n in range(20)]) D-finite sequence a(n): (n + 1)*a(n) + (-1)*a(n+1) = 0 and a(0)=1 sage: D.guess([sum(1/i for i in range(1,n)) for n in range(1,20)], ....: algorithm="rec_sequences", operator_only = True) (1/2*n + 3/2)*Sn^3 + (-n - 5/2)*Sn^2 + (1/2*n + 1)*Sn
- class rec_sequences.DFiniteSequenceRing.DFiniteSequenceRingFunctor
Bases:
rec_sequences.RecurrenceSequenceRing.RecurrenceSequenceRingFunctor- __init__()
Constructs a
DFiniteSequenceRingFunctor.
- _repr_()
Returns a string representation of the functor.
OUTPUT:
The string “DFiniteSequenceRing(*)” .