Local Generic Element#

This file contains a common superclass for \(p\)-adic elements and power series elements.

AUTHORS:

  • David Roe: initial version

  • Julian Rueth (2012-10-15, 2014-06-25, 2017-08-04): added inverse_of_unit(); improved add_bigoh(); added _test_expansion()

class sage.rings.padics.local_generic_element.LocalGenericElement[source]#

Bases: CommutativeRingElement

add_bigoh(absprec)[source]#

Return a copy of this element with absolute precision decreased to absprec.

INPUT:

  • absprec – an integer or positive infinity

EXAMPLES:

sage: K = QpCR(3,4)
sage: o = K(1); o
1 + O(3^4)
sage: o.add_bigoh(2)
1 + O(3^2)
sage: o.add_bigoh(-5)
O(3^-5)
>>> from sage.all import *
>>> K = QpCR(Integer(3),Integer(4))
>>> o = K(Integer(1)); o
1 + O(3^4)
>>> o.add_bigoh(Integer(2))
1 + O(3^2)
>>> o.add_bigoh(-Integer(5))
O(3^-5)

One cannot use add_bigoh to lift to a higher precision; this can be accomplished with lift_to_precision():

sage: o.add_bigoh(5)
1 + O(3^4)
>>> from sage.all import *
>>> o.add_bigoh(Integer(5))
1 + O(3^4)

Negative values of absprec return an element in the fraction field of the element’s parent:

sage: R = ZpCA(3,4)
sage: R(3).add_bigoh(-5)
O(3^-5)
>>> from sage.all import *
>>> R = ZpCA(Integer(3),Integer(4))
>>> R(Integer(3)).add_bigoh(-Integer(5))
O(3^-5)

For fixed-mod elements this method truncates the element:

sage: R = ZpFM(3,4)
sage: R(3).add_bigoh(1)
0
>>> from sage.all import *
>>> R = ZpFM(Integer(3),Integer(4))
>>> R(Integer(3)).add_bigoh(Integer(1))
0

If absprec exceeds the precision of the element, then this method has no effect:

sage: R(3).add_bigoh(5)
3
>>> from sage.all import *
>>> R(Integer(3)).add_bigoh(Integer(5))
3

A negative value for absprec returns an element in the fraction field:

sage: R(3).add_bigoh(-1).parent()
3-adic Field with floating precision 4
>>> from sage.all import *
>>> R(Integer(3)).add_bigoh(-Integer(1)).parent()
3-adic Field with floating precision 4
euclidean_degree()[source]#

Return the degree of this element as an element of an Euclidean domain.

EXAMPLES:

For a field, this is always zero except for the zero element:

sage: K = Qp(2)
sage: K.one().euclidean_degree()
0
sage: K.gen().euclidean_degree()
0
sage: K.zero().euclidean_degree()
Traceback (most recent call last):
...
ValueError: euclidean degree not defined for the zero element
>>> from sage.all import *
>>> K = Qp(Integer(2))
>>> K.one().euclidean_degree()
0
>>> K.gen().euclidean_degree()
0
>>> K.zero().euclidean_degree()
Traceback (most recent call last):
...
ValueError: euclidean degree not defined for the zero element

For a ring which is not a field, this is the valuation of the element:

sage: R = Zp(2)
sage: R.one().euclidean_degree()
0
sage: R.gen().euclidean_degree()
1
sage: R.zero().euclidean_degree()
Traceback (most recent call last):
...
ValueError: euclidean degree not defined for the zero element
>>> from sage.all import *
>>> R = Zp(Integer(2))
>>> R.one().euclidean_degree()
0
>>> R.gen().euclidean_degree()
1
>>> R.zero().euclidean_degree()
Traceback (most recent call last):
...
ValueError: euclidean degree not defined for the zero element
inverse_of_unit()[source]#

Returns the inverse of self if self is a unit.

OUTPUT:

  • an element in the same ring as self

EXAMPLES:

sage: R = ZpCA(3,5)
sage: a = R(2); a
2 + O(3^5)
sage: b = a.inverse_of_unit(); b
2 + 3 + 3^2 + 3^3 + 3^4 + O(3^5)
>>> from sage.all import *
>>> R = ZpCA(Integer(3),Integer(5))
>>> a = R(Integer(2)); a
2 + O(3^5)
>>> b = a.inverse_of_unit(); b
2 + 3 + 3^2 + 3^3 + 3^4 + O(3^5)

