Dense univariate polynomials over \(\RR\), implemented using MPFR#
- class sage.rings.polynomial.polynomial_real_mpfr_dense.PolynomialRealDense[source]#
Bases:
Polynomial
- change_ring(R)[source]#
EXAMPLES:
sage: from sage.rings.polynomial.polynomial_real_mpfr_dense import PolynomialRealDense sage: f = PolynomialRealDense(RR['x'], [-2, 0, 1.5]) sage: f.change_ring(QQ) 3/2*x^2 - 2 sage: f.change_ring(RealField(10)) 1.5*x^2 - 2.0 sage: f.change_ring(RealField(100)) 1.5000000000000000000000000000*x^2 - 2.0000000000000000000000000000
>>> from sage.all import * >>> from sage.rings.polynomial.polynomial_real_mpfr_dense import PolynomialRealDense >>> f = PolynomialRealDense(RR['x'], [-Integer(2), Integer(0), RealNumber('1.5')]) >>> f.change_ring(QQ) 3/2*x^2 - 2 >>> f.change_ring(RealField(Integer(10))) 1.5*x^2 - 2.0 >>> f.change_ring(RealField(Integer(100))) 1.5000000000000000000000000000*x^2 - 2.0000000000000000000000000000
- degree()[source]#
Return the degree of the polynomial.
EXAMPLES:
sage: from sage.rings.polynomial.polynomial_real_mpfr_dense import PolynomialRealDense sage: f = PolynomialRealDense(RR['x'], [1, 2, 3]); f 3.00000000000000*x^2 + 2.00000000000000*x + 1.00000000000000 sage: f.degree() 2
>>> from sage.all import * >>> from sage.rings.polynomial.polynomial_real_mpfr_dense import PolynomialRealDense >>> f = PolynomialRealDense(RR['x'], [Integer(1), Integer(2), Integer(3)]); f 3.00000000000000*x^2 + 2.00000000000000*x + 1.00000000000000 >>> f.degree() 2
- integral()[source]#
EXAMPLES:
sage: from sage.rings.polynomial.polynomial_real_mpfr_dense import PolynomialRealDense sage: f = PolynomialRealDense(RR['x'], [3, pi, 1]) # needs sage.symbolic sage: f.integral() # needs sage.symbolic 0.333333333333333*x^3 + 1.57079632679490*x^2 + 3.00000000000000*x
>>> from sage.all import * >>> from sage.rings.polynomial.polynomial_real_mpfr_dense import PolynomialRealDense >>> f = PolynomialRealDense(RR['x'], [Integer(3), pi, Integer(1)]) # needs sage.symbolic >>> f.integral() # needs sage.symbolic 0.333333333333333*x^3 + 1.57079632679490*x^2 + 3.00000000000000*x
- list(copy=True)[source]#
EXAMPLES:
sage: from sage.rings.polynomial.polynomial_real_mpfr_dense import PolynomialRealDense sage: f = PolynomialRealDense(RR['x'], [1, 0, -2]); f -2.00000000000000*x^2 + 1.00000000000000 sage: f.list() [1.00000000000000, 0.000000000000000, -2.00000000000000]
>>> from sage.all import * >>> from sage.rings.polynomial.polynomial_real_mpfr_dense import PolynomialRealDense >>> f = PolynomialRealDense(RR['x'], [Integer(1), Integer(0), -Integer(2)]); f -2.00000000000000*x^2 + 1.00000000000000 >>> f.list() [1.00000000000000, 0.000000000000000, -2.00000000000000]
- quo_rem(other)[source]#
Return the quotient with remainder of
self
byother
.EXAMPLES:
sage: from sage.rings.polynomial.polynomial_real_mpfr_dense import PolynomialRealDense sage: f = PolynomialRealDense(RR['x'], [-2, 0, 1]) sage: g = PolynomialRealDense(RR['x'], [5, 1]) sage: q, r = f.quo_rem(g) sage: q x - 5.00000000000000 sage: r 23.0000000000000 sage: q*g + r == f True sage: fg = f*g sage: fg.quo_rem(f) (x + 5.00000000000000, 0) sage: fg.quo_rem(g) (x^2 - 2.00000000000000, 0) sage: # needs sage.symbolic sage: f = PolynomialRealDense(RR['x'], range(5)) sage: g = PolynomialRealDense(RR['x'], [pi,3000,4]) sage: q, r = f.quo_rem(g) sage: g*q + r == f True
>>> from sage.all import * >>> from sage.rings.polynomial.polynomial_real_mpfr_dense import PolynomialRealDense >>> f = PolynomialRealDense(RR['x'], [-Integer(2), Integer(0), Integer(1)]) >>> g = PolynomialRealDense(RR['x'], [Integer(5), Integer(1)]) >>> q, r = f.quo_rem(g) >>> q x - 5.00000000000000 >>> r 23.0000000000000 >>> q*g + r == f True >>> fg = f*g >>> fg.quo_rem(f) (x + 5.00000000000000, 0) >>> fg.quo_rem(g) (x^2 - 2.00000000000000, 0) >>> # needs sage.symbolic >>> f = PolynomialRealDense(RR['x'], range(Integer(5))) >>> g = PolynomialRealDense(RR['x'], [pi,Integer(3000),Integer(4)]) >>> q, r = f.quo_rem(g) >>> g*q + r == f True
- reverse(degree=None)[source]#
Return 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: # needs sage.symbolic sage: f = RR['x']([-3, pi, 0, 1]) sage: f.reverse() -3.00000000000000*x^3 + 3.14159265358979*x^2 + 1.00000000000000 sage: f.reverse(2) -3.00000000000000*x^2 + 3.14159265358979*x sage: f.reverse(5) -3.00000000000000*x^5 + 3.14159265358979*x^4 + x^2
>>> from sage.all import * >>> # needs sage.symbolic >>> f = RR['x']([-Integer(3), pi, Integer(0), Integer(1)]) >>> f.reverse() -3.00000000000000*x^3 + 3.14159265358979*x^2 + 1.00000000000000 >>> f.reverse(Integer(2)) -3.00000000000000*x^2 + 3.14159265358979*x >>> f.reverse(Integer(5)) -3.00000000000000*x^5 + 3.14159265358979*x^4 + x^2
- shift(n)[source]#
Returns this polynomial multiplied by the power \(x^n\). If \(n\) is negative, terms below \(x^n\) will be discarded. Does not change this polynomial.
EXAMPLES:
sage: from sage.rings.polynomial.polynomial_real_mpfr_dense import PolynomialRealDense sage: f = PolynomialRealDense(RR['x'], [1, 2, 3]); f 3.00000000000000*x^2 + 2.00000000000000*x + 1.00000000000000 sage: f.shift(10) 3.00000000000000*x^12 + 2.00000000000000*x^11 + x^10 sage: f.shift(-1) 3.00000000000000*x + 2.00000000000000 sage: f.shift(-10) 0
>>> from sage.all import * >>> from sage.rings.polynomial.polynomial_real_mpfr_dense import PolynomialRealDense >>> f = PolynomialRealDense(RR['x'], [Integer(1), Integer(2), Integer(3)]); f 3.00000000000000*x^2 + 2.00000000000000*x + 1.00000000000000 >>> f.shift(Integer(10)) 3.00000000000000*x^12 + 2.00000000000000*x^11 + x^10 >>> f.shift(-Integer(1)) 3.00000000000000*x + 2.00000000000000 >>> f.shift(-Integer(10)) 0
- truncate(n)[source]#
Returns the polynomial of degree \(< n\) which is equivalent to self modulo \(x^n\).
EXAMPLES:
sage: from sage.rings.polynomial.polynomial_real_mpfr_dense import PolynomialRealDense sage: f = PolynomialRealDense(RealField(10)['x'], [1, 2, 4, 8]) sage: f.truncate(3) 4.0*x^2 + 2.0*x + 1.0 sage: f.truncate(100) 8.0*x^3 + 4.0*x^2 + 2.0*x + 1.0 sage: f.truncate(1) 1.0 sage: f.truncate(0) 0
>>> from sage.all import * >>> from sage.rings.polynomial.polynomial_real_mpfr_dense import PolynomialRealDense >>> f = PolynomialRealDense(RealField(Integer(10))['x'], [Integer(1), Integer(2), Integer(4), Integer(8)]) >>> f.truncate(Integer(3)) 4.0*x^2 + 2.0*x + 1.0 >>> f.truncate(Integer(100)) 8.0*x^3 + 4.0*x^2 + 2.0*x + 1.0 >>> f.truncate(Integer(1)) 1.0 >>> f.truncate(Integer(0)) 0
- truncate_abs(bound)[source]#
Truncate all high order coefficients below
bound
.EXAMPLES:
sage: from sage.rings.polynomial.polynomial_real_mpfr_dense import PolynomialRealDense sage: f = PolynomialRealDense(RealField(10)['x'], [10^-k for k in range(10)]) sage: f 1.0e-9*x^9 + 1.0e-8*x^8 + 1.0e-7*x^7 + 1.0e-6*x^6 + 0.000010*x^5 + 0.00010*x^4 + 0.0010*x^3 + 0.010*x^2 + 0.10*x + 1.0 sage: f.truncate_abs(0.5e-6) 1.0e-6*x^6 + 0.000010*x^5 + 0.00010*x^4 + 0.0010*x^3 + 0.010*x^2 + 0.10*x + 1.0 sage: f.truncate_abs(10.0) 0 sage: f.truncate_abs(1e-100) == f True
>>> from sage.all import * >>> from sage.rings.polynomial.polynomial_real_mpfr_dense import PolynomialRealDense >>> f = PolynomialRealDense(RealField(Integer(10))['x'], [Integer(10)**-k for k in range(Integer(10))]) >>> f 1.0e-9*x^9 + 1.0e-8*x^8 + 1.0e-7*x^7 + 1.0e-6*x^6 + 0.000010*x^5 + 0.00010*x^4 + 0.0010*x^3 + 0.010*x^2 + 0.10*x + 1.0 >>> f.truncate_abs(RealNumber('0.5e-6')) 1.0e-6*x^6 + 0.000010*x^5 + 0.00010*x^4 + 0.0010*x^3 + 0.010*x^2 + 0.10*x + 1.0 >>> f.truncate_abs(RealNumber('10.0')) 0 >>> f.truncate_abs(RealNumber('1e-100')) == f True
- sage.rings.polynomial.polynomial_real_mpfr_dense.make_PolynomialRealDense(parent, data)[source]#
EXAMPLES:
sage: from sage.rings.polynomial.polynomial_real_mpfr_dense import make_PolynomialRealDense sage: make_PolynomialRealDense(RR['x'], [1,2,3]) 3.00000000000000*x^2 + 2.00000000000000*x + 1.00000000000000
>>> from sage.all import * >>> from sage.rings.polynomial.polynomial_real_mpfr_dense import make_PolynomialRealDense >>> make_PolynomialRealDense(RR['x'], [Integer(1),Integer(2),Integer(3)]) 3.00000000000000*x^2 + 2.00000000000000*x + 1.00000000000000