Arbitrary precision floating point complex numbers using GNU MPC#
This is a binding for the arbitrary-precision floating point library GNU MPC.
We define a class MPComplexField
, each instance of which
specifies a field of floating-point complex numbers with
a specified precision shared by the real and imaginary part and a rounding
mode stating the rounding mode directions specific to real and imaginary
parts.
Individual floating-point numbers are of class MPComplexNumber
.
For floating-point representation and rounding mode description see the
documentation for the sage.rings.real_mpfr
.
AUTHORS:
Philippe Theveny (2008-10-13): initial version, adapted from
sage.rings.real_mpfr
andsage.rings.complex_mpfr
.Alex Ghitza (2008-11): cache, generators, random element, and many doctests.
Yann Laigle-Chapuy (2010-01): improves compatibility with CC, updates.
Jeroen Demeyer (2012-02): reformat documentation, make MPC a standard package.
Travis Scrimshaw (2012-10-18): Added doctests for full coverage.
Vincent Klein (2017-11-15) : add __mpc__() to class MPComplexNumber. MPComplexNumber constructor support gmpy2.mpz, gmpy2.mpq, gmpy2.mpfr and gmpy2.mpc parameters.
EXAMPLES:
sage: MPC = MPComplexField(42)
sage: a = MPC(12, '15.64E+32'); a
12.0000000000 + 1.56400000000e33*I
sage: a *a *a *a
5.98338564121e132 - 1.83633318912e101*I
sage: a + 1
13.0000000000 + 1.56400000000e33*I
sage: a / 3
4.00000000000 + 5.21333333333e32*I
sage: MPC("infinity + NaN *I")
+infinity + NaN*I
>>> from sage.all import *
>>> MPC = MPComplexField(Integer(42))
>>> a = MPC(Integer(12), '15.64E+32'); a
12.0000000000 + 1.56400000000e33*I
>>> a *a *a *a
5.98338564121e132 - 1.83633318912e101*I
>>> a + Integer(1)
13.0000000000 + 1.56400000000e33*I
>>> a / Integer(3)
4.00000000000 + 5.21333333333e32*I
>>> MPC("infinity + NaN *I")
+infinity + NaN*I
- sage.rings.complex_mpc.MPComplexField(prec=53, rnd='RNDNN', names=None)[source]#
Return the complex field with real and imaginary parts having prec bits of precision.
EXAMPLES:
sage: MPComplexField() Complex Field with 53 bits of precision sage: MPComplexField(100) Complex Field with 100 bits of precision sage: MPComplexField(100).base_ring() Real Field with 100 bits of precision sage: i = MPComplexField(200).gen() sage: i^2 -1.0000000000000000000000000000000000000000000000000000000000
>>> from sage.all import * >>> MPComplexField() Complex Field with 53 bits of precision >>> MPComplexField(Integer(100)) Complex Field with 100 bits of precision >>> MPComplexField(Integer(100)).base_ring() Real Field with 100 bits of precision >>> i = MPComplexField(Integer(200)).gen() >>> i**Integer(2) -1.0000000000000000000000000000000000000000000000000000000000
- class sage.rings.complex_mpc.MPComplexField_class[source]#
Bases:
Field
Initialize
self
.INPUT:
prec
– (integer) precision; default = 53prec is the number of bits used to represent the mantissa of both the real and imaginary part of complex floating-point number.
rnd
– (string) the rounding mode; default ='RNDNN'
Rounding mode is of the form
'RNDxy'
wherex
andy
are the rounding mode for respectively the real and imaginary parts and are one of:'N'
for rounding to nearest'Z'
for rounding towards zero'U'
for rounding towards plus infinity'D'
for rounding towards minus infinity
For example,
'RNDZU'
indicates to round the real part towards zero, and the imaginary part towards plus infinity.
EXAMPLES:
sage: MPComplexField(17) Complex Field with 17 bits of precision sage: MPComplexField() Complex Field with 53 bits of precision sage: MPComplexField(1042,'RNDDZ') Complex Field with 1042 bits of precision and rounding RNDDZ
>>> from sage.all import * >>> MPComplexField(Integer(17)) Complex Field with 17 bits of precision >>> MPComplexField() Complex Field with 53 bits of precision >>> MPComplexField(Integer(1042),'RNDDZ') Complex Field with 1042 bits of precision and rounding RNDDZ
ALGORITHMS: Computations are done using the MPC library.
- characteristic()[source]#
Return 0, since the field of complex numbers has characteristic 0.
EXAMPLES:
sage: MPComplexField(42).characteristic() 0
>>> from sage.all import * >>> MPComplexField(Integer(42)).characteristic() 0
- gen(n=0)[source]#
Return the generator of this complex field over its real subfield.
EXAMPLES:
sage: MPComplexField(34).gen() 1.00000000*I
>>> from sage.all import * >>> MPComplexField(Integer(34)).gen() 1.00000000*I
- is_exact()[source]#
Returns whether or not this field is exact, which is always
False
.EXAMPLES:
sage: MPComplexField(42).is_exact() False
>>> from sage.all import * >>> MPComplexField(Integer(42)).is_exact() False
- name()[source]#
Return the name of the complex field.
EXAMPLES:
sage: C = MPComplexField(10, 'RNDNZ'); C.name() 'MPComplexField10_RNDNZ'
>>> from sage.all import * >>> C = MPComplexField(Integer(10), 'RNDNZ'); C.name() 'MPComplexField10_RNDNZ'
- ngens()[source]#
Return 1, the number of generators of this complex field over its real subfield.
EXAMPLES:
sage: MPComplexField(34).ngens() 1
>>> from sage.all import * >>> MPComplexField(Integer(34)).ngens() 1
- prec()[source]#
Return the precision of this field of complex numbers.
EXAMPLES:
sage: MPComplexField().prec() 53 sage: MPComplexField(22).prec() 22
>>> from sage.all import * >>> MPComplexField().prec() 53 >>> MPComplexField(Integer(22)).prec() 22
- random_element(min=0, max=1)[source]#
Return a random complex number, uniformly distributed with real and imaginary parts between min and max (default 0 to 1).
EXAMPLES:
sage: MPComplexField(100).random_element(-5, 10) # random 1.9305310520925994224072377281 + 0.94745292506956219710477444855*I sage: MPComplexField(10).random_element() # random 0.12 + 0.23*I
>>> from sage.all import * >>> MPComplexField(Integer(100)).random_element(-Integer(5), Integer(10)) # random 1.9305310520925994224072377281 + 0.94745292506956219710477444855*I >>> MPComplexField(Integer(10)).random_element() # random 0.12 + 0.23*I
- rounding_mode()[source]#
Return rounding modes used for each part of a complex number.
EXAMPLES:
sage: MPComplexField().rounding_mode() 'RNDNN' sage: MPComplexField(rnd='RNDZU').rounding_mode() 'RNDZU'
>>> from sage.all import * >>> MPComplexField().rounding_mode() 'RNDNN' >>> MPComplexField(rnd='RNDZU').rounding_mode() 'RNDZU'
- class sage.rings.complex_mpc.MPComplexNumber[source]#
Bases:
FieldElement
A floating point approximation to a complex number using any specified precision common to both real and imaginary part.
- agm(right, algorithm='optimal')[source]#
Return the algebro-geometric mean of
self
andright
.EXAMPLES:
sage: MPC = MPComplexField() sage: u = MPC(1, 4) sage: v = MPC(-2,5) sage: u.agm(v, algorithm="pari") -0.410522769709397 + 4.60061063922097*I sage: u.agm(v, algorithm="principal") 1.24010691168158 - 0.472193567796433*I sage: u.agm(v, algorithm="optimal") -0.410522769709397 + 4.60061063922097*I
>>> from sage.all import * >>> MPC = MPComplexField() >>> u = MPC(Integer(1), Integer(4)) >>> v = MPC(-Integer(2),Integer(5)) >>> u.agm(v, algorithm="pari") -0.410522769709397 + 4.60061063922097*I >>> u.agm(v, algorithm="principal") 1.24010691168158 - 0.472193567796433*I >>> u.agm(v, algorithm="optimal") -0.410522769709397 + 4.60061063922097*I
- algebraic_dependency(n, **kwds)[source]#
Return an irreducible polynomial of degree at most \(n\) which is approximately satisfied by this complex number.
ALGORITHM: Uses the PARI C-library pari:algdep command.
INPUT: Type
algdep?
at the top level prompt. All additional parameters are passed onto the top-level algdep command.EXAMPLES:
sage: MPC = MPComplexField() sage: z = (1/2)*(1 + sqrt(3.0) * MPC.0); z 0.500000000000000 + 0.866025403784439*I sage: p = z.algebraic_dependency(5) sage: p x^2 - x + 1 sage: p(z) 1.11022302462516e-16
>>> from sage.all import * >>> MPC = MPComplexField() >>> z = (Integer(1)/Integer(2))*(Integer(1) + sqrt(RealNumber('3.0')) * MPC.gen(0)); z 0.500000000000000 + 0.866025403784439*I >>> p = z.algebraic_dependency(Integer(5)) >>> p x^2 - x + 1 >>> p(z) 1.11022302462516e-16
- arccos()[source]#
Return the arccosine of this complex number.
EXAMPLES:
sage: MPC = MPComplexField() sage: u = MPC(2, 4) sage: arccos(u) 1.11692611683177 - 2.19857302792094*I
>>> from sage.all import * >>> MPC = MPComplexField() >>> u = MPC(Integer(2), Integer(4)) >>> arccos(u) 1.11692611683177 - 2.19857302792094*I
- arccosh()[source]#
Return the hyperbolic arccos of this complex number.
EXAMPLES:
sage: MPC = MPComplexField() sage: u = MPC(2, 4) sage: arccosh(u) 2.19857302792094 + 1.11692611683177*I
>>> from sage.all import * >>> MPC = MPComplexField() >>> u = MPC(Integer(2), Integer(4)) >>> arccosh(u) 2.19857302792094 + 1.11692611683177*I
- arccoth()[source]#
Return the hyperbolic arccotangent of this complex number.
EXAMPLES:
sage: MPC = MPComplexField(100) sage: MPC(1,1).arccoth() 0.40235947810852509365018983331 - 0.55357435889704525150853273009*I
>>> from sage.all import * >>> MPC = MPComplexField(Integer(100)) >>> MPC(Integer(1),Integer(1)).arccoth() 0.40235947810852509365018983331 - 0.55357435889704525150853273009*I
- arccsch()[source]#
Return the hyperbolic arcsine of this complex number.
EXAMPLES:
sage: MPC = MPComplexField(100) sage: MPC(1,1).arccsch() 0.53063753095251782601650945811 - 0.45227844715119068206365839783*I
>>> from sage.all import * >>> MPC = MPComplexField(Integer(100)) >>> MPC(Integer(1),Integer(1)).arccsch() 0.53063753095251782601650945811 - 0.45227844715119068206365839783*I
- arcsech()[source]#
Return the hyperbolic arcsecant of this complex number.
EXAMPLES:
sage: MPC = MPComplexField(100) sage: MPC(1,1).arcsech() 0.53063753095251782601650945811 - 1.1185178796437059371676632938*I
>>> from sage.all import * >>> MPC = MPComplexField(Integer(100)) >>> MPC(Integer(1),Integer(1)).arcsech() 0.53063753095251782601650945811 - 1.1185178796437059371676632938*I
- arcsin()[source]#
Return the arcsine of this complex number.
EXAMPLES:
sage: MPC = MPComplexField() sage: u = MPC(2, 4) sage: arcsin(u) 0.453870209963122 + 2.19857302792094*I
>>> from sage.all import * >>> MPC = MPComplexField() >>> u = MPC(Integer(2), Integer(4)) >>> arcsin(u) 0.453870209963122 + 2.19857302792094*I
- arcsinh()[source]#
Return the hyperbolic arcsine of this complex number.
EXAMPLES:
sage: MPC = MPComplexField() sage: u = MPC(2, 4) sage: arcsinh(u) 2.18358521656456 + 1.09692154883014*I
>>> from sage.all import * >>> MPC = MPComplexField() >>> u = MPC(Integer(2), Integer(4)) >>> arcsinh(u) 2.18358521656456 + 1.09692154883014*I
- arctan()[source]#
Return the arctangent of this complex number.
EXAMPLES:
sage: MPC = MPComplexField() sage: u = MPC(-2, 4) sage: arctan(u) -1.46704821357730 + 0.200586618131234*I
>>> from sage.all import * >>> MPC = MPComplexField() >>> u = MPC(-Integer(2), Integer(4)) >>> arctan(u) -1.46704821357730 + 0.200586618131234*I
- arctanh()[source]#
Return the hyperbolic arctangent of this complex number.
EXAMPLES:
sage: MPC = MPComplexField() sage: u = MPC(2, 4) sage: arctanh(u) 0.0964156202029962 + 1.37153510396169*I
>>> from sage.all import * >>> MPC = MPComplexField() >>> u = MPC(Integer(2), Integer(4)) >>> arctanh(u) 0.0964156202029962 + 1.37153510396169*I
- argument()[source]#
The argument (angle) of the complex number, normalized so that \(-\pi < \theta \leq \pi\).
EXAMPLES:
sage: MPC = MPComplexField() sage: i = MPC.0 sage: (i^2).argument() 3.14159265358979 sage: (1+i).argument() 0.785398163397448 sage: i.argument() 1.57079632679490 sage: (-i).argument() -1.57079632679490 sage: (RR('-0.001') - i).argument() -1.57179632646156
>>> from sage.all import * >>> MPC = MPComplexField() >>> i = MPC.gen(0) >>> (i**Integer(2)).argument() 3.14159265358979 >>> (Integer(1)+i).argument() 0.785398163397448 >>> i.argument() 1.57079632679490 >>> (-i).argument() -1.57079632679490 >>> (RR('-0.001') - i).argument() -1.57179632646156
- conjugate()[source]#
Return the complex conjugate of this complex number:
\[\mathrm{conjugate}(a + ib) = a - ib.\]EXAMPLES:
sage: MPC = MPComplexField() sage: i = MPC(0, 1) sage: (1+i).conjugate() 1.00000000000000 - 1.00000000000000*I
>>> from sage.all import * >>> MPC = MPComplexField() >>> i = MPC(Integer(0), Integer(1)) >>> (Integer(1)+i).conjugate() 1.00000000000000 - 1.00000000000000*I
- cos()[source]#
Return the cosine of this complex number:
\[\cos(a + ib) = \cos a \cosh b -i \sin a \sinh b.\]EXAMPLES:
sage: MPC = MPComplexField() sage: u = MPC(2, 4) sage: cos(u) -11.3642347064011 - 24.8146514856342*I
>>> from sage.all import * >>> MPC = MPComplexField() >>> u = MPC(Integer(2), Integer(4)) >>> cos(u) -11.3642347064011 - 24.8146514856342*I
- cosh()[source]#
Return the hyperbolic cosine of this complex number:
\[\cosh(a + ib) = \cosh a \cos b + i \sinh a \sin b.\]EXAMPLES:
sage: MPC = MPComplexField() sage: u = MPC(2, 4) sage: cosh(u) -2.45913521391738 - 2.74481700679215*I
>>> from sage.all import * >>> MPC = MPComplexField() >>> u = MPC(Integer(2), Integer(4)) >>> cosh(u) -2.45913521391738 - 2.74481700679215*I
- cot()[source]#
Return the cotangent of this complex number.
EXAMPLES:
sage: MPC = MPComplexField(53) sage: (1+MPC(I)).cot() 0.217621561854403 - 0.868014142895925*I sage: i = MPComplexField(200).0 sage: (1+i).cot() 0.21762156185440268136513424360523807352075436916785404091068 - 0.86801414289592494863584920891627388827343874994609327121115*I sage: i = MPComplexField(220).0 sage: (1+i).cot() 0.21762156185440268136513424360523807352075436916785404091068124239 - 0.86801414289592494863584920891627388827343874994609327121115071646*I
>>> from sage.all import * >>> MPC = MPComplexField(Integer(53)) >>> (Integer(1)+MPC(I)).cot() 0.217621561854403 - 0.868014142895925*I >>> i = MPComplexField(Integer(200)).gen(0) >>> (Integer(1)+i).cot() 0.21762156185440268136513424360523807352075436916785404091068 - 0.86801414289592494863584920891627388827343874994609327121115*I >>> i = MPComplexField(Integer(220)).gen(0) >>> (Integer(1)+i).cot() 0.21762156185440268136513424360523807352075436916785404091068124239 - 0.86801414289592494863584920891627388827343874994609327121115071646*I
- coth()[source]#
Return the hyperbolic cotangent of this complex number.
EXAMPLES:
sage: MPC = MPComplexField(100) sage: MPC(1,1).coth() 0.86801414289592494863584920892 - 0.21762156185440268136513424361*I
>>> from sage.all import * >>> MPC = MPComplexField(Integer(100)) >>> MPC(Integer(1),Integer(1)).coth() 0.86801414289592494863584920892 - 0.21762156185440268136513424361*I
- csc()[source]#
Return the cosecant of this complex number.
EXAMPLES:
sage: MPC = MPComplexField(100) sage: MPC(1,1).csc() 0.62151801717042842123490780586 - 0.30393100162842645033448560451*I
>>> from sage.all import * >>> MPC = MPComplexField(Integer(100)) >>> MPC(Integer(1),Integer(1)).csc() 0.62151801717042842123490780586 - 0.30393100162842645033448560451*I
- csch()[source]#
Return the hyperbolic cosecant of this complex number.
EXAMPLES:
sage: MPC = MPComplexField(100) sage: MPC(1,1).csch() 0.30393100162842645033448560451 - 0.62151801717042842123490780586*I
>>> from sage.all import * >>> MPC = MPComplexField(Integer(100)) >>> MPC(Integer(1),Integer(1)).csch() 0.30393100162842645033448560451 - 0.62151801717042842123490780586*I
- dilog()[source]#
Return the complex dilogarithm of
self
.The complex dilogarithm, or Spence’s function, is defined by
\[Li_2(z) = - \int_0^z \frac{\log|1-\zeta|}{\zeta} d(\zeta) = \sum_{k=1}^\infty \frac{z^k}{k^2}.\]Note that the series definition can only be used for \(|z| < 1\).
EXAMPLES:
sage: MPC = MPComplexField() sage: a = MPC(1,0) sage: a.dilog() 1.64493406684823 sage: float(pi^2/6) # needs sage.symbolic 1.6449340668482262
>>> from sage.all import * >>> MPC = MPComplexField() >>> a = MPC(Integer(1),Integer(0)) >>> a.dilog() 1.64493406684823 >>> float(pi**Integer(2)/Integer(6)) # needs sage.symbolic 1.6449340668482262
sage: b = MPC(0,1) sage: b.dilog() -0.205616758356028 + 0.915965594177219*I
>>> from sage.all import * >>> b = MPC(Integer(0),Integer(1)) >>> b.dilog() -0.205616758356028 + 0.915965594177219*I
sage: c = MPC(0,0) sage: c.dilog() 0
>>> from sage.all import * >>> c = MPC(Integer(0),Integer(0)) >>> c.dilog() 0
- eta(omit_frac=False)[source]#
Return the value of the Dedekind \(\eta\) function on
self
, intelligently computed using \(\mathbb{SL}(2,\ZZ)\) transformations.The \(\eta\) function is
\[\eta(z) = e^{\pi i z / 12} \prod_{n=1}^{\infty}(1-e^{2\pi inz})\]INPUT:
self
– element of the upper half plane (if not, raises aValueError
).omit_frac
– (bool, default:False
), ifTrue
, omit the \(e^{\pi i z / 12}\) factor.
OUTPUT: a complex number
ALGORITHM: Uses the PARI C library.
EXAMPLES:
sage: MPC = MPComplexField() sage: i = MPC.0 sage: z = 1+i; z.eta() 0.742048775836565 + 0.198831370229911*I
>>> from sage.all import * >>> MPC = MPComplexField() >>> i = MPC.gen(0) >>> z = Integer(1)+i; z.eta() 0.742048775836565 + 0.198831370229911*I
- exp()[source]#
Return the exponential of this complex number:
\[\exp(a + ib) = \exp(a) (\cos b + i \sin b).\]EXAMPLES:
sage: MPC = MPComplexField() sage: u = MPC(2, 4) sage: exp(u) -4.82980938326939 - 5.59205609364098*I
>>> from sage.all import * >>> MPC = MPComplexField() >>> u = MPC(Integer(2), Integer(4)) >>> exp(u) -4.82980938326939 - 5.59205609364098*I
- gamma()[source]#
Return the Gamma function evaluated at this complex number.
EXAMPLES:
sage: MPC = MPComplexField(30) sage: i = MPC.0 sage: (1+i).gamma() 0.49801567 - 0.15494983*I
>>> from sage.all import * >>> MPC = MPComplexField(Integer(30)) >>> i = MPC.gen(0) >>> (Integer(1)+i).gamma() 0.49801567 - 0.15494983*I
- gamma_inc(t)[source]#
Return the incomplete Gamma function evaluated at this complex number.
EXAMPLES:
sage: C, i = MPComplexField(30).objgen() sage: (1+i).gamma_inc(2 + 3*i) # abs tol 2e-10 0.0020969149 - 0.059981914*I sage: (1+i).gamma_inc(5) -0.0013781309 + 0.0065198200*I sage: C(2).gamma_inc(1 + i) 0.70709210 - 0.42035364*I
>>> from sage.all import * >>> C, i = MPComplexField(Integer(30)).objgen() >>> (Integer(1)+i).gamma_inc(Integer(2) + Integer(3)*i) # abs tol 2e-10 0.0020969149 - 0.059981914*I >>> (Integer(1)+i).gamma_inc(Integer(5)) -0.0013781309 + 0.0065198200*I >>> C(Integer(2)).gamma_inc(Integer(1) + i) 0.70709210 - 0.42035364*I
- imag()[source]#
Return imaginary part of
self
.EXAMPLES:
sage: C = MPComplexField(100) sage: z = C(2, 3) sage: x = z.imag(); x 3.0000000000000000000000000000 sage: x.parent() Real Field with 100 bits of precision
>>> from sage.all import * >>> C = MPComplexField(Integer(100)) >>> z = C(Integer(2), Integer(3)) >>> x = z.imag(); x 3.0000000000000000000000000000 >>> x.parent() Real Field with 100 bits of precision
- is_imaginary()[source]#
Return
True
ifself
is imaginary, i.e. has real part zero.EXAMPLES:
sage: C200 = MPComplexField(200) sage: C200(1.23*i).is_imaginary() True sage: C200(1+i).is_imaginary() False
>>> from sage.all import * >>> C200 = MPComplexField(Integer(200)) >>> C200(RealNumber('1.23')*i).is_imaginary() True >>> C200(Integer(1)+i).is_imaginary() False
- is_real()[source]#
Return
True
ifself
is real, i.e. has imaginary part zero.EXAMPLES:
sage: C200 = MPComplexField(200) sage: C200(1.23).is_real() True sage: C200(1+i).is_real() False
>>> from sage.all import * >>> C200 = MPComplexField(Integer(200)) >>> C200(RealNumber('1.23')).is_real() True >>> C200(Integer(1)+i).is_real() False
- is_square()[source]#
This function always returns true as \(\CC\) is algebraically closed.
EXAMPLES:
sage: C200 = MPComplexField(200) sage: a = C200(2,1) sage: a.is_square() True
>>> from sage.all import * >>> C200 = MPComplexField(Integer(200)) >>> a = C200(Integer(2),Integer(1)) >>> a.is_square() True
\(\CC\) is algebraically closed, hence every element is a square:
sage: b = C200(5) sage: b.is_square() True
>>> from sage.all import * >>> b = C200(Integer(5)) >>> b.is_square() True
- log()[source]#
Return the logarithm of this complex number with the branch cut on the negative real axis:
\[\log(z) = \log |z| + i \arg(z).\]EXAMPLES:
sage: MPC = MPComplexField() sage: u = MPC(2, 4) sage: log(u) 1.49786613677700 + 1.10714871779409*I
>>> from sage.all import * >>> MPC = MPComplexField() >>> u = MPC(Integer(2), Integer(4)) >>> log(u) 1.49786613677700 + 1.10714871779409*I
- norm()[source]#
Return the norm of a complex number, rounded with the rounding mode of the real part. The norm is the square of the absolute value:
\[\mathrm{norm}(a + ib) = a^2 + b^2.\]OUTPUT:
A floating-point number in the real field of the real part (same precision, same rounding mode).
EXAMPLES:
This indeed acts as the square function when the imaginary component of self is equal to zero:
sage: MPC = MPComplexField() sage: a = MPC(2,1) sage: a.norm() 5.00000000000000 sage: b = MPC(4.2,0) sage: b.norm() 17.6400000000000 sage: b^2 17.6400000000000
>>> from sage.all import * >>> MPC = MPComplexField() >>> a = MPC(Integer(2),Integer(1)) >>> a.norm() 5.00000000000000 >>> b = MPC(RealNumber('4.2'),Integer(0)) >>> b.norm() 17.6400000000000 >>> b**Integer(2) 17.6400000000000
- nth_root(n, all=False)[source]#
The \(n\)-th root function.
INPUT:
all
– bool (default:False
); ifTrue
, return a list of all \(n\)-th roots.
EXAMPLES:
sage: MPC = MPComplexField() sage: a = MPC(27) sage: a.nth_root(3) 3.00000000000000 sage: a.nth_root(3, all=True) [3.00000000000000, -1.50000000000000 + 2.59807621135332*I, -1.50000000000000 - 2.59807621135332*I]
>>> from sage.all import * >>> MPC = MPComplexField() >>> a = MPC(Integer(27)) >>> a.nth_root(Integer(3)) 3.00000000000000 >>> a.nth_root(Integer(3), all=True) [3.00000000000000, -1.50000000000000 + 2.59807621135332*I, -1.50000000000000 - 2.59807621135332*I]
- prec()[source]#
Return precision of this complex number.
EXAMPLES:
sage: i = MPComplexField(2000).0 sage: i.prec() 2000
>>> from sage.all import * >>> i = MPComplexField(Integer(2000)).gen(0) >>> i.prec() 2000
- real()[source]#
Return the real part of
self
.EXAMPLES:
sage: C = MPComplexField(100) sage: z = C(2, 3) sage: x = z.real(); x 2.0000000000000000000000000000 sage: x.parent() Real Field with 100 bits of precision
>>> from sage.all import * >>> C = MPComplexField(Integer(100)) >>> z = C(Integer(2), Integer(3)) >>> x = z.real(); x 2.0000000000000000000000000000 >>> x.parent() Real Field with 100 bits of precision
- sec()[source]#
Return the secant of this complex number.
EXAMPLES:
sage: MPC = MPComplexField(100) sage: MPC(1,1).sec() 0.49833703055518678521380589177 + 0.59108384172104504805039169297*I
>>> from sage.all import * >>> MPC = MPComplexField(Integer(100)) >>> MPC(Integer(1),Integer(1)).sec() 0.49833703055518678521380589177 + 0.59108384172104504805039169297*I
- sech()[source]#
Return the hyperbolic secant of this complex number.
EXAMPLES:
sage: MPC = MPComplexField(100) sage: MPC(1,1).sech() 0.49833703055518678521380589177 - 0.59108384172104504805039169297*I
>>> from sage.all import * >>> MPC = MPComplexField(Integer(100)) >>> MPC(Integer(1),Integer(1)).sech() 0.49833703055518678521380589177 - 0.59108384172104504805039169297*I
- sin()[source]#
Return the sine of this complex number:
\[\sin(a + ib) = \sin a \cosh b + i \cos x \sinh b.\]EXAMPLES:
sage: MPC = MPComplexField() sage: u = MPC(2, 4) sage: sin(u) 24.8313058489464 - 11.3566127112182*I
>>> from sage.all import * >>> MPC = MPComplexField() >>> u = MPC(Integer(2), Integer(4)) >>> sin(u) 24.8313058489464 - 11.3566127112182*I
- sinh()[source]#
Return the hyperbolic sine of this complex number:
\[\sinh(a + ib) = \sinh a \cos b + i \cosh a \sin b.\]EXAMPLES:
sage: MPC = MPComplexField() sage: u = MPC(2, 4) sage: sinh(u) -2.37067416935200 - 2.84723908684883*I
>>> from sage.all import * >>> MPC = MPComplexField() >>> u = MPC(Integer(2), Integer(4)) >>> sinh(u) -2.37067416935200 - 2.84723908684883*I
- sqr()[source]#
Return the square of a complex number:
\[(a + ib)^2 = (a^2 - b^2) + 2iab.\]EXAMPLES:
sage: C = MPComplexField() sage: a = C(5, 1) sage: a.sqr() 24.0000000000000 + 10.0000000000000*I
>>> from sage.all import * >>> C = MPComplexField() >>> a = C(Integer(5), Integer(1)) >>> a.sqr() 24.0000000000000 + 10.0000000000000*I
- sqrt()[source]#
Return the square root, taking the branch cut to be the negative real axis:
\[\sqrt z = \sqrt{|z|}(\cos(\arg(z)/2) + i \sin(\arg(z)/2)).\]EXAMPLES:
sage: C = MPComplexField() sage: a = C(24, 10) sage: a.sqrt() 5.00000000000000 + 1.00000000000000*I
>>> from sage.all import * >>> C = MPComplexField() >>> a = C(Integer(24), Integer(10)) >>> a.sqrt() 5.00000000000000 + 1.00000000000000*I
- str(base=10, **kwds)[source]#
Return a string of
self
.INPUT:
base
– (default: 10) base for output**kwds
– other arguments to pass to thestr()
method of the real numbers in the real and imaginary parts.
EXAMPLES:
sage: MPC = MPComplexField(64) sage: z = MPC(-4, 3)/7 sage: z.str() '-0.571428571428571428564 + 0.428571428571428571436*I' sage: z.str(16) '-0.92492492492492490 + 0.6db6db6db6db6db70*I' sage: z.str(truncate=True) '-0.571428571428571429 + 0.428571428571428571*I' sage: z.str(2) '-0.1001001001001001001001001001001001001001001001001001001001001001 + 0.01101101101101101101101101101101101101101101101101101101101101110*I'
>>> from sage.all import * >>> MPC = MPComplexField(Integer(64)) >>> z = MPC(-Integer(4), Integer(3))/Integer(7) >>> z.str() '-0.571428571428571428564 + 0.428571428571428571436*I' >>> z.str(Integer(16)) '-0.92492492492492490 + 0.6db6db6db6db6db70*I' >>> z.str(truncate=True) '-0.571428571428571429 + 0.428571428571428571*I' >>> z.str(Integer(2)) '-0.1001001001001001001001001001001001001001001001001001001001001001 + 0.01101101101101101101101101101101101101101101101101101101101101110*I'
- tan()[source]#
Return the tangent of this complex number:
\[\tan(a + ib) = (\sin 2a + i \sinh 2b)/(\cos 2a + \cosh 2b).\]EXAMPLES:
sage: MPC = MPComplexField() sage: u = MPC(-2, 4) sage: tan(u) 0.000507980623470039 + 1.00043851320205*I
>>> from sage.all import * >>> MPC = MPComplexField() >>> u = MPC(-Integer(2), Integer(4)) >>> tan(u) 0.000507980623470039 + 1.00043851320205*I
- tanh()[source]#
Return the hyperbolic tangent of this complex number:
\[\tanh(a + ib) = (\sinh 2a + i \sin 2b)/(\cosh 2a + \cos 2b).\]EXAMPLES:
sage: MPC = MPComplexField() sage: u = MPC(2, 4) sage: tanh(u) 1.00468231219024 + 0.0364233692474037*I
>>> from sage.all import * >>> MPC = MPComplexField() >>> u = MPC(Integer(2), Integer(4)) >>> tanh(u) 1.00468231219024 + 0.0364233692474037*I
- zeta()[source]#
Return the Riemann zeta function evaluated at this complex number.
EXAMPLES:
sage: i = MPComplexField(30).gen() sage: z = 1 + i sage: z.zeta() 0.58215806 - 0.92684856*I
>>> from sage.all import * >>> i = MPComplexField(Integer(30)).gen() >>> z = Integer(1) + i >>> z.zeta() 0.58215806 - 0.92684856*I
- class sage.rings.complex_mpc.MPCtoMPC[source]#
Bases:
Map
- section()[source]#
EXAMPLES:
sage: from sage.rings.complex_mpc import * sage: C10 = MPComplexField(10) sage: C100 = MPComplexField(100) sage: f = MPCtoMPC(C100, C10) sage: f.section() Generic map: From: Complex Field with 10 bits of precision To: Complex Field with 100 bits of precision
>>> from sage.all import * >>> from sage.rings.complex_mpc import * >>> C10 = MPComplexField(Integer(10)) >>> C100 = MPComplexField(Integer(100)) >>> f = MPCtoMPC(C100, C10) >>> f.section() Generic map: From: Complex Field with 10 bits of precision To: Complex Field with 100 bits of precision
- sage.rings.complex_mpc.split_complex_string(string, base=10)[source]#
Split and return in that order the real and imaginary parts of a complex in a string.
This is an internal function.
EXAMPLES:
sage: sage.rings.complex_mpc.split_complex_string('123.456e789') ('123.456e789', None) sage: sage.rings.complex_mpc.split_complex_string('123.456e789*I') (None, '123.456e789') sage: sage.rings.complex_mpc.split_complex_string('123.+456e789*I') ('123.', '+456e789') sage: sage.rings.complex_mpc.split_complex_string('123.456e789', base=2) (None, None)
>>> from sage.all import * >>> sage.rings.complex_mpc.split_complex_string('123.456e789') ('123.456e789', None) >>> sage.rings.complex_mpc.split_complex_string('123.456e789*I') (None, '123.456e789') >>> sage.rings.complex_mpc.split_complex_string('123.+456e789*I') ('123.', '+456e789') >>> sage.rings.complex_mpc.split_complex_string('123.456e789', base=Integer(2)) (None, None)