\(p\)-adic ZZ_pX
FM Element#
This file implements elements of Eisenstein and unramified extensions of \(\ZZ_p\) with fixed modulus precision.
For the parent class see padic_extension_leaves.pyx
.
The underlying implementation is through NTL’s ZZ_pX
class. Each
element contains the following data:
value
(ZZ_pX_c
) – An ntlZZ_pX
storing the value. The variable \(x\) is the uniformizer in the case of Eisenstein extensions. ThisZZ_pX
is created with global ntl modulus determined by the parent’s precision cap and shared among all elements.prime_pow
(some subclass ofPowComputer_ZZ_pX
) – a class, identical among all elements with the same parent, holding common data.prime_pow.deg
– the degree of the extensionprime_pow.e
– the ramification indexprime_pow.f
– the inertia degreeprime_pow.prec_cap
– the unramified precision cap: for Eisenstein extensions this is the smallest power of \(p\) that is zeroprime_pow.ram_prec_cap
– the ramified precision cap: for Eisenstein extensions this will be the smallest power of \(x\) that is indistinguishable from zeroprime_pow.pow_ZZ_tmp
, prime_pow.pow_mpz_t_tmp``,prime_pow.pow_Integer
– functions for accessing powers of \(p\). The first two return pointers. Seesage/rings/padics/pow_computer_ext
for examples and important warnings.prime_pow.get_context
,prime_pow.get_context_capdiv
,prime_pow.get_top_context
– obtain anntl_ZZ_pContext_class
corresponding to \(p^n\). The capdiv version divides byprime_pow.e
as appropriate.top_context
corresponds to \(p^{\texttt{prec\_cap}}\).prime_pow.restore_context
,prime_pow.restore_context_capdiv
,prime_pow.restore_top_context
– restores the given contextprime_pow.get_modulus
,get_modulus_capdiv
,get_top_modulus
– Returns aZZ_pX_Modulus_c*
pointing to a polynomial modulus defined modulo \(p^n\) (appropriately divided byprime_pow.e
in the capdiv case).
EXAMPLES:
An Eisenstein extension:
sage: R = ZpFM(5,5)
sage: S.<x> = R[]
sage: f = x^5 + 75*x^3 - 15*x^2 + 125*x - 5
sage: W.<w> = R.ext(f); W
5-adic Eisenstein Extension Ring in w defined by x^5 + 75*x^3 - 15*x^2 + 125*x - 5
sage: z = (1+w)^5; z
1 + w^5 + w^6 + 2*w^7 + 4*w^8 + 3*w^10 + w^12 + 4*w^13 + 4*w^14 + 4*w^15 + 4*w^16 + 4*w^17 + 4*w^20 + w^21 + 4*w^24
sage: y = z >> 1; y
w^4 + w^5 + 2*w^6 + 4*w^7 + 3*w^9 + w^11 + 4*w^12 + 4*w^13 + 4*w^14 + 4*w^15 + 4*w^16 + 4*w^19 + w^20 + 4*w^23 + 4*w^24
sage: y.valuation()
4
sage: y.precision_relative()
21
sage: y.precision_absolute()
25
sage: z - (y << 1)
1
>>> from sage.all import *
>>> R = ZpFM(Integer(5),Integer(5))
>>> S = R['x']; (x,) = S._first_ngens(1)
>>> f = x**Integer(5) + Integer(75)*x**Integer(3) - Integer(15)*x**Integer(2) + Integer(125)*x - Integer(5)
>>> W = R.ext(f, names=('w',)); (w,) = W._first_ngens(1); W
5-adic Eisenstein Extension Ring in w defined by x^5 + 75*x^3 - 15*x^2 + 125*x - 5
>>> z = (Integer(1)+w)**Integer(5); z
1 + w^5 + w^6 + 2*w^7 + 4*w^8 + 3*w^10 + w^12 + 4*w^13 + 4*w^14 + 4*w^15 + 4*w^16 + 4*w^17 + 4*w^20 + w^21 + 4*w^24
>>> y = z >> Integer(1); y
w^4 + w^5 + 2*w^6 + 4*w^7 + 3*w^9 + w^11 + 4*w^12 + 4*w^13 + 4*w^14 + 4*w^15 + 4*w^16 + 4*w^19 + w^20 + 4*w^23 + 4*w^24
>>> y.valuation()
4
>>> y.precision_relative()
21
>>> y.precision_absolute()
25
>>> z - (y << Integer(1))
1
An unramified extension:
sage: # needs sage.libs.flint
sage: g = x^3 + 3*x + 3
sage: A.<a> = R.ext(g)
sage: z = (1+a)^5; z
(2*a^2 + 4*a) + (3*a^2 + 3*a + 1)*5 + (4*a^2 + 3*a + 4)*5^2 + (4*a^2 + 4*a + 4)*5^3 + (4*a^2 + 4*a + 4)*5^4
sage: z - 1 - 5*a - 10*a^2 - 10*a^3 - 5*a^4 - a^5
0
sage: y = z >> 1; y
(3*a^2 + 3*a + 1) + (4*a^2 + 3*a + 4)*5 + (4*a^2 + 4*a + 4)*5^2 + (4*a^2 + 4*a + 4)*5^3
sage: 1/a
(3*a^2 + 4) + (a^2 + 4)*5 + (3*a^2 + 4)*5^2 + (a^2 + 4)*5^3 + (3*a^2 + 4)*5^4
>>> from sage.all import *
>>> # needs sage.libs.flint
>>> g = x**Integer(3) + Integer(3)*x + Integer(3)
>>> A = R.ext(g, names=('a',)); (a,) = A._first_ngens(1)
>>> z = (Integer(1)+a)**Integer(5); z
(2*a^2 + 4*a) + (3*a^2 + 3*a + 1)*5 + (4*a^2 + 3*a + 4)*5^2 + (4*a^2 + 4*a + 4)*5^3 + (4*a^2 + 4*a + 4)*5^4
>>> z - Integer(1) - Integer(5)*a - Integer(10)*a**Integer(2) - Integer(10)*a**Integer(3) - Integer(5)*a**Integer(4) - a**Integer(5)
0
>>> y = z >> Integer(1); y
(3*a^2 + 3*a + 1) + (4*a^2 + 3*a + 4)*5 + (4*a^2 + 4*a + 4)*5^2 + (4*a^2 + 4*a + 4)*5^3
>>> Integer(1)/a
(3*a^2 + 4) + (a^2 + 4)*5 + (3*a^2 + 4)*5^2 + (a^2 + 4)*5^3 + (3*a^2 + 4)*5^4
Different printing modes:
sage: # needs sage.libs.flint
sage: R = ZpFM(5, print_mode='digits'); S.<x> = R[]; f = x^5 + 75*x^3 - 15*x^2 + 125*x -5; W.<w> = R.ext(f)
sage: z = (1+w)^5; repr(z)
'...4110403113210310442221311242000111011201102002023303214332011214403232013144001400444441030421100001'
sage: R = ZpFM(5, print_mode='bars'); S.<x> = R[]; g = x^3 + 3*x + 3; A.<a> = R.ext(g)
sage: z = (1+a)^5; repr(z)
'...[4, 4, 4]|[4, 4, 4]|[4, 4, 4]|[4, 4, 4]|[4, 4, 4]|[4, 4, 4]|[4, 4, 4]|[4, 4, 4]|[4, 4, 4]|[4, 4, 4]|[4, 4, 4]|[4, 4, 4]|[4, 4, 4]|[4, 4, 4]|[4, 4, 4]|[4, 4, 4]|[4, 4, 4]|[4, 3, 4]|[1, 3, 3]|[0, 4, 2]'
sage: R = ZpFM(5, print_mode='terse'); S.<x> = R[]; f = x^5 + 75*x^3 - 15*x^2 + 125*x -5; W.<w> = R.ext(f)
sage: z = (1+w)^5; z
6 + 95367431640505*w + 25*w^2 + 95367431640560*w^3 + 5*w^4
sage: R = ZpFM(5, print_mode='val-unit'); S.<x> = R[]; f = x^5 + 75*x^3 - 15*x^2 + 125*x -5; W.<w> = R.ext(f)
sage: y = (1+w)^5 - 1; y
w^5 * (2090041 + 19073486126901*w + 1258902*w^2 + 57220458985049*w^3 + 16785*w^4)
>>> from sage.all import *
>>> # needs sage.libs.flint
>>> R = ZpFM(Integer(5), print_mode='digits'); S = R['x']; (x,) = S._first_ngens(1); f = x**Integer(5) + Integer(75)*x**Integer(3) - Integer(15)*x**Integer(2) + Integer(125)*x -Integer(5); W = R.ext(f, names=('w',)); (w,) = W._first_ngens(1)
>>> z = (Integer(1)+w)**Integer(5); repr(z)
'...4110403113210310442221311242000111011201102002023303214332011214403232013144001400444441030421100001'
>>> R = ZpFM(Integer(5), print_mode='bars'); S = R['x']; (x,) = S._first_ngens(1); g = x**Integer(3) + Integer(3)*x + Integer(3); A = R.ext(g, names=('a',)); (a,) = A._first_ngens(1)
>>> z = (Integer(1)+a)**Integer(5); repr(z)
'...[4, 4, 4]|[4, 4, 4]|[4, 4, 4]|[4, 4, 4]|[4, 4, 4]|[4, 4, 4]|[4, 4, 4]|[4, 4, 4]|[4, 4, 4]|[4, 4, 4]|[4, 4, 4]|[4, 4, 4]|[4, 4, 4]|[4, 4, 4]|[4, 4, 4]|[4, 4, 4]|[4, 4, 4]|[4, 3, 4]|[1, 3, 3]|[0, 4, 2]'
>>> R = ZpFM(Integer(5), print_mode='terse'); S = R['x']; (x,) = S._first_ngens(1); f = x**Integer(5) + Integer(75)*x**Integer(3) - Integer(15)*x**Integer(2) + Integer(125)*x -Integer(5); W = R.ext(f, names=('w',)); (w,) = W._first_ngens(1)
>>> z = (Integer(1)+w)**Integer(5); z
6 + 95367431640505*w + 25*w^2 + 95367431640560*w^3 + 5*w^4
>>> R = ZpFM(Integer(5), print_mode='val-unit'); S = R['x']; (x,) = S._first_ngens(1); f = x**Integer(5) + Integer(75)*x**Integer(3) - Integer(15)*x**Integer(2) + Integer(125)*x -Integer(5); W = R.ext(f, names=('w',)); (w,) = W._first_ngens(1)
>>> y = (Integer(1)+w)**Integer(5) - Integer(1); y
w^5 * (2090041 + 19073486126901*w + 1258902*w^2 + 57220458985049*w^3 + 16785*w^4)
AUTHORS:
David Roe (2008-01-01) initial version
- sage.rings.padics.padic_ZZ_pX_FM_element.make_ZZpXFMElement(parent, f)[source]#
Create a new
pAdicZZpXFMElement
out of anntl_ZZ_pX
f
, with parentparent
. For use with pickling.EXAMPLES:
sage: R = ZpFM(5,5) sage: S.<x> = R[] sage: f = x^5 + 75*x^3 - 15*x^2 + 125*x - 5 sage: W.<w> = R.ext(f) sage: z = (1 + w)^5 - 1 sage: loads(dumps(z)) == z # indirect doctest True
>>> from sage.all import * >>> R = ZpFM(Integer(5),Integer(5)) >>> S = R['x']; (x,) = S._first_ngens(1) >>> f = x**Integer(5) + Integer(75)*x**Integer(3) - Integer(15)*x**Integer(2) + Integer(125)*x - Integer(5) >>> W = R.ext(f, names=('w',)); (w,) = W._first_ngens(1) >>> z = (Integer(1) + w)**Integer(5) - Integer(1) >>> loads(dumps(z)) == z # indirect doctest True
- class sage.rings.padics.padic_ZZ_pX_FM_element.pAdicZZpXFMElement[source]#
Bases:
pAdicZZpXElement
Creates an element of a fixed modulus, unramified or eisenstein extension of \(\ZZ_p\) or \(\QQ_p\).
INPUT:
parent
– either anEisensteinRingFixedMod
orUnramifiedRingFixedMod
x
– an integer, rational, \(p\)-adic element, polynomial, list, integer_mod, pari int/frac/poly_t/pol_mod, anntl_ZZ_pX
, anntl_ZZX
, anntl_ZZ
, or anntl_ZZ_p
absprec
– not usedrelprec
– not usedempty
– whether to return after initializing to zero (without setting anything)
EXAMPLES:
sage: R = ZpFM(5,5) sage: S.<x> = R[] sage: f = x^5 + 75*x^3 - 15*x^2 + 125*x - 5 sage: W.<w> = R.ext(f) sage: z = (1+w)^5; z # indirect doctest 1 + w^5 + w^6 + 2*w^7 + 4*w^8 + 3*w^10 + w^12 + 4*w^13 + 4*w^14 + 4*w^15 + 4*w^16 + 4*w^17 + 4*w^20 + w^21 + 4*w^24
>>> from sage.all import * >>> R = ZpFM(Integer(5),Integer(5)) >>> S = R['x']; (x,) = S._first_ngens(1) >>> f = x**Integer(5) + Integer(75)*x**Integer(3) - Integer(15)*x**Integer(2) + Integer(125)*x - Integer(5) >>> W = R.ext(f, names=('w',)); (w,) = W._first_ngens(1) >>> z = (Integer(1)+w)**Integer(5); z # indirect doctest 1 + w^5 + w^6 + 2*w^7 + 4*w^8 + 3*w^10 + w^12 + 4*w^13 + 4*w^14 + 4*w^15 + 4*w^16 + 4*w^17 + 4*w^20 + w^21 + 4*w^24
- add_bigoh(absprec)[source]#
Return a new element truncated modulo \(\pi^{\text{absprec}}\).
This is only implemented for unramified extension at this point.
INPUT:
absprec
– an integer
OUTPUT:
A new element truncated modulo \(\pi^{\mbox{absprec}}\).
EXAMPLES:
sage: R = Zp(7,4,'fixed-mod') sage: a = R(1 + 7 + 7^2) sage: a.add_bigoh(1) 1
>>> from sage.all import * >>> R = Zp(Integer(7),Integer(4),'fixed-mod') >>> a = R(Integer(1) + Integer(7) + Integer(7)**Integer(2)) >>> a.add_bigoh(Integer(1)) 1
- expansion(n=None, lift_mode='simple')[source]#
Return a list giving a series representation of this element.
If
lift_mode == 'simple' or 'smallest'
, the returned list will consist ofintegers (in the eisenstein case) or
lists of integers (in the unramified case).
this element can be reconstructed as
a sum of elements of the list times powers of the uniformiser (in the eisenstein case), or
as a sum of powers of the \(p\) times polynomials in the generator (in the unramified case).
If
lift_mode == 'simple'
, all integers will be in the range \([0,p-1]\),If
lift_mode == 'smallest'
they will be in the range \([(1-p)/2, p/2]\).If
lift_mode == 'teichmuller'
, returns a list ofpAdicZZpXCRElements
, all of which are Teichmuller representatives and such that this element is the sum of that list times powers of the uniformizer.
INPUT:
n
– integer (defaultNone
); if given, returns the corresponding entry in the expansion
EXAMPLES:
sage: R = ZpFM(5,5) sage: S.<x> = R[] sage: f = x^5 + 75*x^3 - 15*x^2 + 125*x - 5 sage: W.<w> = R.ext(f) sage: y = W(775); y w^10 + 4*w^12 + 2*w^14 + w^15 + 2*w^16 + 4*w^17 + w^18 + w^20 + 2*w^21 + 3*w^22 + w^23 + w^24 sage: (y>>9).expansion() [0, 1, 0, 4, 0, 2, 1, 2, 4, 1, 0, 1, 2, 3, 1, 1, 4, 1, 2, 4, 1, 0, 0, 3] sage: (y>>9).expansion(lift_mode='smallest') [0, 1, 0, -1, 0, 2, 1, 2, 0, 1, 2, 1, 1, -1, -1, 2, -2, 0, -2, -2, -2, 0, -2, -2, 2] sage: w^10 - w^12 + 2*w^14 + w^15 + 2*w^16 + w^18 + 2*w^19 + w^20 + w^21 - w^22 - w^23 + 2*w^24 w^10 + 4*w^12 + 2*w^14 + w^15 + 2*w^16 + 4*w^17 + w^18 + w^20 + 2*w^21 + 3*w^22 + w^23 + w^24 sage: g = x^3 + 3*x + 3 sage: # needs sage.libs.flint sage: A.<a> = R.ext(g) sage: y = 75 + 45*a + 1200*a^2; y 4*a*5 + (3*a^2 + a + 3)*5^2 + 4*a^2*5^3 + a^2*5^4 sage: E = y.expansion(); E 5-adic expansion of 4*a*5 + (3*a^2 + a + 3)*5^2 + 4*a^2*5^3 + a^2*5^4 sage: list(E) [[], [0, 4], [3, 1, 3], [0, 0, 4], [0, 0, 1]] sage: list(y.expansion(lift_mode='smallest')) [[], [0, -1], [-2, 2, -2], [1], [0, 0, 2]] sage: 5*((-2*5 + 25) + (-1 + 2*5)*a + (-2*5 + 2*125)*a^2) 4*a*5 + (3*a^2 + a + 3)*5^2 + 4*a^2*5^3 + a^2*5^4 sage: W(0).expansion() [] sage: list(A(0,4).expansion()) []
>>> from sage.all import * >>> R = ZpFM(Integer(5),Integer(5)) >>> S = R['x']; (x,) = S._first_ngens(1) >>> f = x**Integer(5) + Integer(75)*x**Integer(3) - Integer(15)*x**Integer(2) + Integer(125)*x - Integer(5) >>> W = R.ext(f, names=('w',)); (w,) = W._first_ngens(1) >>> y = W(Integer(775)); y w^10 + 4*w^12 + 2*w^14 + w^15 + 2*w^16 + 4*w^17 + w^18 + w^20 + 2*w^21 + 3*w^22 + w^23 + w^24 >>> (y>>Integer(9)).expansion() [0, 1, 0, 4, 0, 2, 1, 2, 4, 1, 0, 1, 2, 3, 1, 1, 4, 1, 2, 4, 1, 0, 0, 3] >>> (y>>Integer(9)).expansion(lift_mode='smallest') [0, 1, 0, -1, 0, 2, 1, 2, 0, 1, 2, 1, 1, -1, -1, 2, -2, 0, -2, -2, -2, 0, -2, -2, 2] >>> w**Integer(10) - w**Integer(12) + Integer(2)*w**Integer(14) + w**Integer(15) + Integer(2)*w**Integer(16) + w**Integer(18) + Integer(2)*w**Integer(19) + w**Integer(20) + w**Integer(21) - w**Integer(22) - w**Integer(23) + Integer(2)*w**Integer(24) w^10 + 4*w^12 + 2*w^14 + w^15 + 2*w^16 + 4*w^17 + w^18 + w^20 + 2*w^21 + 3*w^22 + w^23 + w^24 >>> g = x**Integer(3) + Integer(3)*x + Integer(3) >>> # needs sage.libs.flint >>> A = R.ext(g, names=('a',)); (a,) = A._first_ngens(1) >>> y = Integer(75) + Integer(45)*a + Integer(1200)*a**Integer(2); y 4*a*5 + (3*a^2 + a + 3)*5^2 + 4*a^2*5^3 + a^2*5^4 >>> E = y.expansion(); E 5-adic expansion of 4*a*5 + (3*a^2 + a + 3)*5^2 + 4*a^2*5^3 + a^2*5^4 >>> list(E) [[], [0, 4], [3, 1, 3], [0, 0, 4], [0, 0, 1]] >>> list(y.expansion(lift_mode='smallest')) [[], [0, -1], [-2, 2, -2], [1], [0, 0, 2]] >>> Integer(5)*((-Integer(2)*Integer(5) + Integer(25)) + (-Integer(1) + Integer(2)*Integer(5))*a + (-Integer(2)*Integer(5) + Integer(2)*Integer(125))*a**Integer(2)) 4*a*5 + (3*a^2 + a + 3)*5^2 + 4*a^2*5^3 + a^2*5^4 >>> W(Integer(0)).expansion() [] >>> list(A(Integer(0),Integer(4)).expansion()) []
Check that Issue #25879 has been resolved:
sage: K = ZpCA(3,5) sage: R.<a> = K[] sage: L.<a> = K.extension(a^2 - 3) sage: a.residue() 0
>>> from sage.all import * >>> K = ZpCA(Integer(3),Integer(5)) >>> R = K['a']; (a,) = R._first_ngens(1) >>> L = K.extension(a**Integer(2) - Integer(3), names=('a',)); (a,) = L._first_ngens(1) >>> a.residue() 0
- is_equal_to(right, absprec=None)[source]#
Return whether
self
is equal toright
moduloself.uniformizer()^absprec
.If
absprec
isNone
, returns ifself
is equal toright
modulo the precision cap.EXAMPLES:
sage: R = Zp(5,5) sage: S.<x> = R[] sage: f = x^5 + 75*x^3 - 15*x^2 + 125*x - 5 sage: W.<w> = R.ext(f) sage: a = W(47); b = W(47 + 25) sage: a.is_equal_to(b) False sage: a.is_equal_to(b, 7) True
>>> from sage.all import * >>> R = Zp(Integer(5),Integer(5)) >>> S = R['x']; (x,) = S._first_ngens(1) >>> f = x**Integer(5) + Integer(75)*x**Integer(3) - Integer(15)*x**Integer(2) + Integer(125)*x - Integer(5) >>> W = R.ext(f, names=('w',)); (w,) = W._first_ngens(1) >>> a = W(Integer(47)); b = W(Integer(47) + Integer(25)) >>> a.is_equal_to(b) False >>> a.is_equal_to(b, Integer(7)) True
- is_zero(absprec=None)[source]#
Return whether the valuation of
self
is at leastabsprec
; ifabsprec
isNone
, return whetherself
is indistinguishable from zero.EXAMPLES:
sage: R = ZpFM(5,5) sage: S.<x> = R[] sage: f = x^5 + 75*x^3 - 15*x^2 + 125*x - 5 sage: W.<w> = R.ext(f) sage: O(w^189).is_zero() True sage: W(0).is_zero() True sage: a = W(675) sage: a.is_zero() False sage: a.is_zero(7) True sage: a.is_zero(21) False
>>> from sage.all import * >>> R = ZpFM(Integer(5),Integer(5)) >>> S = R['x']; (x,) = S._first_ngens(1) >>> f = x**Integer(5) + Integer(75)*x**Integer(3) - Integer(15)*x**Integer(2) + Integer(125)*x - Integer(5) >>> W = R.ext(f, names=('w',)); (w,) = W._first_ngens(1) >>> O(w**Integer(189)).is_zero() True >>> W(Integer(0)).is_zero() True >>> a = W(Integer(675)) >>> a.is_zero() False >>> a.is_zero(Integer(7)) True >>> a.is_zero(Integer(21)) False
- lift_to_precision(absprec=None)[source]#
Return
self
.EXAMPLES:
sage: R = ZpFM(5,5) sage: S.<x> = R[] sage: f = x^5 + 75*x^3 - 15*x^2 + 125*x - 5 sage: W.<w> = R.ext(f) sage: w.lift_to_precision(10000) w
>>> from sage.all import * >>> R = ZpFM(Integer(5),Integer(5)) >>> S = R['x']; (x,) = S._first_ngens(1) >>> f = x**Integer(5) + Integer(75)*x**Integer(3) - Integer(15)*x**Integer(2) + Integer(125)*x - Integer(5) >>> W = R.ext(f, names=('w',)); (w,) = W._first_ngens(1) >>> w.lift_to_precision(Integer(10000)) w
- matrix_mod_pn()[source]#
Return the matrix of right multiplication by the element on the power basis \(1, x, x^2, \ldots, x^{d-1}\) for this extension field.
The rows of this matrix give the images of each of the \(x^i\). The entries of the matrices are
IntegerMod
elements, defined modulop^(self.absprec() / e)
.Raises an error if
self
has negative valuation.EXAMPLES:
sage: R = ZpFM(5,5) sage: S.<x> = R[] sage: f = x^5 + 75*x^3 - 15*x^2 + 125*x - 5 sage: W.<w> = R.ext(f) sage: a = (3+w)^7 sage: a.matrix_mod_pn() # needs sage.geometry.polyhedron [2757 333 1068 725 2510] [ 50 1507 483 318 725] [ 500 50 3007 2358 318] [1590 1375 1695 1032 2358] [2415 590 2370 2970 1032]
>>> from sage.all import * >>> R = ZpFM(Integer(5),Integer(5)) >>> S = R['x']; (x,) = S._first_ngens(1) >>> f = x**Integer(5) + Integer(75)*x**Integer(3) - Integer(15)*x**Integer(2) + Integer(125)*x - Integer(5) >>> W = R.ext(f, names=('w',)); (w,) = W._first_ngens(1) >>> a = (Integer(3)+w)**Integer(7) >>> a.matrix_mod_pn() # needs sage.geometry.polyhedron [2757 333 1068 725 2510] [ 50 1507 483 318 725] [ 500 50 3007 2358 318] [1590 1375 1695 1032 2358] [2415 590 2370 2970 1032]
- norm(base=None)[source]#
Return the absolute or relative norm of this element.
Note
This is not the \(p\)-adic absolute value. This is a field theoretic norm down to a ground ring.
If you want the \(p\)-adic absolute value, use the
abs()
function instead.If \(K\) is given then \(K\) must be a subfield of the parent \(L\) of
self
, in which case the norm is the relative norm from \(L\) to \(K\). In all other cases, the norm is the absolute norm down to \(\QQ_p\) or \(\ZZ_p\).EXAMPLES:
sage: R = ZpCR(5,5) sage: S.<x> = R[] sage: f = x^5 + 75*x^3 - 15*x^2 + 125*x - 5 sage: W.<w> = R.ext(f) sage: ((1+2*w)^5).norm() 1 + 5^2 + O(5^5) sage: ((1+2*w)).norm()^5 1 + 5^2 + O(5^5)
>>> from sage.all import * >>> R = ZpCR(Integer(5),Integer(5)) >>> S = R['x']; (x,) = S._first_ngens(1) >>> f = x**Integer(5) + Integer(75)*x**Integer(3) - Integer(15)*x**Integer(2) + Integer(125)*x - Integer(5) >>> W = R.ext(f, names=('w',)); (w,) = W._first_ngens(1) >>> ((Integer(1)+Integer(2)*w)**Integer(5)).norm() 1 + 5^2 + O(5^5) >>> ((Integer(1)+Integer(2)*w)).norm()**Integer(5) 1 + 5^2 + O(5^5)
- polynomial(var='x')[source]#
Return a polynomial over the base ring that yields this element when evaluated at the generator of the parent.
INPUT:
var
– string, the variable name for the polynomial
EXAMPLES:
sage: S.<x> = ZZ[] sage: W.<w> = ZpFM(5).extension(x^2 - 5) sage: (w + 5).polynomial() x + 5
>>> from sage.all import * >>> S = ZZ['x']; (x,) = S._first_ngens(1) >>> W = ZpFM(Integer(5)).extension(x**Integer(2) - Integer(5), names=('w',)); (w,) = W._first_ngens(1) >>> (w + Integer(5)).polynomial() x + 5
- precision_absolute()[source]#
Return the absolute precision of
self
, ie the precision cap ofself.parent()
.EXAMPLES:
sage: R = ZpFM(5,5) sage: S.<x> = R[] sage: f = x^5 + 75*x^3 - 15*x^2 + 125*x - 5 sage: W.<w> = R.ext(f) sage: a = W(75); a 3*w^10 + 2*w^12 + w^14 + w^16 + w^17 + 3*w^18 + 3*w^19 + 2*w^21 + 3*w^22 + 3*w^23 sage: a.valuation() 10 sage: a.precision_absolute() 25 sage: a.precision_relative() 15 sage: a.unit_part() 3 + 2*w^2 + w^4 + w^6 + w^7 + 3*w^8 + 3*w^9 + 2*w^11 + 3*w^12 + 3*w^13 + w^15 + 4*w^16 + 2*w^17 + w^18 + 3*w^21 + w^22 + 3*w^24
>>> from sage.all import * >>> R = ZpFM(Integer(5),Integer(5)) >>> S = R['x']; (x,) = S._first_ngens(1) >>> f = x**Integer(5) + Integer(75)*x**Integer(3) - Integer(15)*x**Integer(2) + Integer(125)*x - Integer(5) >>> W = R.ext(f, names=('w',)); (w,) = W._first_ngens(1) >>> a = W(Integer(75)); a 3*w^10 + 2*w^12 + w^14 + w^16 + w^17 + 3*w^18 + 3*w^19 + 2*w^21 + 3*w^22 + 3*w^23 >>> a.valuation() 10 >>> a.precision_absolute() 25 >>> a.precision_relative() 15 >>> a.unit_part() 3 + 2*w^2 + w^4 + w^6 + w^7 + 3*w^8 + 3*w^9 + 2*w^11 + 3*w^12 + 3*w^13 + w^15 + 4*w^16 + 2*w^17 + w^18 + 3*w^21 + w^22 + 3*w^24
- precision_relative()[source]#
Return the relative precision of
self
, ie the precision cap ofself.parent()
minus thevaluation of self
.EXAMPLES:
sage: R = ZpFM(5,5) sage: S.<x> = R[] sage: f = x^5 + 75*x^3 - 15*x^2 + 125*x - 5 sage: W.<w> = R.ext(f) sage: a = W(75); a 3*w^10 + 2*w^12 + w^14 + w^16 + w^17 + 3*w^18 + 3*w^19 + 2*w^21 + 3*w^22 + 3*w^23 sage: a.valuation() 10 sage: a.precision_absolute() 25 sage: a.precision_relative() 15 sage: a.unit_part() 3 + 2*w^2 + w^4 + w^6 + w^7 + 3*w^8 + 3*w^9 + 2*w^11 + 3*w^12 + 3*w^13 + w^15 + 4*w^16 + 2*w^17 + w^18 + 3*w^21 + w^22 + 3*w^24
>>> from sage.all import * >>> R = ZpFM(Integer(5),Integer(5)) >>> S = R['x']; (x,) = S._first_ngens(1) >>> f = x**Integer(5) + Integer(75)*x**Integer(3) - Integer(15)*x**Integer(2) + Integer(125)*x - Integer(5) >>> W = R.ext(f, names=('w',)); (w,) = W._first_ngens(1) >>> a = W(Integer(75)); a 3*w^10 + 2*w^12 + w^14 + w^16 + w^17 + 3*w^18 + 3*w^19 + 2*w^21 + 3*w^22 + 3*w^23 >>> a.valuation() 10 >>> a.precision_absolute() 25 >>> a.precision_relative() 15 >>> a.unit_part() 3 + 2*w^2 + w^4 + w^6 + w^7 + 3*w^8 + 3*w^9 + 2*w^11 + 3*w^12 + 3*w^13 + w^15 + 4*w^16 + 2*w^17 + w^18 + 3*w^21 + w^22 + 3*w^24
- teichmuller_expansion(n=None)[source]#
Return a list \([a_0, a_1, \ldots, a_n]\) such that
\(a_i^q = a_i\)
self.unit_part()
= \(\sum_{i = 0}^n a_i \pi^i\), where \(\pi\) is a uniformizer ofself.parent()
INPUT:
n
– integer (defaultNone
); f given, returns the corresponding entry in the expansion
EXAMPLES:
sage: # needs sage.libs.flint sage: R.<a> = ZqFM(5^4,4) sage: E = a.teichmuller_expansion(); E 5-adic expansion of a (teichmuller) sage: list(E) [a + (2*a^3 + 2*a^2 + 3*a + 4)*5 + (4*a^3 + 3*a^2 + 3*a + 2)*5^2 + (4*a^2 + 2*a + 2)*5^3, (3*a^3 + 3*a^2 + 2*a + 1) + (a^3 + 4*a^2 + 1)*5 + (a^2 + 4*a + 4)*5^2 + (4*a^2 + a + 3)*5^3, (4*a^3 + 2*a^2 + a + 1) + (2*a^3 + 2*a^2 + 2*a + 4)*5 + (3*a^3 + 2*a^2 + a + 1)*5^2 + (a^3 + a^2 + 2)*5^3, (a^3 + a^2 + a + 4) + (3*a^3 + 1)*5 + (3*a^3 + a + 2)*5^2 + (3*a^3 + 3*a^2 + 3*a + 1)*5^3] sage: sum([c * 5^i for i, c in enumerate(E)]) a sage: all(c^625 == c for c in E) True sage: S.<x> = ZZ[] sage: f = x^3 - 98*x + 7 sage: W.<w> = ZpFM(7,3).ext(f) sage: b = (1+w)^5; L = b.teichmuller_expansion(); L [1, 5 + 5*w^3 + w^6 + 4*w^7, 3 + 3*w^3 + w^7, 3 + 3*w^3 + w^7, 0, 4 + 5*w^3 + w^6 + 4*w^7, 3 + 3*w^3 + w^7, 6 + w^3 + 5*w^7, 6 + w^3 + 5*w^7] sage: sum([w^i*L[i] for i in range(len(L))]) == b True sage: all(L[i]^(7^3) == L[i] for i in range(9)) True sage: L = W(3).teichmuller_expansion(); L [3 + 3*w^3 + w^7, 0, 0, 4 + 5*w^3 + w^6 + 4*w^7, 0, 0, 3 + 3*w^3 + w^7, 6 + w^3 + 5*w^7] sage: sum([w^i*L[i] for i in range(len(L))]) 3
>>> from sage.all import * >>> # needs sage.libs.flint >>> R = ZqFM(Integer(5)**Integer(4),Integer(4), names=('a',)); (a,) = R._first_ngens(1) >>> E = a.teichmuller_expansion(); E 5-adic expansion of a (teichmuller) >>> list(E) [a + (2*a^3 + 2*a^2 + 3*a + 4)*5 + (4*a^3 + 3*a^2 + 3*a + 2)*5^2 + (4*a^2 + 2*a + 2)*5^3, (3*a^3 + 3*a^2 + 2*a + 1) + (a^3 + 4*a^2 + 1)*5 + (a^2 + 4*a + 4)*5^2 + (4*a^2 + a + 3)*5^3, (4*a^3 + 2*a^2 + a + 1) + (2*a^3 + 2*a^2 + 2*a + 4)*5 + (3*a^3 + 2*a^2 + a + 1)*5^2 + (a^3 + a^2 + 2)*5^3, (a^3 + a^2 + a + 4) + (3*a^3 + 1)*5 + (3*a^3 + a + 2)*5^2 + (3*a^3 + 3*a^2 + 3*a + 1)*5^3] >>> sum([c * Integer(5)**i for i, c in enumerate(E)]) a >>> all(c**Integer(625) == c for c in E) True >>> S = ZZ['x']; (x,) = S._first_ngens(1) >>> f = x**Integer(3) - Integer(98)*x + Integer(7) >>> W = ZpFM(Integer(7),Integer(3)).ext(f, names=('w',)); (w,) = W._first_ngens(1) >>> b = (Integer(1)+w)**Integer(5); L = b.teichmuller_expansion(); L [1, 5 + 5*w^3 + w^6 + 4*w^7, 3 + 3*w^3 + w^7, 3 + 3*w^3 + w^7, 0, 4 + 5*w^3 + w^6 + 4*w^7, 3 + 3*w^3 + w^7, 6 + w^3 + 5*w^7, 6 + w^3 + 5*w^7] >>> sum([w**i*L[i] for i in range(len(L))]) == b True >>> all(L[i]**(Integer(7)**Integer(3)) == L[i] for i in range(Integer(9))) True >>> L = W(Integer(3)).teichmuller_expansion(); L [3 + 3*w^3 + w^7, 0, 0, 4 + 5*w^3 + w^6 + 4*w^7, 0, 0, 3 + 3*w^3 + w^7, 6 + w^3 + 5*w^7] >>> sum([w**i*L[i] for i in range(len(L))]) 3
- trace(base=None)[source]#
Return the absolute or relative trace of this element.
If \(K\) is given then \(K\) must be a subfield of the parent \(L\) of
self
, in which case the norm is the relative norm from \(L\) to \(K\). In all other cases, the norm is the absolute norm down to \(\QQ_p\) or \(\ZZ_p\).EXAMPLES:
sage: R = ZpCR(5,5) sage: S.<x> = R[] sage: f = x^5 + 75*x^3 - 15*x^2 + 125*x - 5 sage: W.<w> = R.ext(f) sage: a = (2+3*w)^7 sage: b = (6+w^3)^5 sage: a.trace() 3*5 + 2*5^2 + 3*5^3 + 2*5^4 + O(5^5) sage: a.trace() + b.trace() 4*5 + 5^2 + 5^3 + 2*5^4 + O(5^5) sage: (a+b).trace() 4*5 + 5^2 + 5^3 + 2*5^4 + O(5^5)
>>> from sage.all import * >>> R = ZpCR(Integer(5),Integer(5)) >>> S = R['x']; (x,) = S._first_ngens(1) >>> f = x**Integer(5) + Integer(75)*x**Integer(3) - Integer(15)*x**Integer(2) + Integer(125)*x - Integer(5) >>> W = R.ext(f, names=('w',)); (w,) = W._first_ngens(1) >>> a = (Integer(2)+Integer(3)*w)**Integer(7) >>> b = (Integer(6)+w**Integer(3))**Integer(5) >>> a.trace() 3*5 + 2*5^2 + 3*5^3 + 2*5^4 + O(5^5) >>> a.trace() + b.trace() 4*5 + 5^2 + 5^3 + 2*5^4 + O(5^5) >>> (a+b).trace() 4*5 + 5^2 + 5^3 + 2*5^4 + O(5^5)
- unit_part()[source]#
Return the unit part of
self
, ieself / uniformizer^(self.valuation())
Warning
If this element has positive valuation then the unit part is not defined to the full precision of the ring. Asking for the unit part of
ZpFM(5)(0)
will not raise an error, but rather return itself.EXAMPLES:
sage: R = ZpFM(5,5) sage: S.<x> = R[] sage: f = x^5 + 75*x^3 - 15*x^2 + 125*x - 5 sage: W.<w> = R.ext(f) sage: a = W(75); a 3*w^10 + 2*w^12 + w^14 + w^16 + w^17 + 3*w^18 + 3*w^19 + 2*w^21 + 3*w^22 + 3*w^23 sage: a.valuation() 10 sage: a.precision_absolute() 25 sage: a.precision_relative() 15 sage: a.unit_part() 3 + 2*w^2 + w^4 + w^6 + w^7 + 3*w^8 + 3*w^9 + 2*w^11 + 3*w^12 + 3*w^13 + w^15 + 4*w^16 + 2*w^17 + w^18 + 3*w^21 + w^22 + 3*w^24
>>> from sage.all import * >>> R = ZpFM(Integer(5),Integer(5)) >>> S = R['x']; (x,) = S._first_ngens(1) >>> f = x**Integer(5) + Integer(75)*x**Integer(3) - Integer(15)*x**Integer(2) + Integer(125)*x - Integer(5) >>> W = R.ext(f, names=('w',)); (w,) = W._first_ngens(1) >>> a = W(Integer(75)); a 3*w^10 + 2*w^12 + w^14 + w^16 + w^17 + 3*w^18 + 3*w^19 + 2*w^21 + 3*w^22 + 3*w^23 >>> a.valuation() 10 >>> a.precision_absolute() 25 >>> a.precision_relative() 15 >>> a.unit_part() 3 + 2*w^2 + w^4 + w^6 + w^7 + 3*w^8 + 3*w^9 + 2*w^11 + 3*w^12 + 3*w^13 + w^15 + 4*w^16 + 2*w^17 + w^18 + 3*w^21 + w^22 + 3*w^24
The unit part inserts nonsense digits if this element has positive valuation:
sage: (a-a).unit_part() 0
>>> from sage.all import * >>> (a-a).unit_part() 0