Base class for generic \(p\)-adic polynomials#

This provides common functionality for all \(p\)-adic polynomials, such as printing and factoring.

AUTHORS:

  • Jeroen Demeyer (2013-11-22): initial version, split off from other files, made Polynomial_padic the common base class for all p-adic polynomials.

class sage.rings.polynomial.padics.polynomial_padic.Polynomial_padic(parent, x=None, check=True, is_gen=False, construct=False)#

Bases: Polynomial

content()#

Compute the content of this polynomial.

OUTPUT:

If this is the zero polynomial, return the constant coefficient. Otherwise, since the content is only defined up to a unit, return the content as \(\pi^k\) with maximal precision where \(k\) is the minimal valuation of any of the coefficients.

EXAMPLES:

sage: # needs sage.libs.ntl
sage: K = Zp(13,7)
sage: R.<t> = K[]
sage: f = 13^7*t^3 + K(169,4)*t - 13^4
sage: f.content()
13^2 + O(13^9)
sage: R(0).content()
0
sage: f = R(K(0,3)); f
O(13^3)
sage: f.content()
O(13^3)

sage: # needs sage.libs.ntl
sage: P.<x> = ZZ[]
sage: f = x + 2
sage: f.content()
1
sage: fp = f.change_ring(pAdicRing(2, 10))
sage: fp
(1 + O(2^10))*x + 2 + O(2^11)
sage: fp.content()
1 + O(2^10)
sage: (2*fp).content()
2 + O(2^11)

Over a field it would be sufficient to return only zero or one, as the content is only defined up to multiplication with a unit. However, we return \(\pi^k\) where \(k\) is the minimal valuation of any coefficient:

sage: # needs sage.libs.ntl
sage: K = Qp(13,7)
sage: R.<t> = K[]
sage: f = 13^7*t^3 + K(169,4)*t - 13^-4
sage: f.content()
13^-4 + O(13^3)
sage: f = R.zero()
sage: f.content()
0
sage: f = R(K(0,3))
sage: f.content()
O(13^3)
sage: f = 13*t^3 + K(0,1)*t
sage: f.content()
13 + O(13^8)
factor()#

Return the factorization of this polynomial.

EXAMPLES:

sage: # needs sage.libs.ntl
sage: R.<t> = PolynomialRing(Qp(3, 3, print_mode='terse', print_pos=False))
sage: pol = t^8 - 1
sage: for p,e in pol.factor():
....:     print("{} {}".format(e, p))
1 (1 + O(3^3))*t + 1 + O(3^3)
1 (1 + O(3^3))*t - 1 + O(3^3)
1 (1 + O(3^3))*t^2 + (5 + O(3^3))*t - 1 + O(3^3)
1 (1 + O(3^3))*t^2 + (-5 + O(3^3))*t - 1 + O(3^3)
1 (1 + O(3^3))*t^2 + O(3^3)*t + 1 + O(3^3)
sage: R.<t> = PolynomialRing(Qp(5, 6, print_mode='terse', print_pos=False))
sage: pol = 100 * (5*t - 1) * (t - 5); pol
(500 + O(5^9))*t^2 + (-2600 + O(5^8))*t + 500 + O(5^9)
sage: pol.factor()
(500 + O(5^9)) * ((1 + O(5^5))*t - 1/5 + O(5^5)) * ((1 + O(5^6))*t - 5 + O(5^6))
sage: pol.factor().value()
(500 + O(5^8))*t^2 + (-2600 + O(5^8))*t + 500 + O(5^8)

The same factorization over \(\ZZ_p\). In this case, the “unit” part is a \(p\)-adic unit and the power of \(p\) is considered to be a factor:

sage: # needs sage.libs.ntl
sage: R.<t> = PolynomialRing(Zp(5, 6, print_mode='terse', print_pos=False))
sage: pol = 100 * (5*t - 1) * (t - 5); pol
(500 + O(5^9))*t^2 + (-2600 + O(5^8))*t + 500 + O(5^9)
sage: pol.factor()
(4 + O(5^6)) * (5 + O(5^7))^2 * ((1 + O(5^6))*t - 5 + O(5^6)) * ((5 + O(5^6))*t - 1 + O(5^6))
sage: pol.factor().value()
(500 + O(5^8))*t^2 + (-2600 + O(5^8))*t + 500 + O(5^8)

In the following example, the discriminant is zero, so the \(p\)-adic factorization is not well defined:

sage: factor(t^2)                                                           # needs sage.libs.ntl
Traceback (most recent call last):
...
PrecisionError: p-adic factorization not well-defined since
the discriminant is zero up to the requestion p-adic precision

An example of factoring a constant polynomial (see github issue #26669):

sage: R.<x> = Qp(5)[]                                                       # needs sage.libs.ntl
sage: R(2).factor()                                                         # needs sage.libs.ntl
2 + O(5^20)

More examples over \(\ZZ_p\):

sage: R.<w> = PolynomialRing(Zp(5, prec=6, type='capped-abs', print_mode='val-unit'))
sage: f = w^5 - 1
sage: f.factor()                                                            # needs sage.libs.ntl
((1 + O(5^6))*w + 3124 + O(5^6))
* ((1 + O(5^6))*w^4 + (12501 + O(5^6))*w^3 + (9376 + O(5^6))*w^2
     + (6251 + O(5^6))*w + 3126 + O(5^6))

See github issue #4038:

sage: # needs sage.libs.ntl sage.schemes
sage: E = EllipticCurve('37a1')
sage: K = Qp(7,10)
sage: EK = E.base_extend(K)
sage: g = EK.division_polynomial_0(3)
sage: g.factor()
(3 + O(7^10))
* ((1 + O(7^10))*x
     + 1 + 2*7 + 4*7^2 + 2*7^3 + 5*7^4 + 7^5 + 5*7^6 + 3*7^7 + 5*7^8 + 3*7^9 + O(7^10))
* ((1 + O(7^10))*x^3
     + (6 + 4*7 + 2*7^2 + 4*7^3 + 7^4 + 5*7^5
          + 7^6 + 3*7^7 + 7^8 + 3*7^9 + O(7^10))*x^2
     + (6 + 3*7 + 5*7^2 + 2*7^4 + 7^5 + 7^6 + 2*7^8 + 3*7^9 + O(7^10))*x
     + 2 + 5*7 + 4*7^2 + 2*7^3 + 6*7^4 + 3*7^5 + 7^6 + 4*7^7 + O(7^10))
root_field(names, check_irreducible=True, **kwds)#

Return the \(p\)-adic extension field generated by the roots of the irreducible polynomial self.

INPUT:

  • names – name of the generator of the extension

  • check_irreducible – check whether the polynomial is irreducible

  • kwds – see sage.ring.padics.padic_generic.pAdicGeneric.extension()

EXAMPLES:

sage: R.<x> = Qp(3,5,print_mode='digits')[]                                 # needs sage.libs.ntl
sage: f = x^2 - 3                                                           # needs sage.libs.ntl
sage: f.root_field('x')                                                     # needs sage.libs.ntl
3-adic Eisenstein Extension Field in x defined by x^2 - 3
sage: R.<x> = Qp(5,5,print_mode='digits')[]                                 # needs sage.libs.ntl
sage: f = x^2 - 3                                                           # needs sage.libs.ntl
sage: f.root_field('x', print_mode='bars')                                  # needs sage.libs.ntl
5-adic Unramified Extension Field in x defined by x^2 - 3
sage: R.<x> = Qp(11,5,print_mode='digits')[]                                # needs sage.libs.ntl
sage: f = x^2 - 3                                                           # needs sage.libs.ntl
sage: f.root_field('x', print_mode='bars')                                  # needs sage.libs.ntl
Traceback (most recent call last):
...
ValueError: polynomial must be irreducible