A ZeroDivisionError is raised if an element has no inverse in the ring:

sage: R(3).inverse_of_unit()
Traceback (most recent call last):
...
ZeroDivisionError: inverse of 3 + O(3^5) does not exist
>>> from sage.all import *
>>> R(Integer(3)).inverse_of_unit()
Traceback (most recent call last):
...
ZeroDivisionError: inverse of 3 + O(3^5) does not exist

Unlike the usual inverse of an element, the result is in the same ring as self and not just in its fraction field:

sage: c = ~a; c
2 + 3 + 3^2 + 3^3 + 3^4 + O(3^5)
sage: a.parent()
3-adic Ring with capped absolute precision 5
sage: b.parent()
3-adic Ring with capped absolute precision 5
sage: c.parent()
3-adic Field with capped relative precision 5
>>> from sage.all import *
>>> c = ~a; c
2 + 3 + 3^2 + 3^3 + 3^4 + O(3^5)
>>> a.parent()
3-adic Ring with capped absolute precision 5
>>> b.parent()
3-adic Ring with capped absolute precision 5
>>> c.parent()
3-adic Field with capped relative precision 5

For fields this does of course not make any difference:

sage: R = QpCR(3,5)
sage: a = R(2)
sage: b = a.inverse_of_unit()
sage: c = ~a
sage: a.parent()
3-adic Field with capped relative precision 5
sage: b.parent()
3-adic Field with capped relative precision 5
sage: c.parent()
3-adic Field with capped relative precision 5
>>> from sage.all import *
>>> R = QpCR(Integer(3),Integer(5))
>>> a = R(Integer(2))
>>> b = a.inverse_of_unit()
>>> c = ~a
>>> a.parent()
3-adic Field with capped relative precision 5
>>> b.parent()
3-adic Field with capped relative precision 5
>>> c.parent()
3-adic Field with capped relative precision 5
is_integral()[source]#

Returns whether self is an integral element.

INPUT:

  • self – a local ring element

OUTPUT:

  • boolean – whether self is an integral element.

EXAMPLES:

sage: R = Qp(3,20)
sage: a = R(7/3); a.is_integral()
False
sage: b = R(7/5); b.is_integral()
True
>>> from sage.all import *
>>> R = Qp(Integer(3),Integer(20))
>>> a = R(Integer(7)/Integer(3)); a.is_integral()
False
>>> b = R(Integer(7)/Integer(5)); b.is_integral()
True
is_padic_unit()[source]#

Returns whether self is a \(p\)-adic unit. That is, whether it has zero valuation.

INPUT:

  • self – a local ring element

OUTPUT:

  • boolean – whether self is a unit

EXAMPLES:

sage: R = Zp(3,20,'capped-rel'); K = Qp(3,20,'capped-rel')
sage: R(0).is_padic_unit()
False
sage: R(1).is_padic_unit()
True
sage: R(2).is_padic_unit()
True
sage: R(3).is_padic_unit()
False
sage: Qp(5,5)(5).is_padic_unit()
False
>>> from sage.all import *
>>> R = Zp(Integer(3),Integer(20),'capped-rel'); K = Qp(Integer(3),Integer(20),'capped-rel')
>>> R(Integer(0)).is_padic_unit()
False
>>> R(Integer(1)).is_padic_unit()
True
>>> R(Integer(2)).is_padic_unit()
True
>>> R(Integer(3)).is_padic_unit()
False
>>> Qp(Integer(5),Integer(5))(Integer(5)).is_padic_unit()
False
is_unit()[source]#

Returns whether self is a unit

INPUT:

  • self – a local ring element

OUTPUT:

  • boolean – whether self is a unit

Note

For fields all nonzero elements are units. For DVR’s, only those elements of valuation 0 are. An older implementation ignored the case of fields, and returned always the negation of self.valuation()==0. This behavior is now supported with self.is_padic_unit().

EXAMPLES:

