Power series implemented using PARI

EXAMPLES:

This implementation can be selected for any base ring supported by PARI by passing the keyword implementation='pari' to the PowerSeriesRing() constructor:

sage: R.<q> = PowerSeriesRing(ZZ, implementation='pari'); R
Power Series Ring in q over Integer Ring
sage: S.<t> = PowerSeriesRing(CC, implementation='pari'); S                         # needs sage.rings.real_mpfr
Power Series Ring in t over Complex Field with 53 bits of precision
>>> from sage.all import *
>>> R = PowerSeriesRing(ZZ, implementation='pari', names=('q',)); (q,) = R._first_ngens(1); R
Power Series Ring in q over Integer Ring
>>> S = PowerSeriesRing(CC, implementation='pari', names=('t',)); (t,) = S._first_ngens(1); S                         # needs sage.rings.real_mpfr
Power Series Ring in t over Complex Field with 53 bits of precision

Note that only the type of the elements depends on the implementation, not the type of the parents:

sage: type(R)
<class 'sage.rings.power_series_ring.PowerSeriesRing_domain_with_category'>
sage: type(q)
<class 'sage.rings.power_series_pari.PowerSeries_pari'>
sage: type(S)                                                                       # needs sage.rings.real_mpfr
<class 'sage.rings.power_series_ring.PowerSeriesRing_over_field_with_category'>
sage: type(t)                                                                       # needs sage.rings.real_mpfr
<class 'sage.rings.power_series_pari.PowerSeries_pari'>
>>> from sage.all import *
>>> type(R)
<class 'sage.rings.power_series_ring.PowerSeriesRing_domain_with_category'>
>>> type(q)
<class 'sage.rings.power_series_pari.PowerSeries_pari'>
>>> type(S)                                                                       # needs sage.rings.real_mpfr
<class 'sage.rings.power_series_ring.PowerSeriesRing_over_field_with_category'>
>>> type(t)                                                                       # needs sage.rings.real_mpfr
<class 'sage.rings.power_series_pari.PowerSeries_pari'>

If \(k\) is a finite field implemented using PARI, this is the default implementation for power series over \(k\):

sage: k.<c> = GF(5^12)
sage: type(c)
<class 'sage.rings.finite_rings.element_pari_ffelt.FiniteFieldElement_pari_ffelt'>
sage: A.<x> = k[[]]
sage: type(x)
<class 'sage.rings.power_series_pari.PowerSeries_pari'>
>>> from sage.all import *
>>> k = GF(Integer(5)**Integer(12), names=('c',)); (c,) = k._first_ngens(1)
>>> type(c)
<class 'sage.rings.finite_rings.element_pari_ffelt.FiniteFieldElement_pari_ffelt'>
>>> A = k[['x']]; (x,) = A._first_ngens(1)
>>> type(x)
<class 'sage.rings.power_series_pari.PowerSeries_pari'>

Warning

Because this implementation uses the PARI interface, the PARI variable ordering must be respected in the sense that the variable name of the power series ring must have higher priority than any variable names occurring in the base ring:

sage: R.<y> = QQ[]
sage: S.<x> = PowerSeriesRing(R, implementation='pari'); S
Power Series Ring in x over Univariate Polynomial Ring in y over Rational Field
>>> from sage.all import *
>>> R = QQ['y']; (y,) = R._first_ngens(1)
>>> S = PowerSeriesRing(R, implementation='pari', names=('x',)); (x,) = S._first_ngens(1); S
Power Series Ring in x over Univariate Polynomial Ring in y over Rational Field

Reversing the variable ordering leads to errors:

sage: R.<x> = QQ[]
sage: S.<y> = PowerSeriesRing(R, implementation='pari')
Traceback (most recent call last):
...
PariError: incorrect priority in gtopoly: variable x <= y
>>> from sage.all import *
>>> R = QQ['x']; (x,) = R._first_ngens(1)
>>> S = PowerSeriesRing(R, implementation='pari', names=('y',)); (y,) = S._first_ngens(1)
Traceback (most recent call last):
...
PariError: incorrect priority in gtopoly: variable x <= y

AUTHORS:

  • Peter Bruin (December 2013): initial version

class sage.rings.power_series_pari.PowerSeries_pari[source]

Bases: PowerSeries

A power series implemented using PARI.

INPUT:

  • parent – the power series ring to use as the parent

  • f – object from which to construct a power series

  • prec – (default: infinity) precision of the element to be constructed

  • check – ignored, but accepted for compatibility with PowerSeries_poly

dict()[source]

Return a dictionary of coefficients for self.

This is simply a dict for the underlying polynomial; it need not have keys corresponding to every number smaller than self.prec().

EXAMPLES:

sage: R.<t> = PowerSeriesRing(ZZ, implementation='pari')
sage: f = 1 + t^10 + O(t^12)
sage: f.monomial_coefficients()
{0: 1, 10: 1}
>>> from sage.all import *
>>> R = PowerSeriesRing(ZZ, implementation='pari', names=('t',)); (t,) = R._first_ngens(1)
>>> f = Integer(1) + t**Integer(10) + O(t**Integer(12))
>>> f.monomial_coefficients()
{0: 1, 10: 1}

dict is an alias:

sage: f.dict()
{0: 1, 10: 1}
>>> from sage.all import *
>>> f.dict()
{0: 1, 10: 1}
integral(var=None)[source]

Return the formal integral of self.

By default, the integration variable is the variable of the power series. Otherwise, the integration variable is the optional parameter var.

Note

The integral is always chosen so the constant term is 0.

EXAMPLES:

sage: k.<w> = PowerSeriesRing(QQ, implementation='pari')
sage: (1+17*w+15*w^3+O(w^5)).integral()
w + 17/2*w^2 + 15/4*w^4 + O(w^6)
sage: (w^3 + 4*w^4 + O(w^7)).integral()
1/4*w^4 + 4/5*w^5 + O(w^8)
sage: (3*w^2).integral()
w^3
>>> from sage.all import *
>>> k = PowerSeriesRing(QQ, implementation='pari', names=('w',)); (w,) = k._first_ngens(1)
>>> (Integer(1)+Integer(17)*w+Integer(15)*w**Integer(3)+O(w**Integer(5))).integral()
w + 17/2*w^2 + 15/4*w^4 + O(w^6)
>>> (w**Integer(3) + Integer(4)*w**Integer(4) + O(w**Integer(7))).integral()
1/4*w^4 + 4/5*w^5 + O(w^8)
>>> (Integer(3)*w**Integer(2)).integral()
w^3
list()[source]

Return the list of known coefficients for self.

This is just the list of coefficients of the underlying polynomial; it need not have length equal to self.prec().

EXAMPLES:

sage: R.<t> = PowerSeriesRing(ZZ, implementation='pari')
sage: f = 1 - 5*t^3 + t^5 + O(t^7)
sage: f.list()
[1, 0, 0, -5, 0, 1]

sage: # needs sage.rings.padics
sage: S.<u> = PowerSeriesRing(pAdicRing(5), implementation='pari')
sage: (2 + u).list()
[2 + O(5^20), 1 + O(5^20)]
>>> from sage.all import *
>>> R = PowerSeriesRing(ZZ, implementation='pari', names=('t',)); (t,) = R._first_ngens(1)
>>> f = Integer(1) - Integer(5)*t**Integer(3) + t**Integer(5) + O(t**Integer(7))
>>> f.list()
[1, 0, 0, -5, 0, 1]

>>> # needs sage.rings.padics
>>> S = PowerSeriesRing(pAdicRing(Integer(5)), implementation='pari', names=('u',)); (u,) = S._first_ngens(1)
>>> (Integer(2) + u).list()
[2 + O(5^20), 1 + O(5^20)]
monomial_coefficients()[source]

Return a dictionary of coefficients for self.

This is simply a dict for the underlying polynomial; it need not have keys corresponding to every number smaller than self.prec().

EXAMPLES:

sage: R.<t> = PowerSeriesRing(ZZ, implementation='pari')
sage: f = 1 + t^10 + O(t^12)
sage: f.monomial_coefficients()
{0: 1, 10: 1}
>>> from sage.all import *
>>> R = PowerSeriesRing(ZZ, implementation='pari', names=('t',)); (t,) = R._first_ngens(1)
>>> f = Integer(1) + t**Integer(10) + O(t**Integer(12))
>>> f.monomial_coefficients()
{0: 1, 10: 1}

dict is an alias:

sage: f.dict()
{0: 1, 10: 1}
>>> from sage.all import *
>>> f.dict()
{0: 1, 10: 1}
padded_list(n=None)[source]

Return a list of coefficients of self up to (but not including) \(q^n\).

The list is padded with zeroes on the right so that it has length \(n\).

INPUT:

  • n – nonnegative integer (optional); if \(n\) is not given, it will be taken to be the precision of self, unless this is +Infinity, in which case we just return self.list()

EXAMPLES:

sage: R.<q> = PowerSeriesRing(QQ, implementation='pari')
sage: f = 1 - 17*q + 13*q^2 + 10*q^4 + O(q^7)
sage: f.list()
[1, -17, 13, 0, 10]
sage: f.padded_list(7)
[1, -17, 13, 0, 10, 0, 0]
sage: f.padded_list(10)
[1, -17, 13, 0, 10, 0, 0, 0, 0, 0]
sage: f.padded_list(3)
[1, -17, 13]
sage: f.padded_list()
[1, -17, 13, 0, 10, 0, 0]
sage: g = 1 - 17*q + 13*q^2 + 10*q^4
sage: g.list()
[1, -17, 13, 0, 10]
sage: g.padded_list()
[1, -17, 13, 0, 10]
sage: g.padded_list(10)
[1, -17, 13, 0, 10, 0, 0, 0, 0, 0]
>>> from sage.all import *
>>> R = PowerSeriesRing(QQ, implementation='pari', names=('q',)); (q,) = R._first_ngens(1)
>>> f = Integer(1) - Integer(17)*q + Integer(13)*q**Integer(2) + Integer(10)*q**Integer(4) + O(q**Integer(7))
>>> f.list()
[1, -17, 13, 0, 10]
>>> f.padded_list(Integer(7))
[1, -17, 13, 0, 10, 0, 0]
>>> f.padded_list(Integer(10))
[1, -17, 13, 0, 10, 0, 0, 0, 0, 0]
>>> f.padded_list(Integer(3))
[1, -17, 13]
>>> f.padded_list()
[1, -17, 13, 0, 10, 0, 0]
>>> g = Integer(1) - Integer(17)*q + Integer(13)*q**Integer(2) + Integer(10)*q**Integer(4)
>>> g.list()
[1, -17, 13, 0, 10]
>>> g.padded_list()
[1, -17, 13, 0, 10]
>>> g.padded_list(Integer(10))
[1, -17, 13, 0, 10, 0, 0, 0, 0, 0]
polynomial()[source]

Convert self to a polynomial.

EXAMPLES:

sage: R.<t> = PowerSeriesRing(GF(7), implementation='pari')
sage: f = 3 - t^3 + O(t^5)
sage: f.polynomial()
6*t^3 + 3
>>> from sage.all import *
>>> R = PowerSeriesRing(GF(Integer(7)), implementation='pari', names=('t',)); (t,) = R._first_ngens(1)
>>> f = Integer(3) - t**Integer(3) + O(t**Integer(5))
>>> f.polynomial()
6*t^3 + 3
reverse(precision=None)[source]

Return the reverse of self.

The reverse of a power series \(f\) is the power series \(g\) such that \(g(f(x)) = x\). This exists if and only if the valuation of self is exactly 1 and the coefficient of \(x\) is a unit.

If the optional argument precision is given, the reverse is returned with this precision. If f has infinite precision and the argument precision is not given, then the reverse is returned with the default precision of f.parent().

EXAMPLES:

sage: R.<x> = PowerSeriesRing(QQ, implementation='pari')
sage: f = 2*x + 3*x^2 - x^4 + O(x^5)
sage: g = f.reverse()
sage: g
1/2*x - 3/8*x^2 + 9/16*x^3 - 131/128*x^4 + O(x^5)
sage: f(g)
x + O(x^5)
sage: g(f)
x + O(x^5)

sage: A.<t> = PowerSeriesRing(ZZ, implementation='pari')
sage: a = t - t^2 - 2*t^4 + t^5 + O(t^6)
sage: b = a.reverse(); b
t + t^2 + 2*t^3 + 7*t^4 + 25*t^5 + O(t^6)
sage: a(b)
t + O(t^6)
sage: b(a)
t + O(t^6)

sage: B.<b,c> = PolynomialRing(ZZ)
sage: A.<t> = PowerSeriesRing(B, implementation='pari')
sage: f = t + b*t^2 + c*t^3 + O(t^4)
sage: g = f.reverse(); g
t - b*t^2 + (2*b^2 - c)*t^3 + O(t^4)
sage: f(g)
t + O(t^4)
sage: g(f)
t + O(t^4)

sage: A.<t> = PowerSeriesRing(ZZ, implementation='pari')
sage: B.<x> = PowerSeriesRing(A, implementation='pari')
sage: f = (1 - 3*t + 4*t^3 + O(t^4))*x + (2 + t + t^2 + O(t^3))*x^2 + O(x^3)
sage: g = f.reverse(); g
(1 + 3*t + 9*t^2 + 23*t^3 + O(t^4))*x + (-2 - 19*t - 118*t^2 + O(t^3))*x^2 + O(x^3)
>>> from sage.all import *
>>> R = PowerSeriesRing(QQ, implementation='pari', names=('x',)); (x,) = R._first_ngens(1)
>>> f = Integer(2)*x + Integer(3)*x**Integer(2) - x**Integer(4) + O(x**Integer(5))
>>> g = f.reverse()
>>> g
1/2*x - 3/8*x^2 + 9/16*x^3 - 131/128*x^4 + O(x^5)
>>> f(g)
x + O(x^5)
>>> g(f)
x + O(x^5)

