\(p\)-adic Base Leaves

Implementations of \(\ZZ_p\) and \(\QQ_p\)

AUTHORS:

  • David Roe

  • Genya Zaytman: documentation

  • David Harvey: doctests

  • William Stein: doctest updates

EXAMPLES:

\(p\)-adic rings and fields are examples of inexact structures, as the reals are. That means that elements cannot generally be stored exactly: to do so would take an infinite amount of storage. Instead, we store an approximation to the elements with varying precision.

There are two types of precision for a \(p\)-adic element. The first is relative precision, which gives the number of known \(p\)-adic digits:

sage: R = Qp(5, 20, 'capped-rel', 'series'); a = R(675); a
2*5^2 + 5^4 + O(5^22)
sage: a.precision_relative()
20
>>> from sage.all import *
>>> R = Qp(Integer(5), Integer(20), 'capped-rel', 'series'); a = R(Integer(675)); a
2*5^2 + 5^4 + O(5^22)
>>> a.precision_relative()
20
R = Qp(5, 20, 'capped-rel', 'series'); a = R(675); a
a.precision_relative()

The second type of precision is absolute precision, which gives the power of \(p\) that this element is stored modulo:

sage: a.precision_absolute()
22
>>> from sage.all import *
>>> a.precision_absolute()
22
a.precision_absolute()

The number of times that \(p\) divides the element is called the valuation, and can be accessed with the methods valuation() and ordp():

sage: a.valuation() 2

The following relationship holds:

self.valuation() + self.precision_relative() == self.precision_absolute().

sage: a.valuation() + a.precision_relative() == a.precision_absolute() True

In the capped relative case, the relative precision of an element is restricted to be at most a certain value, specified at the creation of the field. Individual elements also store their own precision, so the effect of various arithmetic operations on precision is tracked. When you cast an exact element into a capped relative field, it truncates it to the precision cap of the field.:

sage: R = Qp(5, 5); a = R(4006); a
1 + 5 + 2*5^3 + 5^4 + O(5^5)
sage: b = R(17/3); b
4 + 2*5 + 3*5^2 + 5^3 + 3*5^4 + O(5^5)
sage: c = R(4025); c
5^2 + 2*5^3 + 5^4 + 5^5 + O(5^7)
sage: a + b
4*5 + 3*5^2 + 3*5^3 + 4*5^4 + O(5^5)
sage: a + b + c
4*5 + 4*5^2 + 5^4 + O(5^5)
>>> from sage.all import *
>>> R = Qp(Integer(5), Integer(5)); a = R(Integer(4006)); a
1 + 5 + 2*5^3 + 5^4 + O(5^5)
>>> b = R(Integer(17)/Integer(3)); b
4 + 2*5 + 3*5^2 + 5^3 + 3*5^4 + O(5^5)
>>> c = R(Integer(4025)); c
5^2 + 2*5^3 + 5^4 + 5^5 + O(5^7)
>>> a + b
4*5 + 3*5^2 + 3*5^3 + 4*5^4 + O(5^5)
>>> a + b + c
4*5 + 4*5^2 + 5^4 + O(5^5)
R = Qp(5, 5); a = R(4006); a
b = R(17/3); b
c = R(4025); c
a + b
a + b + c

sage: R = Zp(5, 5, 'capped-rel', 'series'); a = R(4006); a
1 + 5 + 2*5^3 + 5^4 + O(5^5)
sage: b = R(17/3); b
4 + 2*5 + 3*5^2 + 5^3 + 3*5^4 + O(5^5)
sage: c = R(4025); c
5^2 + 2*5^3 + 5^4 + 5^5 + O(5^7)
sage: a + b
4*5 + 3*5^2 + 3*5^3 + 4*5^4 + O(5^5)
sage: a + b + c
4*5 + 4*5^2 + 5^4 + O(5^5)
>>> from sage.all import *
>>> R = Zp(Integer(5), Integer(5), 'capped-rel', 'series'); a = R(Integer(4006)); a
1 + 5 + 2*5^3 + 5^4 + O(5^5)
>>> b = R(Integer(17)/Integer(3)); b
4 + 2*5 + 3*5^2 + 5^3 + 3*5^4 + O(5^5)
>>> c = R(Integer(4025)); c
5^2 + 2*5^3 + 5^4 + 5^5 + O(5^7)
>>> a + b
4*5 + 3*5^2 + 3*5^3 + 4*5^4 + O(5^5)
>>> a + b + c
4*5 + 4*5^2 + 5^4 + O(5^5)
R = Zp(5, 5, 'capped-rel', 'series'); a = R(4006); a
b = R(17/3); b
c = R(4025); c
a + b
a + b + c

