p-adic Capped Relative Dense Polynomials#
- class sage.rings.polynomial.padics.polynomial_padic_capped_relative_dense.Polynomial_padic_capped_relative_dense(parent, x=None, check=True, is_gen=False, construct=False, absprec=+Infinity, relprec=+Infinity)#
Bases:
Polynomial_generic_cdv
,Polynomial_padic
- degree(secure=False)#
Return the degree of
self
.INPUT:
secure – a boolean (default:
False
)
If
secure
isTrue
and the degree of this polynomial is not determined (because the leading coefficient is indistinguishable from 0), an error is raised.If
secure
isFalse
, the returned value is the largest \(n\) so that the coefficient of \(x^n\) does not compare equal to \(0\).EXAMPLES:
sage: K = Qp(3,10) sage: R.<T> = K[] sage: f = T + 2; f (1 + O(3^10))*T + 2 + O(3^10) sage: f.degree() 1 sage: (f-T).degree() 0 sage: (f-T).degree(secure=True) Traceback (most recent call last): ... PrecisionError: the leading coefficient is indistinguishable from 0 sage: x = O(3^5) sage: li = [3^i * x for i in range(0,5)]; li [O(3^5), O(3^6), O(3^7), O(3^8), O(3^9)] sage: f = R(li); f O(3^9)*T^4 + O(3^8)*T^3 + O(3^7)*T^2 + O(3^6)*T + O(3^5) sage: f.degree() -1 sage: f.degree(secure=True) Traceback (most recent call last): ... PrecisionError: the leading coefficient is indistinguishable from 0
- disc()#
- factor_mod()#
Return the factorization of
self
modulo \(p\).
- is_eisenstein(secure=False)#
Return
True
if this polynomial is an Eisenstein polynomial.EXAMPLES:
sage: K = Qp(5) sage: R.<t> = K[] sage: f = 5 + 5*t + t^4 sage: f.is_eisenstein() True
AUTHOR:
Xavier Caruso (2013-03)
- lift()#
Return an integer polynomial congruent to this one modulo the precision of each coefficient.
Note
The lift that is returned will not necessarily be the same for polynomials with the same coefficients (i.e. same values and precisions): it will depend on how the polynomials are created.
EXAMPLES:
sage: K = Qp(13,7) sage: R.<t> = K[] sage: a = 13^7*t^3 + K(169,4)*t - 13^4 sage: a.lift() 62748517*t^3 + 169*t - 28561
- list(copy=True)#
Return a list of coefficients of
self
.Note
The length of the list returned may be greater than expected since it includes any leading zeros that have finite absolute precision.
EXAMPLES:
sage: K = Qp(13,7) sage: R.<t> = K[] sage: a = 2*t^3 + 169*t - 1 sage: a (2 + O(13^7))*t^3 + (13^2 + O(13^9))*t + 12 + 12*13 + 12*13^2 + 12*13^3 + 12*13^4 + 12*13^5 + 12*13^6 + O(13^7) sage: a.list() [12 + 12*13 + 12*13^2 + 12*13^3 + 12*13^4 + 12*13^5 + 12*13^6 + O(13^7), 13^2 + O(13^9), 0, 2 + O(13^7)]
- lshift_coeffs(shift, no_list=False)#
Return a new polynomials whose coefficients are multiplied by p^shift.
EXAMPLES:
sage: K = Qp(13, 4) sage: R.<t> = K[] sage: a = t + 52 sage: a.lshift_coeffs(3) (13^3 + O(13^7))*t + 4*13^4 + O(13^8)
- newton_polygon()#
Return the Newton polygon of this polynomial.
Note
If some coefficients have not enough precision an error is raised.
OUTPUT:
a Newton polygon
EXAMPLES:
sage: K = Qp(2, prec=5) sage: P.<x> = K[] sage: f = x^4 + 2^3*x^3 + 2^13*x^2 + 2^21*x + 2^37 sage: f.newton_polygon() Finite Newton polygon with 4 vertices: (0, 37), (1, 21), (3, 3), (4, 0) sage: K = Qp(5) sage: R.<t> = K[] sage: f = 5 + 3*t + t^4 + 25*t^10 sage: f.newton_polygon() Finite Newton polygon with 4 vertices: (0, 1), (1, 0), (4, 0), (10, 2)
Here is an example where the computation fails because precision is not sufficient:
sage: g = f + K(0,0)*t^4; g (5^2 + O(5^22))*t^10 + O(5^0)*t^4 + (3 + O(5^20))*t + 5 + O(5^21) sage: g.newton_polygon() Traceback (most recent call last): ... PrecisionError: The coefficient of t^4 has not enough precision
AUTHOR:
Xavier Caruso (2013-03-20)
- newton_slopes(repetition=True)#
Return a list of the Newton slopes of this polynomial.
These are the valuations of the roots of this polynomial.
If
repetition
isTrue
, each slope is repeated a number of times equal to its multiplicity. Otherwise it appears only one time.INPUT:
repetition
– boolean (defaultTrue
)
OUTPUT:
a list of rationals
EXAMPLES:
sage: K = Qp(5) sage: R.<t> = K[] sage: f = 5 + 3*t + t^4 + 25*t^10 sage: f.newton_polygon() Finite Newton polygon with 4 vertices: (0, 1), (1, 0), (4, 0), (10, 2) sage: f.newton_slopes() [1, 0, 0, 0, -1/3, -1/3, -1/3, -1/3, -1/3, -1/3] sage: f.newton_slopes(repetition=False) [1, 0, -1/3]
AUTHOR:
Xavier Caruso (2013-03-20)
- prec_degree()#
Return the largest \(n\) so that precision information is stored about the coefficient of \(x^n\).
Always greater than or equal to degree.
EXAMPLES:
sage: K = Qp(3,10) sage: R.<T> = K[] sage: f = T + 2; f (1 + O(3^10))*T + 2 + O(3^10) sage: f.prec_degree() 1
- precision_absolute(n=None)#
Return absolute precision information about
self
.INPUT:
self
– a p-adic polynomialn –
None
or an integer (defaultNone
).OUTPUT:
If n == None, returns a list of absolute precisions of coefficients. Otherwise, returns the absolute precision of the coefficient of x^n.
EXAMPLES:
sage: K = Qp(3,10) sage: R.<T> = K[] sage: f = T + 2; f (1 + O(3^10))*T + 2 + O(3^10) sage: f.precision_absolute() [10, 10]
- precision_relative(n=None)#
Return relative precision information about
self
.INPUT:
self
– a p-adic polynomialn –
None
or an integer (defaultNone
).OUTPUT:
If n == None, returns a list of relative precisions of coefficients. Otherwise, returns the relative precision of the coefficient of x^n.
EXAMPLES:
sage: K = Qp(3,10) sage: R.<T> = K[] sage: f = T + 2; f (1 + O(3^10))*T + 2 + O(3^10) sage: f.precision_relative() [10, 10]
- quo_rem(right, secure=False)#
Return the quotient and remainder in division of
self
byright
.EXAMPLES:
sage: K = Qp(3,10) sage: R.<T> = K[] sage: f = T + 2 sage: g = T**4 + 3*T+22 sage: g.quo_rem(f) ((1 + O(3^10))*T^3 + (1 + 2*3 + 2*3^2 + 2*3^3 + 2*3^4 + 2*3^5 + 2*3^6 + 2*3^7 + 2*3^8 + 2*3^9 + O(3^10))*T^2 + (1 + 3 + O(3^10))*T + 1 + 3 + 2*3^2 + 2*3^3 + 2*3^4 + 2*3^5 + 2*3^6 + 2*3^7 + 2*3^8 + 2*3^9 + O(3^10), 2 + 3 + 3^3 + O(3^10))
- rescale(a)#
Return f(a*X)
Todo
Need to write this function for integer polynomials before this works.
EXAMPLES:
sage: K = Zp(13, 5) sage: R.<t> = K[] sage: f = t^3 + K(13, 3) * t sage: f.rescale(2) # not implemented
- reverse(degree=None)#
Return the reverse of the input polynomial, thought as a polynomial of degree
degree
.If \(f\) is a degree-\(d\) polynomial, its reverse is \(x^d f(1/x)\).
INPUT:
degree
(None
or an integer) - if specified, truncate or zero pad the list of coefficients to this degree before reversing it.
EXAMPLES:
sage: K = Qp(13,7) sage: R.<t> = K[] sage: f = t^3 + 4*t; f (1 + O(13^7))*t^3 + (4 + O(13^7))*t sage: f.reverse() 0*t^3 + (4 + O(13^7))*t^2 + 1 + O(13^7) sage: f.reverse(3) 0*t^3 + (4 + O(13^7))*t^2 + 1 + O(13^7) sage: f.reverse(2) 0*t^2 + (4 + O(13^7))*t sage: f.reverse(4) 0*t^4 + (4 + O(13^7))*t^3 + (1 + O(13^7))*t sage: f.reverse(6) 0*t^6 + (4 + O(13^7))*t^5 + (1 + O(13^7))*t^3
- rshift_coeffs(shift, no_list=False)#
Return a new polynomial whose coefficients are p-adically shifted to the right by
shift
.Note
Type
Qp(5)(0).__rshift__?
for more information.EXAMPLES:
sage: K = Zp(13, 4) sage: R.<t> = K[] sage: a = t^2 + K(13,3)*t + 169; a (1 + O(13^4))*t^2 + (13 + O(13^3))*t + 13^2 + O(13^6) sage: b = a.rshift_coeffs(1); b O(13^3)*t^2 + (1 + O(13^2))*t + 13 + O(13^5) sage: b.list() [13 + O(13^5), 1 + O(13^2), O(13^3)] sage: b = a.rshift_coeffs(2); b O(13^2)*t^2 + O(13)*t + 1 + O(13^4) sage: b.list() [1 + O(13^4), O(13), O(13^2)]
- valuation(val_of_var=None)#
Return the valuation of
self
.INPUT:
self
– a p-adic polynomialval_of_var –
None
or a rational (defaultNone
).OUTPUT:
If val_of_var == None, returns the largest power of the variable dividing self. Otherwise, returns the valuation of
self
where the variable is assigned valuation val_of_varEXAMPLES:
sage: K = Qp(3,10) sage: R.<T> = K[] sage: f = T + 2; f (1 + O(3^10))*T + 2 + O(3^10) sage: f.valuation() 0
- valuation_of_coefficient(n=None)#
Return valuation information about
self
’s coefficients.INPUT:
self
– a p-adic polynomialn –
None
or an integer (defaultNone
).OUTPUT:
If n == None, returns a list of valuations of coefficients. Otherwise, returns the valuation of the coefficient of x^n.
EXAMPLES:
sage: K = Qp(3,10) sage: R.<T> = K[] sage: f = T + 2; f (1 + O(3^10))*T + 2 + O(3^10) sage: f.valuation_of_coefficient(1) 0
- sage.rings.polynomial.padics.polynomial_padic_capped_relative_dense.make_padic_poly(parent, x, version)#