p-Adic Generic Element
Elements of
-Adic Rings and Fields
AUTHORS:
Bases: sage.rings.padics.local_generic_element.LocalGenericElement
INPUT:
Return the
-adic absolute value of self.
This is normalized so that the absolute value of
is
.
INPUT:
EXAMPLES:
sage: a = Qp(5)(15); a.abs()
1/5
sage: a.abs(53)
0.200000000000000
sage: Qp(7)(0).abs()
0
sage: Qp(7)(0).abs(prec=20)
0.00000
An unramified extension:
sage: R = Zp(5,5)
sage: P.<x> = PolynomialRing(R)
sage: Z25.<u> = R.ext(x^2 - 3)
sage: u.abs()
1
sage: (u^24-1).abs()
1/5
A ramified extension:
sage: W.<w> = R.ext(x^5 + 75*x^3 - 15*x^2 + 125*x - 5)
sage: w.abs()
0.724779663677696
sage: W(0).abs()
0.000000000000000
Returns the additive order of self, where self is considered
to be zero if it is zero modulo
.
INPUT:
OUTPUT:
integer – the additive order of self
EXAMPLES:
sage: R = Zp(7, 4, 'capped-rel', 'series'); a = R(7^3); a.additive_order(3)
1
sage: a.additive_order(4)
+Infinity
sage: R = Zp(7, 4, 'fixed-mod', 'series'); a = R(7^5); a.additive_order(6)
1
Returns a polynomial of degree at most
which is approximately
satisfied by this number. Note that the returned polynomial need not be
irreducible, and indeed usually won’t be if this number is a good
approximation to an algebraic number of degree less than
.
ALGORITHM: Uses the PARI C-library algdep command.
INPUT:
OUTPUT:
polynomial – degree n polynomial approximately satisfied by self
EXAMPLES:
sage: K = Qp(3,20,'capped-rel','series'); R = Zp(3,20,'capped-rel','series')
sage: a = K(7/19); a
1 + 2*3 + 3^2 + 3^3 + 2*3^4 + 2*3^5 + 3^8 + 2*3^9 + 3^11 + 3^12 + 2*3^15 + 2*3^16 + 3^17 + 2*3^19 + O(3^20)
sage: a.algdep(1)
19*x - 7
sage: K2 = Qp(7,20,'capped-rel')
sage: b = K2.zeta(); b.algdep(2)
x^2 - x + 1
sage: K2 = Qp(11,20,'capped-rel')
sage: b = K2.zeta(); b.algdep(4)
x^4 - x^3 + x^2 - x + 1
sage: a = R(7/19); a
1 + 2*3 + 3^2 + 3^3 + 2*3^4 + 2*3^5 + 3^8 + 2*3^9 + 3^11 + 3^12 + 2*3^15 + 2*3^16 + 3^17 + 2*3^19 + O(3^20)
sage: a.algdep(1)
19*x - 7
sage: R2 = Zp(7,20,'capped-rel')
sage: b = R2.zeta(); b.algdep(2)
x^2 - x + 1
sage: R2 = Zp(11,20,'capped-rel')
sage: b = R2.zeta(); b.algdep(4)
x^4 - x^3 + x^2 - x + 1
Returns a polynomial of degree at most
which is approximately
satisfied by this number. Note that the returned polynomial need not
be irreducible, and indeed usually won’t be if this number is a good
approximation to an algebraic number of degree less than
.
ALGORITHM: Uses the PARI C-library algdep command.
INPUT:
OUTPUT:
polynomial – degree n polynomial approximately satisfied by self
EXAMPLES:
sage: K = Qp(3,20,'capped-rel','series'); R = Zp(3,20,'capped-rel','series')
sage: a = K(7/19); a
1 + 2*3 + 3^2 + 3^3 + 2*3^4 + 2*3^5 + 3^8 + 2*3^9 + 3^11 + 3^12 + 2*3^15 + 2*3^16 + 3^17 + 2*3^19 + O(3^20)
sage: a.algebraic_dependency(1)
19*x - 7
sage: K2 = Qp(7,20,'capped-rel')
sage: b = K2.zeta(); b.algebraic_dependency(2)
x^2 - x + 1
sage: K2 = Qp(11,20,'capped-rel')
sage: b = K2.zeta(); b.algebraic_dependency(4)
x^4 - x^3 + x^2 - x + 1
sage: a = R(7/19); a
1 + 2*3 + 3^2 + 3^3 + 2*3^4 + 2*3^5 + 3^8 + 2*3^9 + 3^11 + 3^12 + 2*3^15 + 2*3^16 + 3^17 + 2*3^19 + O(3^20)
sage: a.algebraic_dependency(1)
19*x - 7
sage: R2 = Zp(7,20,'capped-rel')
sage: b = R2.zeta(); b.algebraic_dependency(2)
x^2 - x + 1
sage: R2 = Zp(11,20,'capped-rel')
sage: b = R2.zeta(); b.algebraic_dependency(4)
x^4 - x^3 + x^2 - x + 1
Compute the
-adic exponential of this element if the exponential
series converges.
INPUT:
ALGORITHM: If self has a lift method (which should happen for
elements of
and
), then one uses the rule:
modulo the precision. The value of
is
precomputed. Otherwise, use the power series expansion of
,
evaluating a certain number of terms which does about
multiplications.
EXAMPLES:
log() and exp() are inverse to each other:
sage: Z13 = Zp(13, 10)
sage: a = Z13(14); a
1 + 13 + O(13^10)
sage: a.log().exp()
1 + 13 + O(13^10)
An error occurs if this is called with an element for which the exponential series does not converge:
sage: Z13.one().exp()
Traceback (most recent call last):
...
ValueError: Exponential does not converge for that input.
The next few examples illustrate precision when computing
-adic
exponentials:
sage: R = Zp(5,10)
sage: e = R(2*5 + 2*5**2 + 4*5**3 + 3*5**4 + 5**5 + 3*5**7 + 2*5**8 + 4*5**9).add_bigoh(10); e
2*5 + 2*5^2 + 4*5^3 + 3*5^4 + 5^5 + 3*5^7 + 2*5^8 + 4*5^9 + O(5^10)
sage: e.exp()*R.teichmuller(4)
4 + 2*5 + 3*5^3 + O(5^10)
sage: K = Qp(5,10)
sage: e = K(2*5 + 2*5**2 + 4*5**3 + 3*5**4 + 5**5 + 3*5**7 + 2*5**8 + 4*5**9).add_bigoh(10); e
2*5 + 2*5^2 + 4*5^3 + 3*5^4 + 5^5 + 3*5^7 + 2*5^8 + 4*5^9 + O(5^10)
sage: e.exp()*K.teichmuller(4)
4 + 2*5 + 3*5^3 + O(5^10)
Logarithms and exponentials in extension fields. First, in an Eisenstein extension:
sage: R = Zp(5,5)
sage: S.<x> = R[]
sage: f = x^4 + 15*x^2 + 625*x - 5
sage: W.<w> = R.ext(f)
sage: z = 1 + w^2 + 4*w^7; z
1 + w^2 + 4*w^7 + O(w^20)
sage: z.log().exp()
1 + w^2 + 4*w^7 + O(w^20)
Now an unramified example:
sage: R = Zp(5,5)
sage: S.<x> = R[]
sage: g = x^3 + 3*x + 3
sage: A.<a> = R.ext(g)
sage: b = 1 + 5*(1 + a^2) + 5^3*(3 + 2*a); b
1 + (a^2 + 1)*5 + (2*a + 3)*5^3 + O(5^5)
sage: b.log().exp()
1 + (a^2 + 1)*5 + (2*a + 3)*5^3 + O(5^5)
TESTS:
Check that results are consistent over a range of precision:
sage: max_prec = 40
sage: p = 3
sage: K = Zp(p, max_prec)
sage: full_exp = (K(p)).exp()
sage: for prec in range(2, max_prec):
... ll = (K(p).add_bigoh(prec)).exp()
... assert ll == full_exp
... assert ll.precision_absolute() == prec
sage: K = Qp(p, max_prec)
sage: full_exp = (K(p)).exp()
sage: for prec in range(2, max_prec):
... ll = (K(p).add_bigoh(prec)).exp()
... assert ll == full_exp
... assert ll.precision_absolute() == prec
Check that this also works for capped-absolute implementations:
sage: Z13 = ZpCA(13, 10)
sage: a = Z13(14); a
1 + 13 + O(13^10)
sage: a.log().exp()
1 + 13 + O(13^10)
sage: R = ZpCA(5,5)
sage: S.<x> = R[]
sage: f = x^4 + 15*x^2 + 625*x - 5
sage: W.<w> = R.ext(f)
sage: z = 1 + w^2 + 4*w^7; z
1 + w^2 + 4*w^7 + O(w^16)
sage: z.log().exp()
1 + w^2 + 4*w^7 + O(w^16)
Check that this also works for fixed-mod implementations:
sage: Z13 = ZpFM(13, 10)
sage: a = Z13(14); a
1 + 13 + O(13^10)
sage: a.log().exp()
1 + 13 + O(13^10)
sage: R = ZpFM(5,5)
sage: S.<x> = R[]
sage: f = x^4 + 15*x^2 + 625*x - 5
sage: W.<w> = R.ext(f)
sage: z = 1 + w^2 + 4*w^7; z
1 + w^2 + 4*w^7 + O(w^20)
sage: z.log().exp()
1 + w^2 + 4*w^7 + O(w^20)
Some corner cases:
sage: Z2 = Zp(2, 5)
sage: Z2(2).exp()
Traceback (most recent call last):
...
ValueError: Exponential does not converge for that input.
sage: S.<x> = Z2[]
sage: W.<w> = Z2.ext(x^3-2)
sage: (w^2).exp()
Traceback (most recent call last):
...
ValueError: Exponential does not converge for that input.
sage: (w^3).exp()
Traceback (most recent call last):
...
ValueError: Exponential does not converge for that input.
sage: (w^4).exp()
1 + w^4 + w^5 + w^7 + w^9 + w^10 + w^14 + O(w^15)
AUTHORS:
Return a greatest common divisor of this element and other.
INPUT:
OUTPUT:
The gcd of self and other.
AUTHORS:
Note
Since the elements are only given with finite precision,
their greatest common divisor is in general not unique (not even up
to units). For example
is a representative for the elements
0 and 3 in the 3-adic ring
. The greatest common
divisior of
and
could be (among others) 3 or 0 which
have different valuation. The algorithm implemented here, will
return an element of minimal valuation among the possible greatest
common divisors.
EXAMPLES:
The greatest common divisor is either zero or a power of the uniformizing parameter:
sage: R = Zp(3)
sage: R.zero().gcd(R.zero())
0
sage: R(3).gcd(9)
3 + O(3^21)
A non-zero result is always lifted to the maximal precision possible in the ring:
sage: a = R(3,2); a
3 + O(3^2)
sage: b = R(9,3); b
3^2 + O(3^3)
sage: a.gcd(b)
3 + O(3^21)
sage: a.gcd(0)
3 + O(3^21)
If both elements are zero, then the result is zero with the precision set to the smallest of their precisions:
sage: a = R.zero(); a
0
sage: b = R(0,2); b
O(3^2)
sage: a.gcd(b)
O(3^2)
One could argue that it is mathematically correct to return 9 + O(3^22) instead. However, this would lead to some confusing behaviour:
sage: alternative_gcd = R(9,22); alternative_gcd
3^2 + O(3^22)
sage: a.is_zero()
True
sage: b.is_zero()
True
sage: alternative_gcd.is_zero()
False
If exactly one element is zero, then the result depends on the valuation of the other element:
sage: R(0,3).gcd(3^4)
O(3^3)
sage: R(0,4).gcd(3^4)
O(3^4)
sage: R(0,5).gcd(3^4)
3^4 + O(3^24)
Over a field, the greatest common divisor is either zero (possibly with finite precision) or one:
sage: K = Qp(3)
sage: K(3).gcd(0)
1 + O(3^20)
sage: K.zero().gcd(0)
0
sage: K.zero().gcd(K(0,2))
O(3^2)
sage: K(3).gcd(4)
1 + O(3^20)
TESTS:
The implementation also works over extensions:
sage: K = Qp(3)
sage: R.<a> = K[]
sage: L.<a> = K.extension(a^3-3)
sage: (a+3).gcd(3)
1 + O(a^60)
sage: R = Zp(3)
sage: S.<a> = R[]
sage: S.<a> = R.extension(a^3-3)
sage: (a+3).gcd(3)
a + O(a^61)
sage: K = Qp(3)
sage: R.<a> = K[]
sage: L.<a> = K.extension(a^2-2)
sage: (a+3).gcd(3)
1 + O(3^20)
sage: R = Zp(3)
sage: S.<a> = R[]
sage: S.<a> = R.extension(a^2-2)
sage: (a+3).gcd(3)
1 + O(3^20)
For elements with a fixed modulus:
sage: R = ZpFM(3)
sage: R(3).gcd(9)
3 + O(3^20)
And elements with a capped absolute precision:
sage: R = ZpCA(3)
sage: R(3).gcd(9)
3 + O(3^20)
Returns whether self is a square
INPUT:
OUTPUT:
boolean – whether self is a square
EXAMPLES:
sage: R = Zp(3,20,'capped-rel')
sage: R(0).is_square()
True
sage: R(1).is_square()
True
sage: R(2).is_square()
False
TESTS:
sage: R(3).is_square()
False
sage: R(4).is_square()
True
sage: R(6).is_square()
False
sage: R(9).is_square()
True
sage: R2 = Zp(2,20,'capped-rel')
sage: R2(0).is_square()
True
sage: R2(1).is_square()
True
sage: R2(2).is_square()
False
sage: R2(3).is_square()
False
sage: R2(4).is_square()
True
sage: R2(5).is_square()
False
sage: R2(6).is_square()
False
sage: R2(7).is_square()
False
sage: R2(8).is_square()
False
sage: R2(9).is_square()
True
sage: K = Qp(3,20,'capped-rel')
sage: K(0).is_square()
True
sage: K(1).is_square()
True
sage: K(2).is_square()
False
sage: K(3).is_square()
False
sage: K(4).is_square()
True
sage: K(6).is_square()
False
sage: K(9).is_square()
True
sage: K(1/3).is_square()
False
sage: K(1/9).is_square()
True
sage: K2 = Qp(2,20,'capped-rel')
sage: K2(0).is_square()
True
sage: K2(1).is_square()
True
sage: K2(2).is_square()
False
sage: K2(3).is_square()
False
sage: K2(4).is_square()
True
sage: K2(5).is_square()
False
sage: K2(6).is_square()
False
sage: K2(7).is_square()
False
sage: K2(8).is_square()
False
sage: K2(9).is_square()
True
sage: K2(1/2).is_square()
False
sage: K2(1/4).is_square()
True
Compute the
-adic logarithm of this element.
The usual power series for the logarithm with values in the additive
group of a
-adic ring only converges for 1-units (units congruent to
1 modulo
). However, there is a unique extension of the logarithm
to a homomorphism defined on all the units: If
is a
unit with
and
a Teichmuller representative,
then we define
. This is the correct extension
because the units
split as a product
, where
is the subgroup of 1-units and
is a fundamental
root of unity. The
factor is torsion, so must go
to 0 under any homomorphism to the fraction field, which is a torsion
free group.
INPUTS:
to branch.
-adic field, however, for most neighborhoods of 1, it lies
in the ring of integers. This flag decides if the codomain should be
the same as the input (default) or if it should change to the
fraction field of the input.NOTES:
What some other systems do:
Todo
There is a soft-linear time algorith for logarithm described by Dan Berstein at http://cr.yp.to/lineartime/multapps-20041007.pdf
ALGORITHM:
of the input.2. Raise
to
where
is the inertia degree of the ring
extension, to obtain a 1-unit.

