C-Finite Sequences#

C-finite infinite sequences satisfy homogeneous linear recurrences with constant coefficients:

\[a_{n+d} = c_0a_n + c_1a_{n+1} + \cdots + c_{d-1}a_{n+d-1}, \quad d>0.\]

CFiniteSequences are completely defined by their ordinary generating function (o.g.f., which is always a fraction of polynomials over \(\ZZ\) or \(\QQ\) ).

EXAMPLES:

sage: fibo = CFiniteSequence(x/(1-x-x^2))        # the Fibonacci sequence
sage: fibo
C-finite sequence, generated by -x/(x^2 + x - 1)
sage: fibo.parent()
The ring of C-Finite sequences in x over Rational Field
sage: fibo.parent().category()
Category of commutative rings
sage: C.<x> = CFiniteSequences(QQ)
sage: fibo.parent() == C
True
sage: C
The ring of C-Finite sequences in x over Rational Field
sage: C(x/(1-x-x^2))
C-finite sequence, generated by -x/(x^2 + x - 1)
sage: C(x/(1-x-x^2)) == fibo
True
sage: var('y')
y
sage: CFiniteSequence(y/(1-y-y^2))
C-finite sequence, generated by -y/(y^2 + y - 1)
sage: CFiniteSequence(y/(1-y-y^2)) == fibo
False

Finite subsets of the sequence are accessible via python slices:

sage: fibo[137]                 #the 137th term of the Fibonacci sequence
19134702400093278081449423917
sage: fibo[137] == fibonacci(137)
True
sage: fibo[0:12]
[0, 1, 1, 2, 3, 5, 8, 13, 21, 34, 55, 89]
sage: fibo[14:4:-2]
[377, 144, 55, 21, 8]

They can be created also from the coefficients and start values of a recurrence:

sage: r = C.from_recurrence([1,1],[0,1])
sage: r == fibo
True

Given enough values, the o.g.f. of a C-finite sequence can be guessed:

sage: r = C.guess([0,1,1,2,3,5,8])
sage: r == fibo
True

AUTHORS:

  • Ralf Stephan (2014): initial version

REFERENCES:

class sage.rings.cfinite_sequence.CFiniteSequence(parent, ogf)#

Bases: sage.structure.element.FieldElement

Create a C-finite sequence given its ordinary generating function.

INPUT:

  • ogf – a rational function, the ordinary generating function (can be an element from the symbolic ring, fraction field or polynomial ring)

OUTPUT:

  • A CFiniteSequence object

EXAMPLES:

sage: CFiniteSequence((2-x)/(1-x-x^2))     # the Lucas sequence
C-finite sequence, generated by (x - 2)/(x^2 + x - 1)
sage: CFiniteSequence(x/(1-x)^3)           # triangular numbers
C-finite sequence, generated by -x/(x^3 - 3*x^2 + 3*x - 1)

Polynomials are interpreted as finite sequences, or recurrences of degree 0:

sage: CFiniteSequence(x^2-4*x^5)
Finite sequence [1, 0, 0, -4], offset = 2
sage: CFiniteSequence(1)
Finite sequence [1], offset = 0

This implementation allows any polynomial fraction as o.g.f. by interpreting any power of \(x\) dividing the o.g.f. numerator or denominator as a right or left shift of the sequence offset:

sage: CFiniteSequence(x^2+3/x)
Finite sequence [3, 0, 0, 1], offset = -1
sage: CFiniteSequence(1/x+4/x^3)
Finite sequence [4, 0, 1], offset = -3
sage: P = LaurentPolynomialRing(QQ.fraction_field(), 'X')
sage: X=P.gen()
sage: CFiniteSequence(1/(1-X))
C-finite sequence, generated by -1/(X - 1)

The o.g.f. is always normalized to get a denominator constant coefficient of \(+1\):

sage: CFiniteSequence(1/(x-2))
C-finite sequence, generated by 1/(x - 2)