>>> A = PowerSeriesRing(ZZ, implementation='pari', names=('t',)); (t,) = A._first_ngens(1)
>>> a = t - t**Integer(2) - Integer(2)*t**Integer(4) + t**Integer(5) + O(t**Integer(6))
>>> b = a.reverse(); b
t + t^2 + 2*t^3 + 7*t^4 + 25*t^5 + O(t^6)
>>> a(b)
t + O(t^6)
>>> b(a)
t + O(t^6)

>>> B = PolynomialRing(ZZ, names=('b', 'c',)); (b, c,) = B._first_ngens(2)
>>> A = PowerSeriesRing(B, implementation='pari', names=('t',)); (t,) = A._first_ngens(1)
>>> f = t + b*t**Integer(2) + c*t**Integer(3) + O(t**Integer(4))
>>> g = f.reverse(); g
t - b*t^2 + (2*b^2 - c)*t^3 + O(t^4)
>>> f(g)
t + O(t^4)
>>> g(f)
t + O(t^4)

>>> A = PowerSeriesRing(ZZ, implementation='pari', names=('t',)); (t,) = A._first_ngens(1)
>>> B = PowerSeriesRing(A, implementation='pari', names=('x',)); (x,) = B._first_ngens(1)
>>> f = (Integer(1) - Integer(3)*t + Integer(4)*t**Integer(3) + O(t**Integer(4)))*x + (Integer(2) + t + t**Integer(2) + O(t**Integer(3)))*x**Integer(2) + O(x**Integer(3))
>>> g = f.reverse(); g
(1 + 3*t + 9*t^2 + 23*t^3 + O(t^4))*x + (-2 - 19*t - 118*t^2 + O(t^3))*x^2 + O(x^3)

The optional argument precision sets the precision of the output:

sage: R.<x> = PowerSeriesRing(QQ, implementation='pari')
sage: f = 2*x + 3*x^2 - 7*x^3 + x^4 + O(x^5)
sage: g = f.reverse(precision=3); g
1/2*x - 3/8*x^2 + O(x^3)
sage: f(g)
x + O(x^3)
sage: g(f)
x + O(x^3)
>>> from sage.all import *
>>> R = PowerSeriesRing(QQ, implementation='pari', names=('x',)); (x,) = R._first_ngens(1)
>>> f = Integer(2)*x + Integer(3)*x**Integer(2) - Integer(7)*x**Integer(3) + x**Integer(4) + O(x**Integer(5))
>>> g = f.reverse(precision=Integer(3)); g
1/2*x - 3/8*x^2 + O(x^3)
>>> f(g)
x + O(x^3)
>>> g(f)
x + O(x^3)

If the input series has infinite precision, the precision of the output is automatically set to the default precision of the parent ring:

sage: R.<x> = PowerSeriesRing(QQ, default_prec=20, implementation='pari')
sage: (x - x^2).reverse()  # get some Catalan numbers
x + x^2 + 2*x^3 + 5*x^4 + 14*x^5 + 42*x^6 + 132*x^7 + 429*x^8
 + 1430*x^9 + 4862*x^10 + 16796*x^11 + 58786*x^12 + 208012*x^13
 + 742900*x^14 + 2674440*x^15 + 9694845*x^16 + 35357670*x^17
 + 129644790*x^18 + 477638700*x^19 + O(x^20)
sage: (x - x^2).reverse(precision=3)
x + x^2 + O(x^3)
>>> from sage.all import *
>>> R = PowerSeriesRing(QQ, default_prec=Integer(20), implementation='pari', names=('x',)); (x,) = R._first_ngens(1)
>>> (x - x**Integer(2)).reverse()  # get some Catalan numbers
x + x^2 + 2*x^3 + 5*x^4 + 14*x^5 + 42*x^6 + 132*x^7 + 429*x^8
 + 1430*x^9 + 4862*x^10 + 16796*x^11 + 58786*x^12 + 208012*x^13
 + 742900*x^14 + 2674440*x^15 + 9694845*x^16 + 35357670*x^17
 + 129644790*x^18 + 477638700*x^19 + O(x^20)
>>> (x - x**Integer(2)).reverse(precision=Integer(3))
x + x^2 + O(x^3)
valuation()[source]

Return the valuation of self.

EXAMPLES:

sage: R.<t> = PowerSeriesRing(QQ, implementation='pari')
sage: (5 - t^8 + O(t^11)).valuation()
0
sage: (-t^8 + O(t^11)).valuation()
8
sage: O(t^7).valuation()
7
sage: R(0).valuation()
+Infinity
>>> from sage.all import *
>>> R = PowerSeriesRing(QQ, implementation='pari', names=('t',)); (t,) = R._first_ngens(1)
>>> (Integer(5) - t**Integer(8) + O(t**Integer(11))).valuation()
0
>>> (-t**Integer(8) + O(t**Integer(11))).valuation()
8
>>> O(t**Integer(7)).valuation()
7
>>> R(Integer(0)).valuation()
+Infinity