sage: R = Zp(3,20,'capped-rel'); K = Qp(3,20,'capped-rel')
sage: R(0).is_unit()
False
sage: R(1).is_unit()
True
sage: R(2).is_unit()
True
sage: R(3).is_unit()
False
sage: Qp(5,5)(5).is_unit()  # Note that 5 is invertible in `QQ_5`, even if it has positive valuation!
True
sage: Qp(5,5)(5).is_padic_unit()
False
>>> from sage.all import *
>>> R = Zp(Integer(3),Integer(20),'capped-rel'); K = Qp(Integer(3),Integer(20),'capped-rel')
>>> R(Integer(0)).is_unit()
False
>>> R(Integer(1)).is_unit()
True
>>> R(Integer(2)).is_unit()
True
>>> R(Integer(3)).is_unit()
False
>>> Qp(Integer(5),Integer(5))(Integer(5)).is_unit()  # Note that 5 is invertible in `QQ_5`, even if it has positive valuation!
True
>>> Qp(Integer(5),Integer(5))(Integer(5)).is_padic_unit()
False
normalized_valuation()[source]#

Returns the normalized valuation of this local ring element, i.e., the valuation divided by the absolute ramification index.

INPUT:

  • self – a local ring element.

OUTPUT:

rational – the normalized valuation of self.

EXAMPLES:

sage: Q7 = Qp(7)
sage: R.<x> = Q7[]                                                          # needs sage.libs.ntl
sage: F.<z> = Q7.ext(x^3+7*x+7)                                             # needs sage.libs.ntl
sage: z.normalized_valuation()                                              # needs sage.libs.ntl
1/3
>>> from sage.all import *
>>> Q7 = Qp(Integer(7))
>>> R = Q7['x']; (x,) = R._first_ngens(1)# needs sage.libs.ntl
>>> F = Q7.ext(x**Integer(3)+Integer(7)*x+Integer(7), names=('z',)); (z,) = F._first_ngens(1)# needs sage.libs.ntl
>>> z.normalized_valuation()                                              # needs sage.libs.ntl
1/3
quo_rem(other, integral=False)[source]#

Return the quotient with remainder of the division of this element by other.

INPUT:

  • other – an element in the same ring

  • integral – if True, use integral-style remainders even when the parent is a field. Namely, the remainder will have no terms in its p-adic expansion above the valuation of other.

EXAMPLES:

sage: R = Zp(3, 5)
sage: R(12).quo_rem(R(2))
(2*3 + O(3^6), 0)
sage: R(2).quo_rem(R(12))
(O(3^4), 2 + O(3^5))

sage: K = Qp(3, 5)
sage: K(12).quo_rem(K(2))
(2*3 + O(3^6), 0)
sage: K(2).quo_rem(K(12))
(2*3^-1 + 1 + 3 + 3^2 + 3^3 + O(3^4), 0)
>>> from sage.all import *
>>> R = Zp(Integer(3), Integer(5))
>>> R(Integer(12)).quo_rem(R(Integer(2)))
(2*3 + O(3^6), 0)
>>> R(Integer(2)).quo_rem(R(Integer(12)))
(O(3^4), 2 + O(3^5))

>>> K = Qp(Integer(3), Integer(5))
>>> K(Integer(12)).quo_rem(K(Integer(2)))
(2*3 + O(3^6), 0)
>>> K(Integer(2)).quo_rem(K(Integer(12)))
(2*3^-1 + 1 + 3 + 3^2 + 3^3 + O(3^4), 0)

You can get the same behavior for fields as for rings by using integral=True:

sage: K(12).quo_rem(K(2), integral=True)
(2*3 + O(3^6), 0)
sage: K(2).quo_rem(K(12), integral=True)
(O(3^4), 2 + O(3^5))
>>> from sage.all import *
>>> K(Integer(12)).quo_rem(K(Integer(2)), integral=True)
(2*3 + O(3^6), 0)
>>> K(Integer(2)).quo_rem(K(Integer(12)), integral=True)
(O(3^4), 2 + O(3^5))
slice(i, j, k=1, lift_mode='simple')[source]#

Returns the sum of the \(pi^{i + l \cdot k}\) terms of the series expansion of this element, where pi is the uniformizer, for \(i + l \cdot k\) between i and j-1 inclusive, and nonnegative integers \(l\). Behaves analogously to the slice function for lists.

INPUT:

  • i – an integer; if set to None, the sum will start with the first non-zero term of the series.

  • j – an integer; if set to None or \(\infty\), this method behaves as if it was set to the absolute precision of this element.

  • k – (default: 1) a positive integer

EXAMPLES:

sage: R = Zp(5, 6, 'capped-rel')
sage: a = R(1/2); a
3 + 2*5 + 2*5^2 + 2*5^3 + 2*5^4 + 2*5^5 + O(5^6)
sage: a.slice(2, 4)
2*5^2 + 2*5^3 + O(5^4)
sage: a.slice(1, 6, 2)
2*5 + 2*5^3 + 2*5^5 + O(5^6)
>>> from sage.all import *
>>> R = Zp(Integer(5), Integer(6), 'capped-rel')
>>> a = R(Integer(1)/Integer(2)); a
3 + 2*5 + 2*5^2 + 2*5^3 + 2*5^4 + 2*5^5 + O(5^6)
>>> a.slice(Integer(2), Integer(4))
2*5^2 + 2*5^3 + O(5^4)
>>> a.slice(Integer(1), Integer(6), Integer(2))
2*5 + 2*5^3 + 2*5^5 + O(5^6)

The step size k has to be positive:

sage: a.slice(0, 3, 0)
Traceback (most recent call last):
...
ValueError: slice step must be positive
sage: a.slice(0, 3, -1)
Traceback (most recent call last):
...
ValueError: slice step must be positive
>>> from sage.all import *
>>> a.slice(Integer(0), Integer(3), Integer(0))
Traceback (most recent call last):
...
ValueError: slice step must be positive
>>> a.slice(Integer(0), Integer(3), -Integer(1))
Traceback (most recent call last):
...
ValueError: slice step must be positive

If i exceeds j, then the result will be zero, with the precision given by j:

sage: a.slice(5, 4)
O(5^4)
sage: a.slice(6, 5)
O(5^5)
>>> from sage.all import *
>>> a.slice(Integer(5), Integer(4))
O(5^4)
>>> a.slice(Integer(6), Integer(5))
O(5^5)

However, the precision cannot exceed the precision of the element:

sage: a.slice(101,100)
O(5^6)
sage: a.slice(0,5,2)
3 + 2*5^2 + 2*5^4 + O(5^5)
sage: a.slice(0,6,2)
3 + 2*5^2 + 2*5^4 + O(5^6)
sage: a.slice(0,7,2)
3 + 2*5^2 + 2*5^4 + O(5^6)
>>> from sage.all import *
>>> a.slice(Integer(101),Integer(100))
O(5^6)
>>> a.slice(Integer(0),Integer(5),Integer(2))
3 + 2*5^2 + 2*5^4 + O(5^5)
>>> a.slice(Integer(0),Integer(6),Integer(2))
3 + 2*5^2 + 2*5^4 + O(5^6)
>>> a.slice(Integer(0),Integer(7),Integer(2))
3 + 2*5^2 + 2*5^4 + O(5^6)

If start is left blank, it is set to the valuation:

sage: K = Qp(5, 6)
sage: x = K(1/25 + 5); x
5^-2 + 5 + O(5^4)
sage: x.slice(None, 3)
5^-2 + 5 + O(5^3)
sage: x[:3]
doctest:warning
...
DeprecationWarning: __getitem__ is changing to match the behavior of number fields. Please use expansion instead.
See https://github.com/sagemath/sage/issues/14825 for details.
5^-2 + 5 + O(5^3)
>>> from sage.all import *
>>> K = Qp(Integer(5), Integer(6))
>>> x = K(Integer(1)/Integer(25) + Integer(5)); x
5^-2 + 5 + O(5^4)
>>> x.slice(None, Integer(3))
5^-2 + 5 + O(5^3)
>>> x[:Integer(3)]
doctest:warning
...
DeprecationWarning: __getitem__ is changing to match the behavior of number fields. Please use expansion instead.
See https://github.com/sagemath/sage/issues/14825 for details.
5^-2 + 5 + O(5^3)
sqrt(extend=True, all=False, algorithm=None)[source]#

Return the square root of this element.

INPUT:

  • self – a \(p\)-adic element.

  • extend – a boolean (default: True); if True, return a square root in an extension if necessary; if False and no root exists in the given ring or field, raise a ValueError.

  • all – a boolean (default: False); if True, return a list of all square roots.

  • algorithm"pari", "sage" or None (default: None); Sage provides an implementation for any extension of \(Q_p\) whereas only square roots over \(Q_p\) is implemented in Pari; the default is "pari" if the ground field is \(Q_p\), "sage" otherwise.

OUTPUT:

The square root or the list of all square roots of this element.

Note

The square root is chosen (resp. the square roots are ordered) in a deterministic way, which is compatible with change of precision.

EXAMPLES:

sage: R = Zp(3, 20)
sage: sqrt(R(0))
0

sage: sqrt(R(1))
1 + O(3^20)

sage: R(2).sqrt(extend=False)
Traceback (most recent call last):
...
ValueError: element is not a square