The given ogf is used to create an appropriate parent: it can be a symbolic expression, a polynomial , or a fraction field element as long as it can be coerced into a proper fraction field over the rationals:

sage: var('x')
x
sage: f1 = CFiniteSequence((2-x)/(1-x-x^2))
sage: P.<x> = QQ[]
sage: f2 = CFiniteSequence((2-x)/(1-x-x^2))
sage: f1 == f2
True
sage: f1.parent()
The ring of C-Finite sequences in x over Rational Field
sage: f1.ogf().parent()
Fraction Field of Univariate Polynomial Ring in x over Rational Field
sage: CFiniteSequence(log(x))
Traceback (most recent call last):
...
TypeError: unable to convert log(x) to a rational
closed_form(n='n')#

Return a symbolic expression in n, which equals the n-th term of the sequence.

It is a well-known property of C-finite sequences a_n that they have a closed form of the type:

\[a_n = \sum_{i=1}^d c_i(n) \cdot r_i^n,\]

where r_i are the roots of the characteristic equation and c_i(n) is a polynomial (whose degree equals the multiplicity of r_i minus one). This is a natural generalization of Binet’s formula for Fibonacci numbers. See, for instance, [KP2011, Theorem 4.1].

Note that if the o.g.f. has a polynomial part, that is, if the numerator degree is not strictly less than the denominator degree, then this closed form holds only when n exceeds the degree of that polynomial part. In that case, the returned expression will differ from the sequence for small n.

EXAMPLES:

sage: CFiniteSequence(1/(1-x)).closed_form()
1
sage: CFiniteSequence(x^2/(1-x)).closed_form()
1
sage: CFiniteSequence(1/(1-x^2)).closed_form()
1/2*(-1)^n + 1/2
sage: CFiniteSequence(1/(1+x^3)).closed_form()
1/3*(-1)^n + 1/3*(1/2*I*sqrt(3) + 1/2)^n + 1/3*(-1/2*I*sqrt(3) + 1/2)^n
sage: CFiniteSequence(1/(1-x)/(1-2*x)/(1-3*x)).closed_form()
9/2*3^n - 4*2^n + 1/2

Binet’s formula for the Fibonacci numbers:

sage: CFiniteSequence(x/(1-x-x^2)).closed_form()
sqrt(1/5)*(1/2*sqrt(5) + 1/2)^n - sqrt(1/5)*(-1/2*sqrt(5) + 1/2)^n
sage: [_.subs(n=k).full_simplify() for k in range(6)]
[0, 1, 1, 2, 3, 5]

sage: CFiniteSequence((4*x+3)/(1-2*x-5*x^2)).closed_form()
1/2*(sqrt(6) + 1)^n*(7*sqrt(1/6) + 3) - 1/2*(-sqrt(6) + 1)^n*(7*sqrt(1/6) - 3)

Examples with multiple roots:

sage: CFiniteSequence(x*(x^2+4*x+1)/(1-x)^5).closed_form()
1/4*n^4 + 1/2*n^3 + 1/4*n^2
sage: CFiniteSequence((1+2*x-x^2)/(1-x)^4/(1+x)^2).closed_form()
1/12*n^3 - 1/8*(-1)^n*(n + 1) + 3/4*n^2 + 43/24*n + 9/8
sage: CFiniteSequence(1/(1-x)^3/(1-2*x)^4).closed_form()
4/3*(n^3 - 3*n^2 + 20*n - 36)*2^n + 1/2*n^2 + 19/2*n + 49
sage: CFiniteSequence((x/(1-x-x^2))^2).closed_form()
1/5*(n - sqrt(1/5))*(1/2*sqrt(5) + 1/2)^n + 1/5*(n + sqrt(1/5))*(-1/2*sqrt(5) + 1/2)^n
coefficients()#

Return the coefficients of the recurrence representation of the C-finite sequence.

OUTPUT:

  • A list of values

EXAMPLES:

sage: C.<x> = CFiniteSequences(QQ)
sage: lucas = C((2-x)/(1-x-x^2))   # the Lucas sequence
sage: lucas.coefficients()
[1, 1]
denominator()#

Return the numerator of the o.g.f of self.

EXAMPLES:

sage: f = CFiniteSequence((2-x)/(1-x-x^2)); f
C-finite sequence, generated by (x - 2)/(x^2 + x - 1)
sage: f.denominator()
x^2 + x - 1
numerator()#

Return the numerator of the o.g.f of self.

EXAMPLES:

sage: f = CFiniteSequence((2-x)/(1-x-x^2)); f
C-finite sequence, generated by (x - 2)/(x^2 + x - 1)
sage: f.numerator()
x - 2
ogf()#

Return the ordinary generating function associated with the CFiniteSequence.

This is always a fraction of polynomials in the base ring.

EXAMPLES:

sage: C.<x> = CFiniteSequences(QQ)
sage: r = C.from_recurrence([2],[1])
sage: r.ogf()
-1/2/(x - 1/2)
sage: C(0).ogf()
0
recurrence_repr()#

Return a string with the recurrence representation of the C-finite sequence.

OUTPUT:

  • A string

EXAMPLES:

sage: C.<x> = CFiniteSequences(QQ)
sage: C((2-x)/(1-x-x^2)).recurrence_repr()
'homogeneous linear recurrence with constant coefficients of degree 2: a(n+2) = a(n+1) + a(n), starting a(0...) = [2, 1]'
sage: C(x/(1-x)^3).recurrence_repr()
'homogeneous linear recurrence with constant coefficients of degree 3: a(n+3) = 3*a(n+2) - 3*a(n+1) + a(n), starting a(1...) = [1, 3, 6]'
sage: C(1).recurrence_repr()
'Finite sequence [1], offset 0'
sage: r = C((-2*x^3 + x^2 - x + 1)/(2*x^2 - 3*x + 1))
sage: r.recurrence_repr()
'homogeneous linear recurrence with constant coefficients of degree 2: a(n+2) = 3*a(n+1) - 2*a(n), starting a(0...) = [1, 2, 5, 9]'
sage: r = CFiniteSequence(x^3/(1-x-x^2))
sage: r.recurrence_repr()
'homogeneous linear recurrence with constant coefficients of degree 2: a(n+2) = a(n+1) + a(n), starting a(3...) = [1, 1, 2, 3]'
series(n)#

Return the Laurent power series associated with the CFiniteSequence, with precision \(n\).

INPUT:

  • \(n\) – a nonnegative integer

EXAMPLES:

sage: C.<x> = CFiniteSequences(QQ)
sage: r = C.from_recurrence([-1,2],[0,1])
sage: s = r.series(4); s
x + 2*x^2 + 3*x^3 + 4*x^4 + O(x^5)
sage: type(s)
<class 'sage.rings.laurent_series_ring_element.LaurentSeries'>
sage.rings.cfinite_sequence.CFiniteSequences(base_ring, names=None, category=None)#

Return the ring of C-Finite sequences.

The ring is defined over a base ring (\(\ZZ\) or \(\QQ\) ) and each element is represented by its ordinary generating function (ogf) which is a rational function over the base ring.

INPUT:

  • base_ring – the base ring to construct the fraction field representing the C-Finite sequences

  • names – (optional) the list of variables.

EXAMPLES:

sage: C.<x> = CFiniteSequences(QQ)
sage: C
The ring of C-Finite sequences in x over Rational Field
sage: C.an_element()
C-finite sequence, generated by (x - 2)/(x^2 + x - 1)
sage: C.category()
Category of commutative rings
sage: C.one()
Finite sequence [1], offset = 0
sage: C.zero()
Constant infinite sequence 0.
sage: C(x)
Finite sequence [1], offset = 1
sage: C(1/x)
Finite sequence [1], offset = -1
sage: C((-x + 2)/(-x^2 - x + 1))
C-finite sequence, generated by (x - 2)/(x^2 + x - 1)
class sage.rings.cfinite_sequence.CFiniteSequences_generic(polynomial_ring, category)#

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