In the capped absolute type, instead of having a cap on the relative precision of an element there is instead a cap on the absolute precision. Elements still store their own precisions, and as with the capped relative case, exact elements are truncated when cast into the ring.:

sage: R = ZpCA(5, 5); a = R(4005); a
5 + 2*5^3 + 5^4 + O(5^5)
sage: b = R(4025); b
5^2 + 2*5^3 + 5^4 + O(5^5)
sage: a * b
5^3 + 2*5^4 + O(5^5)
sage: (a * b) // 5^3
1 + 2*5 + O(5^2)
sage: type((a * b) // 5^3)
<class 'sage.rings.padics.padic_capped_absolute_element.pAdicCappedAbsoluteElement'>
sage: (a * b) / 5^3
1 + 2*5 + O(5^2)
sage: type((a * b) / 5^3)
<class 'sage.rings.padics.padic_capped_relative_element.pAdicCappedRelativeElement'>
>>> from sage.all import *
>>> R = ZpCA(Integer(5), Integer(5)); a = R(Integer(4005)); a
5 + 2*5^3 + 5^4 + O(5^5)
>>> b = R(Integer(4025)); b
5^2 + 2*5^3 + 5^4 + O(5^5)
>>> a * b
5^3 + 2*5^4 + O(5^5)
>>> (a * b) // Integer(5)**Integer(3)
1 + 2*5 + O(5^2)
>>> type((a * b) // Integer(5)**Integer(3))
<class 'sage.rings.padics.padic_capped_absolute_element.pAdicCappedAbsoluteElement'>
>>> (a * b) / Integer(5)**Integer(3)
1 + 2*5 + O(5^2)
>>> type((a * b) / Integer(5)**Integer(3))
<class 'sage.rings.padics.padic_capped_relative_element.pAdicCappedRelativeElement'>
R = ZpCA(5, 5); a = R(4005); a
b = R(4025); b
a * b
(a * b) // 5^3
type((a * b) // 5^3)
(a * b) / 5^3
type((a * b) / 5^3)

The fixed modulus type is the leanest of the \(p\)-adic rings: it is basically just a wrapper around \(\ZZ / p^n \ZZ\) providing a unified interface with the rest of the \(p\)-adics. This is the type you should use if your primary interest is in speed (though it’s not all that much faster than other \(p\)-adic types). It does not track precision of elements.:

sage: R = ZpFM(5, 5); a = R(4005); a
5 + 2*5^3 + 5^4
sage: a // 5
1 + 2*5^2 + 5^3
>>> from sage.all import *
>>> R = ZpFM(Integer(5), Integer(5)); a = R(Integer(4005)); a
5 + 2*5^3 + 5^4
>>> a // Integer(5)
1 + 2*5^2 + 5^3
R = ZpFM(5, 5); a = R(4005); a
a // 5

\(p\)-adic rings and fields should be created using the creation functions Zp() and Qp() as above. This will ensure that there is only one instance of \(\ZZ_p\) and \(\QQ_p\) of a given type, \(p\), print mode and precision. It also saves typing very long class names.:

sage: Qp(17,10)
17-adic Field with capped relative precision 10
sage: R = Qp(7, prec = 20, print_mode = 'val-unit'); S = Qp(7, prec = 20, print_mode = 'val-unit'); R is S
True
sage: Qp(2)
2-adic Field with capped relative precision 20
>>> from sage.all import *
>>> Qp(Integer(17),Integer(10))
17-adic Field with capped relative precision 10
>>> R = Qp(Integer(7), prec = Integer(20), print_mode = 'val-unit'); S = Qp(Integer(7), prec = Integer(20), print_mode = 'val-unit'); R is S
True
>>> Qp(Integer(2))
2-adic Field with capped relative precision 20
Qp(17,10)
R = Qp(7, prec = 20, print_mode = 'val-unit'); S = Qp(7, prec = 20, print_mode = 'val-unit'); R is S
Qp(2)

Once one has a \(p\)-adic ring or field, one can cast elements into it in the standard way. Integers, ints, longs, Rationals, other \(p\)-adic types, pari \(p\)-adics and elements of \(\ZZ / p^n \ZZ\) can all be cast into a \(p\)-adic field.:

sage: R = Qp(5, 5, 'capped-rel','series'); a = R(16); a
1 + 3*5 + O(5^5)
sage: b = R(23/15); b
5^-1 + 3 + 3*5 + 5^2 + 3*5^3 + O(5^4)
sage: S = Zp(5, 5, 'fixed-mod','val-unit'); c = S(Mod(75,125)); c
5^2 * 3
sage: R(c)
3*5^2 + O(5^5)
>>> from sage.all import *
>>> R = Qp(Integer(5), Integer(5), 'capped-rel','series'); a = R(Integer(16)); a
1 + 3*5 + O(5^5)
>>> b = R(Integer(23)/Integer(15)); b
5^-1 + 3 + 3*5 + 5^2 + 3*5^3 + O(5^4)
>>> S = Zp(Integer(5), Integer(5), 'fixed-mod','val-unit'); c = S(Mod(Integer(75),Integer(125))); c
5^2 * 3
>>> R(c)
3*5^2 + O(5^5)
R = Qp(5, 5, 'capped-rel','series'); a = R(16); a
b = R(23/15); b
S = Zp(5, 5, 'fixed-mod','val-unit'); c = S(Mod(75,125)); c
R(c)

In the previous example, since fixed-mod elements don’t keep track of their precision, we assume that it has the full precision of the ring. This is why you have to cast manually here.

While you can cast explicitly as above, the chains of automatic coercion are more restricted. As always in Sage, the following arrows are transitive and the diagram is commutative.:

int -> long -> Integer -> Zp capped-rel -> Zp capped_abs -> IntegerMod
Integer -> Zp fixed-mod -> IntegerMod
Integer -> Zp capped-abs -> Qp capped-rel

In addition, there are arrows within each type. For capped relative and capped absolute rings and fields, these arrows go from lower precision cap to higher precision cap. This works since elements track their own precision: choosing the parent with higher precision cap means that precision is less likely to be truncated unnecessarily. For fixed modulus parents, the arrow goes from higher precision cap to lower. The fact that elements do not track precision necessitates this choice in order to not produce incorrect results.

class sage.rings.padics.padic_base_leaves.pAdicFieldCappedRelative(p, prec, print_mode, names)[source]

Bases: pAdicFieldBaseGeneric, pAdicCappedRelativeFieldGeneric

An implementation of \(p\)-adic fields with capped relative precision.

EXAMPLES:

sage: K = Qp(17, 1000000)  # indirect doctest
sage: K = Qp(101)  # indirect doctest
>>> from sage.all import *
>>> K = Qp(Integer(17), Integer(1000000))  # indirect doctest
>>> K = Qp(Integer(101))  # indirect doctest
K = Qp(17, 1000000)  # indirect doctest
K = Qp(101)  # indirect doctest
random_element(algorithm='default')[source]

Return a random element of self, optionally using the algorithm argument to decide how it generates the element. Algorithms currently implemented:

  • 'default': Choose an integer \(k\) using the standard distribution on the integers. Then choose an integer \(a\) uniformly in the range \(0 \le a < p^N\) where \(N\) is the precision cap of self. Return self(p^k * a, absprec = k + self.precision_cap()).

EXAMPLES:

sage: Qp(17,6).random_element().parent() is Qp(17,6)
True
>>> from sage.all import *
>>> Qp(Integer(17),Integer(6)).random_element().parent() is Qp(Integer(17),Integer(6))
True
Qp(17,6).random_element().parent() is Qp(17,6)
class sage.rings.padics.padic_base_leaves.pAdicFieldFloatingPoint(p, prec, print_mode, names)[source]

Bases: pAdicFieldBaseGeneric, pAdicFloatingPointFieldGeneric

An implementation of the \(p\)-adic rationals with floating point precision.

class sage.rings.padics.padic_base_leaves.pAdicFieldLattice(p, prec, subtype, print_mode, names, label=None)[source]

Bases: pAdicLatticeGeneric, pAdicFieldBaseGeneric

An implementation of the \(p\)-adic numbers with lattice precision.

INPUT:

  • p – prime

  • prec – precision cap, given as a pair (relative_cap, absolute_cap)

  • subtype – either 'cap' or 'float'

  • print_mode – dictionary with print options

  • names – how to print the prime

  • label – the label of this ring

See also

label()

EXAMPLES:

sage: R = QpLC(next_prime(10^60))  # indirect doctest
doctest:...: FutureWarning: This class/method/function is marked as experimental.
It, its functionality or its interface might change without a formal deprecation.
See https://github.com/sagemath/sage/issues/23505 for details.
sage: type(R)
<class 'sage.rings.padics.padic_base_leaves.pAdicFieldLattice_with_category'>

sage: R = QpLC(2,label='init')  # indirect doctest
sage: R
2-adic Field with lattice-cap precision (label: init)
>>> from sage.all import *
>>> R = QpLC(next_prime(Integer(10)**Integer(60)))  # indirect doctest
doctest:...: FutureWarning: This class/method/function is marked as experimental.
It, its functionality or its interface might change without a formal deprecation.
See https://github.com/sagemath/sage/issues/23505 for details.
>>> type(R)
<class 'sage.rings.padics.padic_base_leaves.pAdicFieldLattice_with_category'>

>>> R = QpLC(Integer(2),label='init')  # indirect doctest
>>> R
2-adic Field with lattice-cap precision (label: init)
R = QpLC(next_prime(10^60))  # indirect doctest
type(R)
R = QpLC(2,label='init')  # indirect doctest
R
random_element(prec=None, integral=False)[source]

Return a random element of this ring.

INPUT:

  • prec – integer or None (default); the absolute precision of the generated random element

  • integral – boolean (default: False); if True, return an element in the ring of integers

EXAMPLES:

sage: K = QpLC(2)
sage: K.random_element()   # not tested, known bug (see :issue:`32126`)
2^-8 + 2^-7 + 2^-6 + 2^-5 + 2^-3 + 1 + 2^2 + 2^3 + 2^5 + O(2^12)
sage: K.random_element(integral=True)    # random
2^3 + 2^4 + 2^5 + 2^6 + 2^7 + 2^10 + 2^11 + 2^14 + 2^15 + 2^16
 + 2^17 + 2^18 + 2^19 + O(2^20)

sage: K.random_element(prec=10)    # random
2^(-3) + 1 + 2 + 2^4 + 2^8 + O(2^10)
>>> from sage.all import *
>>> K = QpLC(Integer(2))
>>> K.random_element()   # not tested, known bug (see :issue:`32126`)
2^-8 + 2^-7 + 2^-6 + 2^-5 + 2^-3 + 1 + 2^2 + 2^3 + 2^5 + O(2^12)
>>> K.random_element(integral=True)    # random
2^3 + 2^4 + 2^5 + 2^6 + 2^7 + 2^10 + 2^11 + 2^14 + 2^15 + 2^16
 + 2^17 + 2^18 + 2^19 + O(2^20)

>>> K.random_element(prec=Integer(10))    # random
2^(-3) + 1 + 2 + 2^4 + 2^8 + O(2^10)
K = QpLC(2)
K.random_element()   # not tested, known bug (see :issue:`32126`)
K.random_element(integral=True)    # random
K.random_element(prec=10)    # random

If the given precision is higher than the internal cap of the parent, then the cap is used:

sage: K.precision_cap_relative()
20
sage: K.random_element(prec=100)    # random
2^5 + 2^8 + 2^11 + 2^12 + 2^14 + 2^18 + 2^20 + 2^24 + O(2^25)
>>> from sage.all import *
>>> K.precision_cap_relative()
20
>>> K.random_element(prec=Integer(100))    # random
2^5 + 2^8 + 2^11 + 2^12 + 2^14 + 2^18 + 2^20 + 2^24 + O(2^25)
K.precision_cap_relative()
K.random_element(prec=100)    # random
class sage.rings.padics.padic_base_leaves.pAdicFieldRelaxed(p, prec, print_mode, names)[source]

Bases: pAdicRelaxedGeneric, pAdicFieldBaseGeneric

An implementation of relaxed arithmetics over \(\QQ_p\).

INPUT:

  • p – prime

  • prec – default precision

  • print_mode – dictionary with print options

  • names – how to print the prime

EXAMPLES:

sage: R = QpER(5)  # indirect doctest                                           # needs sage.libs.flint
sage: type(R)                                                                   # needs sage.libs.flint
<class 'sage.rings.padics.padic_base_leaves.pAdicFieldRelaxed_with_category'>
>>> from sage.all import *
>>> R = QpER(Integer(5))  # indirect doctest                                           # needs sage.libs.flint
>>> type(R)                                                                   # needs sage.libs.flint
<class 'sage.rings.padics.padic_base_leaves.pAdicFieldRelaxed_with_category'>
R = QpER(5)  # indirect doctest                                           # needs sage.libs.flint
type(R)                                                                   # needs sage.libs.flint
class sage.rings.padics.padic_base_leaves.pAdicRingCappedAbsolute(p, prec, print_mode, names)[source]

Bases: pAdicRingBaseGeneric, pAdicCappedAbsoluteRingGeneric

An implementation of the \(p\)-adic integers with capped absolute precision.

class sage.rings.padics.padic_base_leaves.pAdicRingCappedRelative(p, prec, print_mode, names)[source]

Bases: pAdicRingBaseGeneric, pAdicCappedRelativeRingGeneric

An implementation of the \(p\)-adic integers with capped relative precision.

class sage.rings.padics.padic_base_leaves.pAdicRingFixedMod(p, prec, print_mode, names)[source]

Bases: pAdicRingBaseGeneric, pAdicFixedModRingGeneric

An implementation of the \(p\)-adic integers using fixed modulus.

class sage.rings.padics.padic_base_leaves.pAdicRingFloatingPoint(p, prec, print_mode, names)[source]

Bases: pAdicRingBaseGeneric, pAdicFloatingPointRingGeneric

An implementation of the \(p\)-adic integers with floating point precision.

class sage.rings.padics.padic_base_leaves.pAdicRingLattice(p, prec, subtype, print_mode, names, label=None)[source]

Bases: pAdicLatticeGeneric, pAdicRingBaseGeneric

An implementation of the \(p\)-adic integers with lattice precision.

INPUT:

  • p – prime

  • prec – precision cap, given as a pair (relative_cap, absolute_cap)

  • subtype – either 'cap' or 'float'

  • print_mode – dictionary with print options

  • names – how to print the prime

  • label – the label of this ring

See also

label()

EXAMPLES:

sage: R = ZpLC(next_prime(10^60))  # indirect doctest
doctest:...: FutureWarning: This class/method/function is marked as experimental.
It, its functionality or its interface might change without a formal deprecation.
See https://github.com/sagemath/sage/issues/23505 for details.
sage: type(R)
<class 'sage.rings.padics.padic_base_leaves.pAdicRingLattice_with_category'>

sage: R = ZpLC(2, label='init')  # indirect doctest
sage: R
2-adic Ring with lattice-cap precision (label: init)
>>> from sage.all import *
>>> R = ZpLC(next_prime(Integer(10)**Integer(60)))  # indirect doctest
doctest:...: FutureWarning: This class/method/function is marked as experimental.
It, its functionality or its interface might change without a formal deprecation.
See https://github.com/sagemath/sage/issues/23505 for details.
>>> type(R)
<class 'sage.rings.padics.padic_base_leaves.pAdicRingLattice_with_category'>

>>> R = ZpLC(Integer(2), label='init')  # indirect doctest
>>> R
2-adic Ring with lattice-cap precision (label: init)
R = ZpLC(next_prime(10^60))  # indirect doctest
type(R)
R = ZpLC(2, label='init')  # indirect doctest
R
random_element(prec=None)[source]

Return a random element of this ring.

INPUT:

  • prec – integer or None (default); the absolute precision of the generated random element

EXAMPLES:

sage: R = ZpLC(2)
sage: R.random_element()    # random
2^3 + 2^4 + 2^5 + 2^6 + 2^7 + 2^10 + 2^11 + 2^14 + 2^15 + 2^16
 + 2^17 + 2^18 + 2^19 + 2^21 + O(2^23)

sage: R.random_element(prec=10)    # random
1 + 2^3 + 2^4 + 2^7 + O(2^10)
>>> from sage.all import *
>>> R = ZpLC(Integer(2))
>>> R.random_element()    # random
2^3 + 2^4 + 2^5 + 2^6 + 2^7 + 2^10 + 2^11 + 2^14 + 2^15 + 2^16
 + 2^17 + 2^18 + 2^19 + 2^21 + O(2^23)

>>> R.random_element(prec=Integer(10))    # random
1 + 2^3 + 2^4 + 2^7 + O(2^10)
R = ZpLC(2)
R.random_element()    # random
R.random_element(prec=10)    # random
class sage.rings.padics.padic_base_leaves.pAdicRingRelaxed(p, prec, print_mode, names)[source]

Bases: pAdicRelaxedGeneric, pAdicRingBaseGeneric

An implementation of relaxed arithmetics over \(\ZZ_p\).

INPUT:

  • p – prime

  • prec – default precision

  • print_mode – dictionary with print options

  • names – how to print the prime

EXAMPLES:

sage: R = ZpER(5)  # indirect doctest                                           # needs sage.libs.flint
sage: type(R)                                                                   # needs sage.libs.flint
<class 'sage.rings.padics.padic_base_leaves.pAdicRingRelaxed_with_category'>
>>> from sage.all import *
>>> R = ZpER(Integer(5))  # indirect doctest                                           # needs sage.libs.flint
>>> type(R)                                                                   # needs sage.libs.flint
<class 'sage.rings.padics.padic_base_leaves.pAdicRingRelaxed_with_category'>
R = ZpER(5)  # indirect doctest                                           # needs sage.libs.flint
type(R)                                                                   # needs sage.libs.flint