sage: s = sqrt(R(4)); -s
2 + O(3^20)

sage: s = sqrt(R(9)); s
3 + O(3^21)
>>> from sage.all import *
>>> R = Zp(Integer(3), Integer(20))
>>> sqrt(R(Integer(0)))
0

>>> sqrt(R(Integer(1)))
1 + O(3^20)

>>> R(Integer(2)).sqrt(extend=False)
Traceback (most recent call last):
...
ValueError: element is not a square

>>> s = sqrt(R(Integer(4))); -s
2 + O(3^20)

>>> s = sqrt(R(Integer(9))); s
3 + O(3^21)

Over the \(2\)-adics, the precision of the square root is less than the input:

sage: R2 = Zp(2, 20)
sage: sqrt(R2(1))
1 + O(2^19)
sage: sqrt(R2(4))
2 + O(2^20)

sage: R.<t> = Zq(2^10, 10)                                                  # needs sage.libs.ntl
sage: u = 1 + 8*t                                                           # needs sage.libs.ntl
sage: sqrt(u)                                                               # needs sage.libs.ntl
1 + t*2^2 + t^2*2^3 + t^2*2^4 + (t^4 + t^3 + t^2)*2^5 + (t^4 + t^2)*2^6
  + (t^5 + t^2)*2^7 + (t^6 + t^5 + t^4 + t^2)*2^8 + O(2^9)

sage: R.<a> = Zp(2).extension(x^3 - 2)
sage: u = R(1 + a^4 + a^5 + a^7 + a^8, 10); u
1 + a^4 + a^5 + a^7 + a^8 + O(a^10)
sage: v = sqrt(u); v                                                        # needs sage.libs.ntl
1 + a^2 + a^4 + a^6 + O(a^7)
>>> from sage.all import *
>>> R2 = Zp(Integer(2), Integer(20))
>>> sqrt(R2(Integer(1)))
1 + O(2^19)
>>> sqrt(R2(Integer(4)))
2 + O(2^20)

>>> R = Zq(Integer(2)**Integer(10), Integer(10), names=('t',)); (t,) = R._first_ngens(1)# needs sage.libs.ntl
>>> u = Integer(1) + Integer(8)*t                                                           # needs sage.libs.ntl
>>> sqrt(u)                                                               # needs sage.libs.ntl
1 + t*2^2 + t^2*2^3 + t^2*2^4 + (t^4 + t^3 + t^2)*2^5 + (t^4 + t^2)*2^6
  + (t^5 + t^2)*2^7 + (t^6 + t^5 + t^4 + t^2)*2^8 + O(2^9)

>>> R = Zp(Integer(2)).extension(x**Integer(3) - Integer(2), names=('a',)); (a,) = R._first_ngens(1)
>>> u = R(Integer(1) + a**Integer(4) + a**Integer(5) + a**Integer(7) + a**Integer(8), Integer(10)); u
1 + a^4 + a^5 + a^7 + a^8 + O(a^10)
>>> v = sqrt(u); v                                                        # needs sage.libs.ntl
1 + a^2 + a^4 + a^6 + O(a^7)

However, observe that the precision increases to its original value when we recompute the square of the square root:

sage: v^2                                                                   # needs sage.libs.ntl
1 + a^4 + a^5 + a^7 + a^8 + O(a^10)
>>> from sage.all import *
>>> v**Integer(2)                                                                   # needs sage.libs.ntl
1 + a^4 + a^5 + a^7 + a^8 + O(a^10)

If the input does not have enough precision in order to determine if the given element has a square root in the ground field, an error is raised:

sage: R(1, 6).sqrt()
Traceback (most recent call last):
...
PrecisionError: not enough precision to be sure that this element has a square root

sage: R(1, 7).sqrt()
1 + O(a^4)

sage: R(1+a^6, 7).sqrt(extend=False)
Traceback (most recent call last):
...
ValueError: element is not a square
>>> from sage.all import *
>>> R(Integer(1), Integer(6)).sqrt()
Traceback (most recent call last):
...
PrecisionError: not enough precision to be sure that this element has a square root

>>> R(Integer(1), Integer(7)).sqrt()
1 + O(a^4)

>>> R(Integer(1)+a**Integer(6), Integer(7)).sqrt(extend=False)
Traceback (most recent call last):
...
ValueError: element is not a square

In particular, an error is raised when we try to compute the square root of an inexact