to compute the logarithm
.
EXAMPLES:
sage: Z13 = Zp(13, 10)
sage: a = Z13(14); a
1 + 13 + O(13^10)
Note that the relative precision decreases when we take log – it is the absolute precision that is preserved:
sage: a.log()
13 + 6*13^2 + 2*13^3 + 5*13^4 + 10*13^6 + 13^7 + 11*13^8 + 8*13^9 + O(13^10)
sage: Q13 = Qp(13, 10)
sage: a = Q13(14); a
1 + 13 + O(13^10)
sage: a.log()
13 + 6*13^2 + 2*13^3 + 5*13^4 + 10*13^6 + 13^7 + 11*13^8 + 8*13^9 + O(13^10)
The next few examples illustrate precision when computing
-adic
logarithms:
sage: R = Zp(5,10)
sage: e = R(389); e
4 + 2*5 + 3*5^3 + O(5^10)
sage: e.log()
2*5 + 2*5^2 + 4*5^3 + 3*5^4 + 5^5 + 3*5^7 + 2*5^8 + 4*5^9 + O(5^10)
sage: K = Qp(5,10)
sage: e = K(389); e
4 + 2*5 + 3*5^3 + O(5^10)
sage: e.log()
2*5 + 2*5^2 + 4*5^3 + 3*5^4 + 5^5 + 3*5^7 + 2*5^8 + 4*5^9 + O(5^10)
The logarithm is not only defined for 1-units:
sage: R = Zp(5,10)
sage: a = R(2)
sage: a.log()
2*5 + 3*5^2 + 2*5^3 + 4*5^4 + 2*5^6 + 2*5^7 + 4*5^8 + 2*5^9 + O(5^10)
If you want to take the logarithm of a non-unit you must specify either p_branch or pi_branch:
sage: b = R(5)
sage: b.log()
Traceback (most recent call last):
...
ValueError: You must specify a branch of the logarithm for non-units
sage: b.log(p_branch=4)
4 + O(5^10)
sage: c = R(10)
sage: c.log(p_branch=4)
4 + 2*5 + 3*5^2 + 2*5^3 + 4*5^4 + 2*5^6 + 2*5^7 + 4*5^8 + 2*5^9 + O(5^10)
The branch parameters are only relevant for elements of non-zero valuation:
sage: a.log(p_branch=0)
2*5 + 3*5^2 + 2*5^3 + 4*5^4 + 2*5^6 + 2*5^7 + 4*5^8 + 2*5^9 + O(5^10)
sage: a.log(p_branch=1)
2*5 + 3*5^2 + 2*5^3 + 4*5^4 + 2*5^6 + 2*5^7 + 4*5^8 + 2*5^9 + O(5^10)
Logarithms can also be computed in extension fields. First, in an Eisenstein extension:
sage: R = Zp(5,5)
sage: S.<x> = ZZ[]
sage: f = x^4 + 15*x^2 + 625*x - 5
sage: W.<w> = R.ext(f)
sage: z = 1 + w^2 + 4*w^7; z
1 + w^2 + 4*w^7 + O(w^20)
sage: z.log()
w^2 + 2*w^4 + 3*w^6 + 4*w^7 + w^9 + 4*w^10 + 4*w^11 + 4*w^12 + 3*w^14 + w^15 + w^17 + 3*w^18 + 3*w^19 + O(w^20)
In an extension, there will usually be a difference between specifying p_branch and pi_branch:
sage: b = W(5)
sage: b.log()
Traceback (most recent call last):
...
ValueError: You must specify a branch of the logarithm for non-units
sage: b.log(p_branch=0)
O(w^20)
sage: b.log(p_branch=w)
w + O(w^20)
sage: b.log(pi_branch=0)
3*w^2 + 2*w^4 + 2*w^6 + 3*w^8 + 4*w^10 + w^13 + w^14 + 2*w^15 + 2*w^16 + w^18 + 4*w^19 + O(w^20)
sage: b.unit_part().log()
3*w^2 + 2*w^4 + 2*w^6 + 3*w^8 + 4*w^10 + w^13 + w^14 + 2*w^15 + 2*w^16 + w^18 + 4*w^19 + O(w^20)
sage: y = w^2 * 4*w^7; y
4*w^9 + O(w^29)
sage: y.log(p_branch=0)
2*w^2 + 2*w^4 + 2*w^6 + 2*w^8 + w^10 + w^12 + 4*w^13 + 4*w^14 + 3*w^15 + 4*w^16 + 4*w^17 + w^18 + 4*w^19 + O(w^20)
sage: y.log(p_branch=w)
w + 2*w^2 + 2*w^4 + 4*w^5 + 2*w^6 + 2*w^7 + 2*w^8 + 4*w^9 + w^10 + 3*w^11 + w^12 + 4*w^14 + 4*w^16 + 2*w^17 + w^19 + O(w^20)
Check that log is multiplicative:
sage: y.log(p_branch=0) + z.log() - (y*z).log(p_branch=0)
O(w^20)
Now an unramified example:
sage: g = x^3 + 3*x + 3
sage: A.<a> = R.ext(g)
sage: b = 1 + 5*(1 + a^2) + 5^3*(3 + 2*a)
sage: b.log()
(a^2 + 1)*5 + (3*a^2 + 4*a + 2)*5^2 + (3*a^2 + 2*a)*5^3 + (3*a^2 + 2*a + 2)*5^4 + O(5^5)
Check that log is multiplicative:
sage: c = 3 + 5^2*(2 + 4*a)
sage: b.log() + c.log() - (b*c).log()
O(5^5)
We illustrate the effect of the precision argument:
sage: R = ZpCA(7,10)
sage: x = R(41152263); x
5 + 3*7^2 + 4*7^3 + 3*7^4 + 5*7^5 + 6*7^6 + 7^9 + O(7^10)
sage: x.log(aprec = 5)
7 + 3*7^2 + 4*7^3 + 3*7^4 + O(7^5)
sage: x.log(aprec = 7)
7 + 3*7^2 + 4*7^3 + 3*7^4 + 7^5 + 3*7^6 + O(7^7)
sage: x.log()
7 + 3*7^2 + 4*7^3 + 3*7^4 + 7^5 + 3*7^6 + 7^7 + 3*7^8 + 4*7^9 + O(7^10)
The logarithm is not defined for zero:
sage: R.zero().log()
Traceback (most recent call last):
...
ValueError: logarithm is not defined at zero
For elements in a
-adic ring, the logarithm will be returned in the
same ring:
sage: x = R(2)
sage: x.log().parent()
7-adic Ring with capped absolute precision 10
sage: x = R(14)
sage: x.log(p_branch=0).parent()
7-adic Ring with capped absolute precision 10
This is not possible if the logarithm has negative valuation:
sage: R = ZpCA(3,10)
sage: S.<x> = R[]
sage: f = x^3 - 3
sage: W.<w> = R.ext(f)
sage: w.log(p_branch=2)
Traceback (most recent call last):
...
ValueError: logarithm is not integral, use change_frac=True to obtain a result in the fraction field
sage: w.log(p_branch=2, change_frac=True)
2*w^-3 + O(w^21)
TESTS:
Check that results are consistent over a range of precision:
sage: max_prec = 40
sage: p = 3
sage: K = Zp(p, max_prec)
sage: full_log = (K(1 + p)).log()
sage: for prec in range(2, max_prec):
... ll1 = (K(1+p).add_bigoh(prec)).log()
... ll2 = K(1+p).log(prec)
... assert ll1 == full_log
... assert ll2 == full_log
... assert ll1.precision_absolute() == prec
Check that aprec works for fixed-mod elements:
sage: R = ZpFM(7,10)
sage: x = R(41152263); x
5 + 3*7^2 + 4*7^3 + 3*7^4 + 5*7^5 + 6*7^6 + 7^9 + O(7^10)
sage: x.log(aprec = 5)
7 + 3*7^2 + 4*7^3 + 3*7^4 + O(7^10)
sage: x.log(aprec = 7)
7 + 3*7^2 + 4*7^3 + 3*7^4 + 7^5 + 3*7^6 + O(7^10)
sage: x.log()
7 + 3*7^2 + 4*7^3 + 3*7^4 + 7^5 + 3*7^6 + 7^7 + 3*7^8 + 4*7^9 + O(7^10)
Check that precision is computed correctly in highly ramified extensions:
sage: S.<x> = ZZ[]
sage: K = Qp(5,5)
sage: f = x^625 - 5*x - 5
sage: W.<w> = K.extension(f)
sage: z = 1 - w^2 + O(w^11)
sage: x = 1 - z
sage: z.log().precision_absolute()
-975
sage: (x^5/5).precision_absolute()
-570
sage: (x^25/25).precision_absolute()
-975
sage: (x^125/125).precision_absolute()
-775
sage: z = 1 - w + O(w^2)
sage: x = 1 - z
sage: z.log().precision_absolute()
-1625
sage: (x^5/5).precision_absolute()
-615
sage: (x^25/25).precision_absolute()
-1200
sage: (x^125/125).precision_absolute()
-1625
sage: (x^625/625).precision_absolute()
-1250
sage: z.log().precision_relative()
250
AUTHORS:
-adic class
-adic rings.Returns a minimal polynomial of this
-adic element, i.e., x - self
INPUT:
-adic elementEXAMPLES:
sage: Zp(5,5)(1/3).minimal_polynomial('x')
(1 + O(5^5))*x + (3 + 5 + 3*5^2 + 5^3 + 3*5^4 + O(5^5))
Returns the multiplicative order of self, where self is
considered to be one if it is one modulo
.
INPUT:
OUTPUT:
EXAMPLES:
sage: K = Qp(5,20,'capped-rel')
sage: K(-1).multiplicative_order(20)
2
sage: K(1).multiplicative_order(20)
1
sage: K(2).multiplicative_order(20)
+Infinity
sage: K(3).multiplicative_order(20)
+Infinity
sage: K(4).multiplicative_order(20)
+Infinity
sage: K(5).multiplicative_order(20)
+Infinity
sage: K(25).multiplicative_order(20)
+Infinity
sage: K(1/5).multiplicative_order(20)
+Infinity
sage: K(1/25).multiplicative_order(20)
+Infinity
sage: K.zeta().multiplicative_order(20)
4
sage: R = Zp(5,20,'capped-rel')
sage: R(-1).multiplicative_order(20)
2
sage: R(1).multiplicative_order(20)
1
sage: R(2).multiplicative_order(20)
+Infinity
sage: R(3).multiplicative_order(20)
+Infinity
sage: R(4).multiplicative_order(20)
+Infinity
sage: R(5).multiplicative_order(20)
+Infinity
sage: R(25).multiplicative_order(20)
+Infinity
sage: R.zeta().multiplicative_order(20)
4
Returns the norm of this
-adic element over the ground ring.
Warning
This is not the
-adic absolute value. This is a field
theoretic norm down to a ground ring. If you want the
-adic absolute value, use the abs() function
instead.
INPUT:
EXAMPLES:
sage: Zp(5)(5).norm()
5 + O(5^21)
Returns the valuation of self, normalized so that the valuation of
is 1
INPUT:
NOTE: The optional argument p is used for consistency with the valuation methods on integer and rational.
OUTPUT:
integer – the valuation of self, normalized so that the valuation of
is 1
EXAMPLES:
sage: R = Zp(5,20,'capped-rel')
sage: R(0).ordp()
+Infinity
sage: R(1).ordp()
0
sage: R(2).ordp()
0
sage: R(5).ordp()
1
sage: R(10).ordp()
1
sage: R(25).ordp()
2
sage: R(50).ordp()
2
sage: R(1/2).ordp()
0
Returns a rational approximation to this p-adic number
INPUT:
OUTPUT:
rational – an approximation to self
EXAMPLES:
sage: R = Zp(5,20,'capped-rel')
sage: for i in range(11):
... for j in range(1,10):
... if j == 5:
... continue
... assert i/j == R(i/j).rational_reconstruction()
Returns the square root of this p-adic number
INPUT:
OUTPUT:
p-adic element – the square root of this p-adic number
If all=False, the square root chosen is the one whose
reduction mod
is in the range
.
EXAMPLES:
sage: R = Zp(3,20,'capped-rel', 'val-unit')
sage: R(0).square_root()
0
sage: R(1).square_root()
1 + O(3^20)
sage: R(2).square_root(extend = False)
Traceback (most recent call last):
...
ValueError: element is not a square
sage: R(4).square_root() == R(-2)
True
sage: R(9).square_root()
3 * 1 + O(3^21)
When p = 2, the precision of the square root is one less than the input:
sage: R2 = Zp(2,20,'capped-rel')
sage: R2(0).square_root()
0
sage: R2(1).square_root()
1 + O(2^19)
sage: R2(4).square_root()
2 + O(2^20)
sage: R2(9).square_root() == R2(3, 19) or R2(9).square_root() == R2(-3, 19)
True
sage: R2(17).square_root()
1 + 2^3 + 2^5 + 2^6 + 2^7 + 2^9 + 2^10 + 2^13 + 2^16 + 2^17 + O(2^19)
sage: R3 = Zp(5,20,'capped-rel')
sage: R3(0).square_root()
0
sage: R3(1).square_root()
1 + O(5^20)
sage: R3(-1).square_root() == R3.teichmuller(2) or R3(-1).square_root() == R3.teichmuller(3)
True
TESTS:
sage: R = Qp(3,20,'capped-rel')
sage: R(0).square_root()
0
sage: R(1).square_root()
1 + O(3^20)
sage: R(4).square_root() == R(-2)
True
sage: R(9).square_root()
3 + O(3^21)
sage: R(1/9).square_root()
3^-1 + O(3^19)
sage: R2 = Qp(2,20,'capped-rel')
sage: R2(0).square_root()
0
sage: R2(1).square_root()
1 + O(2^19)
sage: R2(4).square_root()
2 + O(2^20)
sage: R2(9).square_root() == R2(3,19) or R2(9).square_root() == R2(-3,19)
True
sage: R2(17).square_root()
1 + 2^3 + 2^5 + 2^6 + 2^7 + 2^9 + 2^10 + 2^13 + 2^16 + 2^17 + O(2^19)
sage: R3 = Qp(5,20,'capped-rel')
sage: R3(0).square_root()
0
sage: R3(1).square_root()
1 + O(5^20)
sage: R3(-1).square_root() == R3.teichmuller(2) or R3(-1).square_root() == R3.teichmuller(3)
True
sage: R = Zp(3,20,'capped-abs')
sage: R(1).square_root()
1 + O(3^20)
sage: R(4).square_root() == R(-2)
True
sage: R(9).square_root()
3 + O(3^19)
sage: R2 = Zp(2,20,'capped-abs')
sage: R2(1).square_root()
1 + O(2^19)
sage: R2(4).square_root()
2 + O(2^18)
sage: R2(9).square_root() == R2(3) or R2(9).square_root() == R2(-3)
True
sage: R2(17).square_root()
1 + 2^3 + 2^5 + 2^6 + 2^7 + 2^9 + 2^10 + 2^13 + 2^16 + 2^17 + O(2^19)
sage: R3 = Zp(5,20,'capped-abs')
sage: R3(1).square_root()
1 + O(5^20)
sage: R3(-1).square_root() == R3.teichmuller(2) or R3(-1).square_root() == R3.teichmuller(3)
True
Returns a string representation of self.
EXAMPLES:
sage: Zp(5,5,print_mode='bars')(1/3).str()[3:]
'1|3|1|3|2'
Returns the trace of this
-adic element over the ground ring
INPUT:
OUTPUT:
-adic element over the
ground ringEXAMPLES:
sage: Zp(5,5)(5).trace()
5 + O(5^6)
Return (self.valuation(), self.unit_part()). To be overridden in derived classes.
EXAMPLES:
sage: Zp(5,5)(5).val_unit()
(1, 1 + O(5^5))
Returns the valuation of this element.
INPUT:
NOTE: The optional argument p is used for consistency with the valuation methods on integer and rational.
OUTPUT:
integer – the valuation of self
EXAMPLES:
sage: R = Zp(17, 4,'capped-rel')
sage: a = R(2*17^2)
sage: a.valuation()
2
sage: R = Zp(5, 4,'capped-rel')
sage: R(0).valuation()
+Infinity
TESTS:
sage: R(1).valuation()
0
sage: R(2).valuation()
0
sage: R(5).valuation()
1
sage: R(10).valuation()
1
sage: R(25).valuation()
2
sage: R(50).valuation()
2
sage: R = Qp(17, 4)
sage: a = R(2*17^2)
sage: a.valuation()
2
sage: R = Qp(5, 4)
sage: R(0).valuation()
+Infinity
sage: R(1).valuation()
0
sage: R(2).valuation()
0
sage: R(5).valuation()
1
sage: R(10).valuation()
1
sage: R(25).valuation()
2
sage: R(50).valuation()
2
sage: R(1/2).valuation()
0
sage: R(1/5).valuation()
-1
sage: R(1/10).valuation()
-1
sage: R(1/25).valuation()
-2
sage: R(1/50).valuation()
-2
sage: K.<a> = Qq(25)
sage: K(0).valuation()
+Infinity
sage: R(1/50).valuation(5)
-2
sage: R(1/50).valuation(3)
Traceback (most recent call last):
...
ValueError: Ring (5-adic Field with capped relative precision 4) residue field of the wrong characteristic.