The class representing the ring of C-Finite Sequences

Element#

alias of CFiniteSequence

an_element()#

Return an element of C-Finite Sequences.

OUTPUT:

The Lucas sequence.

EXAMPLES:

sage: C.<x> = CFiniteSequences(QQ)
sage: C.an_element()
C-finite sequence, generated by (x - 2)/(x^2 + x - 1)
fraction_field()#

Return the fraction field used to represent the elements of self.

EXAMPLES:

sage: C.<x> = CFiniteSequences(QQ)
sage: C.fraction_field()
Fraction Field of Univariate Polynomial Ring in x over Rational Field
from_recurrence(coefficients, values)#

Create a C-finite sequence given the coefficients \(c\) and starting values \(a\) of a homogeneous linear recurrence.

\[a_{n+d} = c_0a_n + c_1a_{n+1} + \cdots + c_{d-1}a_{n+d-1}, \quad d\ge0.\]

INPUT:

  • coefficients – a list of rationals

  • values – start values, a list of rationals

OUTPUT:

  • A CFiniteSequence object

EXAMPLES:

sage: C.<x> = CFiniteSequences(QQ)
sage: C.from_recurrence([1,1],[0,1])   # Fibonacci numbers
C-finite sequence, generated by -x/(x^2 + x - 1)
sage: C.from_recurrence([-1,2],[0,1])    # natural numbers
C-finite sequence, generated by x/(x^2 - 2*x + 1)
sage: r = C.from_recurrence([-1],[1])
sage: s = C.from_recurrence([-1],[1,-1])
sage: r == s
True
sage: r = C(x^3/(1-x-x^2))
sage: s = C.from_recurrence([1,1],[0,0,0,1,1])
sage: r == s
True
sage: C.from_recurrence(1,1)
Traceback (most recent call last):
...
ValueError: Wrong type for recurrence coefficient list.
gen(i=0)#

Return the i-th generator of self.

INPUT:

  • i – an integer (default:0)

EXAMPLES:

sage: C.<x> = CFiniteSequences(QQ)
sage: C.gen()
x
sage: x == C.gen()
True
guess(sequence, algorithm='sage')#

Return the minimal CFiniteSequence that generates the sequence.

Assume the first value has index 0.

INPUT:

  • sequence – list of integers

  • algorithm – string
    • ‘sage’ - the default is to use Sage’s matrix kernel function

    • ‘pari’ - use Pari’s implementation of LLL

    • ‘bm’ - use Sage’s Berlekamp-Massey algorithm

OUTPUT:

  • a CFiniteSequence, or 0 if none could be found

With the default kernel method, trailing zeroes are chopped off before a guessing attempt. This may reduce the data below the accepted length of six values.

EXAMPLES:

sage: C.<x> = CFiniteSequences(QQ)
sage: C.guess([1,2,4,8,16,32])
C-finite sequence, generated by -1/2/(x - 1/2)
sage: r = C.guess([1,2,3,4,5])
Traceback (most recent call last):
...
ValueError: Sequence too short for guessing.

With Berlekamp-Massey, if an odd number of values is given, the last one is dropped. So with an odd number of values the result may not generate the last value:

sage: r = C.guess([1,2,4,8,9], algorithm='bm'); r
C-finite sequence, generated by -1/2/(x - 1/2)
sage: r[0:5]
[1, 2, 4, 8, 16]
ngens()#

Return the number of generators of self

EXAMPLES:

sage: from sage.rings.cfinite_sequence import CFiniteSequences
sage: C.<x> = CFiniteSequences(QQ)
sage: C.ngens()
1
polynomial_ring()#

Return the polynomial ring used to represent the elements of self.

EXAMPLES:

sage: C.<x> = CFiniteSequences(QQ)
sage: C.polynomial_ring()
Univariate Polynomial Ring in x over Rational Field