Double precision floating point complex numbers#

Sage supports arithmetic using double-precision complex numbers. A double-precision complex number is a complex number x + I*y with \(x\), \(y\) 64-bit (8 byte) floating point numbers (double precision).

The field ComplexDoubleField implements the field of all double-precision complex numbers. You can refer to this field by the shorthand CDF. Elements of this field are of type ComplexDoubleElement. If \(x\) and \(y\) are coercible to doubles, you can create a complex double element using ComplexDoubleElement(x,y). You can coerce more general objects \(z\) to complex doubles by typing either ComplexDoubleField(x) or CDF(x).

EXAMPLES:

sage: ComplexDoubleField()
Complex Double Field
sage: CDF
Complex Double Field
sage: type(CDF.0)
<class 'sage.rings.complex_double.ComplexDoubleElement'>
sage: ComplexDoubleElement(sqrt(2), 3)                                              # needs sage.symbolic
1.4142135623730951 + 3.0*I
sage: parent(CDF(-2))
Complex Double Field
>>> from sage.all import *
>>> ComplexDoubleField()
Complex Double Field
>>> CDF
Complex Double Field
>>> type(CDF.gen(0))
<class 'sage.rings.complex_double.ComplexDoubleElement'>
>>> ComplexDoubleElement(sqrt(Integer(2)), Integer(3))                                              # needs sage.symbolic
1.4142135623730951 + 3.0*I
>>> parent(CDF(-Integer(2)))
Complex Double Field
sage: CC == CDF
False
sage: CDF is ComplexDoubleField()     # CDF is the shorthand
True
sage: CDF == ComplexDoubleField()
True
>>> from sage.all import *
>>> CC == CDF
False
>>> CDF is ComplexDoubleField()     # CDF is the shorthand
True
>>> CDF == ComplexDoubleField()
True

The underlying arithmetic of complex numbers is implemented using functions and macros in GSL (the GNU Scientific Library), and should be very fast. Also, all standard complex trig functions, log, exponents, etc., are implemented using GSL, and are also robust and fast. Several other special functions, e.g. eta, gamma, incomplete gamma, etc., are implemented using the PARI C library.

AUTHORS:

  • William Stein (2006-09): first version

  • Travis Scrimshaw (2012-10-18): Added doctests to get full coverage

  • Jeroen Demeyer (2013-02-27): fixed all PARI calls (Issue #14082)

  • Vincent Klein (2017-11-15) : add __mpc__() to class ComplexDoubleElement. ComplexDoubleElement constructor support and gmpy2.mpc parameter.

class sage.rings.complex_double.ComplexDoubleElement[source]#

Bases: FieldElement

An approximation to a complex number using double precision floating point numbers. Answers derived from calculations with such approximations may differ from what they would be if those calculations were performed with true complex numbers. This is due to the rounding errors inherent to finite precision calculations.

abs()[source]#

This function returns the magnitude \(|z|\) of the complex number \(z\).

See also

EXAMPLES:

sage: CDF(2,3).abs()
3.605551275463989
>>> from sage.all import *
>>> CDF(Integer(2),Integer(3)).abs()
3.605551275463989
abs2()[source]#

This function returns the squared magnitude \(|z|^2\) of the complex number \(z\), otherwise known as the complex norm.

See also

EXAMPLES:

sage: CDF(2,3).abs2()
13.0
>>> from sage.all import *
>>> CDF(Integer(2),Integer(3)).abs2()
13.0
agm(right, algorithm='optimal')[source]#

Return the Arithmetic-Geometric Mean (AGM) of self and right.

INPUT:

  • right (complex) – another complex number

  • algorithm (string, default "optimal") – the algorithm to use (see below).

OUTPUT:

(complex) A value of the AGM of self and right. Note that this is a multi-valued function, and the algorithm used affects the value returned, as follows:

  • 'pari': Call the pari:agm function from the pari library.

  • 'optimal': Use the AGM sequence such that at each stage \((a,b)\) is replaced by \((a_1,b_1)=((a+b)/2,\pm\sqrt{ab})\) where the sign is chosen so that \(|a_1-b_1| \leq |a_1+b_1|\), or equivalently \(\Re(b_1/a_1) \geq 0\). The resulting limit is maximal among all possible values.

  • 'principal': Use the AGM sequence such that at each stage \((a,b)\) is replaced by \((a_1,b_1)=((a+b)/2,\pm\sqrt{ab})\) where the sign is chosen so that \(\Re(b_1/a_1) \geq 0\) (the so-called principal branch of the square root).

See Wikipedia article Arithmetic-geometric mean

EXAMPLES:

sage: i = CDF(I)                                                            # needs sage.symbolic
sage: (1+i).agm(2-i)  # rel tol 1e-15                                       # needs sage.symbolic
1.6278054848727064 + 0.1368275483973686*I
>>> from sage.all import *
>>> i = CDF(I)                                                            # needs sage.symbolic
>>> (Integer(1)+i).agm(Integer(2)-i)  # rel tol 1e-15                                       # needs sage.symbolic
1.6278054848727064 + 0.1368275483973686*I

An example to show that the returned value depends on the algorithm parameter:

sage: a = CDF(-0.95,-0.65)
sage: b = CDF(0.683,0.747)
sage: a.agm(b, algorithm='optimal')
-0.3715916523517613 + 0.31989466020683*I
sage: a.agm(b, algorithm='principal')  # rel tol 1e-15
0.33817546298618006 - 0.013532696956540503*I
sage: a.agm(b, algorithm='pari')
-0.37159165235176134 + 0.31989466020683005*I
>>> from sage.all import *
>>> a = CDF(-RealNumber('0.95'),-RealNumber('0.65'))
>>> b = CDF(RealNumber('0.683'),RealNumber('0.747'))
>>> a.agm(b, algorithm='optimal')
-0.3715916523517613 + 0.31989466020683*I
>>> a.agm(b, algorithm='principal')  # rel tol 1e-15
0.33817546298618006 - 0.013532696956540503*I
>>> a.agm(b, algorithm='pari')
-0.37159165235176134 + 0.31989466020683005*I

Some degenerate cases:

sage: CDF(0).agm(a)
0.0
sage: a.agm(0)
0.0
sage: a.agm(-a)
0.0
>>> from sage.all import *
>>> CDF(Integer(0)).agm(a)
0.0
>>> a.agm(Integer(0))
0.0
>>> a.agm(-a)
0.0
algdep(n)[source]#

Returns a polynomial of degree at most \(n\) which is approximately satisfied by this complex number. Note that the returned polynomial need not be irreducible, and indeed usually won’t be if \(z\) is a good approximation to an algebraic number of degree less than \(n\).

ALGORITHM: Uses the PARI C-library algdep command.

EXAMPLES:

sage: z = (1/2)*(1 + RDF(sqrt(3)) * CDF.0); z   # abs tol 1e-16             # needs sage.symbolic
0.5 + 0.8660254037844387*I
sage: p = z.algdep(5); p                                                    # needs sage.libs.pari sage.symbolic
x^2 - x + 1
sage: abs(z^2 - z + 1) < 1e-14                                              # needs sage.symbolic
True
>>> from sage.all import *
>>> z = (Integer(1)/Integer(2))*(Integer(1) + RDF(sqrt(Integer(3))) * CDF.gen(0)); z   # abs tol 1e-16             # needs sage.symbolic
0.5 + 0.8660254037844387*I
>>> p = z.algdep(Integer(5)); p                                                    # needs sage.libs.pari sage.symbolic
x^2 - x + 1
>>> abs(z**Integer(2) - z + Integer(1)) < RealNumber('1e-14')                                              # needs sage.symbolic
True
sage: CDF(0,2).algdep(10)                                                   # needs sage.libs.pari
x^2 + 4
sage: CDF(1,5).algdep(2)                                                    # needs sage.libs.pari
x^2 - 2*x + 26
>>> from sage.all import *
>>> CDF(Integer(0),Integer(2)).algdep(Integer(10))                                                   # needs sage.libs.pari
x^2 + 4
>>> CDF(Integer(1),Integer(5)).algdep(Integer(2))                                                    # needs sage.libs.pari
x^2 - 2*x + 26
arccos()[source]#

This function returns the complex arccosine of the complex number \(z\), \({\rm arccos}(z)\). The branch cuts are on the real axis, less than -1 and greater than 1.

EXAMPLES:

sage: CDF(1,1).arccos()
0.9045568943023814 - 1.0612750619050357*I
>>> from sage.all import *
>>> CDF(Integer(1),Integer(1)).arccos()
0.9045568943023814 - 1.0612750619050357*I
arccosh()[source]#

This function returns the complex hyperbolic arccosine of the complex number \(z\), \({\rm arccosh}(z)\). The branch cut is on the real axis, less than 1.

EXAMPLES:

sage: CDF(1,1).arccosh()
1.0612750619050357 + 0.9045568943023814*I
>>> from sage.all import *
>>> CDF(Integer(1),Integer(1)).arccosh()
1.0612750619050357 + 0.9045568943023814*I
arccot()[source]#

This function returns the complex arccotangent of the complex number \(z\), \({\rm arccot}(z) = {\rm arctan}(1/z).\)

EXAMPLES:

sage: CDF(1,1).arccot()  # rel tol 1e-15
0.5535743588970452 - 0.4023594781085251*I
>>> from sage.all import *
>>> CDF(Integer(1),Integer(1)).arccot()  # rel tol 1e-15
0.5535743588970452 - 0.4023594781085251*I
arccoth()[source]#

This function returns the complex hyperbolic arccotangent of the complex number \(z\), \({\rm arccoth}(z) = {\rm arctanh(1/z)}\).

EXAMPLES:

sage: CDF(1,1).arccoth()  # rel tol 1e-15
0.4023594781085251 - 0.5535743588970452*I
>>> from sage.all import *
>>> CDF(Integer(1),Integer(1)).arccoth()  # rel tol 1e-15
0.4023594781085251 - 0.5535743588970452*I
arccsc()[source]#

This function returns the complex arccosecant of the complex number \(z\), \({\rm arccsc}(z) = {\rm arcsin}(1/z)\).

EXAMPLES:

sage: CDF(1,1).arccsc()  # rel tol 1e-15
0.45227844715119064 - 0.5306375309525178*I
>>> from sage.all import *
>>> CDF(Integer(1),Integer(1)).arccsc()  # rel tol 1e-15
0.45227844715119064 - 0.5306375309525178*I
arccsch()[source]#

This function returns the complex hyperbolic arccosecant of the complex number \(z\), \({\rm arccsch}(z) = {\rm arcsin}(1/z)\).

EXAMPLES:

sage: CDF(1,1).arccsch()  # rel tol 1e-15
0.5306375309525178 - 0.45227844715119064*I
>>> from sage.all import *
>>> CDF(Integer(1),Integer(1)).arccsch()  # rel tol 1e-15
0.5306375309525178 - 0.45227844715119064*I
arcsec()[source]#

This function returns the complex arcsecant of the complex number \(z\), \({\rm arcsec}(z) = {\rm arccos}(1/z)\).

EXAMPLES:

sage: CDF(1,1).arcsec()  # rel tol 1e-15
1.118517879643706 + 0.5306375309525178*I
>>> from sage.all import *
>>> CDF(Integer(1),Integer(1)).arcsec()  # rel tol 1e-15
1.118517879643706 + 0.5306375309525178*I
arcsech()[source]#

This function returns the complex hyperbolic arcsecant of the complex number \(z\), \({\rm arcsech}(z) = {\rm arccosh}(1/z)\).

EXAMPLES:

sage: CDF(1,1).arcsech()  # rel tol 1e-15
0.5306375309525176 - 1.118517879643706*I
>>> from sage.all import *
>>> CDF(Integer(1),Integer(1)).arcsech()  # rel tol 1e-15
0.5306375309525176 - 1.118517879643706*I
arcsin()[source]#

This function returns the complex arcsine of the complex number \(z\), \({\rm arcsin}(z)\). The branch cuts are on the real axis, less than -1 and greater than 1.

EXAMPLES:

sage: CDF(1,1).arcsin()
0.6662394324925152 + 1.0612750619050357*I
>>> from sage.all import *
>>> CDF(Integer(1),Integer(1)).arcsin()
0.6662394324925152 + 1.0612750619050357*I
arcsinh()[source]#

This function returns the complex hyperbolic arcsine of the complex number \(z\), \({\rm arcsinh}(z)\). The branch cuts are on the imaginary axis, below \(-i\) and above \(i\).

EXAMPLES:

sage: CDF(1,1).arcsinh()
1.0612750619050357 + 0.6662394324925152*I
>>> from sage.all import *
>>> CDF(Integer(1),Integer(1)).arcsinh()
1.0612750619050357 + 0.6662394324925152*I
arctan()[source]#

This function returns the complex arctangent of the complex number \(z\), \({\rm arctan}(z)\). The branch cuts are on the imaginary axis, below \(-i\) and above \(i\).

EXAMPLES:

sage: CDF(1,1).arctan()
1.0172219678978514 + 0.4023594781085251*I
>>> from sage.all import *
>>> CDF(Integer(1),Integer(1)).arctan()
1.0172219678978514 + 0.4023594781085251*I
arctanh()[source]#

This function returns the complex hyperbolic arctangent of the complex number \(z\), \({\rm arctanh} (z)\). The branch cuts are on the real axis, less than -1 and greater than 1.

EXAMPLES:

sage: CDF(1,1).arctanh()
0.4023594781085251 + 1.0172219678978514*I
>>> from sage.all import *
>>> CDF(Integer(1),Integer(1)).arctanh()
0.4023594781085251 + 1.0172219678978514*I
arg()[source]#

This function returns the argument of self, the complex number \(z\), denoted by \(\arg(z)\), where \(-\pi < \arg(z) <= \pi\).

EXAMPLES:

sage: CDF(1,0).arg()
0.0
sage: CDF(0,1).arg()
1.5707963267948966
sage: CDF(0,-1).arg()
-1.5707963267948966
sage: CDF(-1,0).arg()
3.141592653589793
>>> from sage.all import *
>>> CDF(Integer(1),Integer(0)).arg()
0.0
>>> CDF(Integer(0),Integer(1)).arg()
1.5707963267948966
>>> CDF(Integer(0),-Integer(1)).arg()
-1.5707963267948966
>>> CDF(-Integer(1),Integer(0)).arg()
3.141592653589793
argument()[source]#

This function returns the argument of the self, the complex number \(z\), in the interval \(-\pi < arg(z) \leq \pi\).

EXAMPLES:

sage: CDF(6).argument()
0.0
sage: CDF(i).argument()                                                     # needs sage.symbolic
1.5707963267948966
sage: CDF(-1).argument()
3.141592653589793
sage: CDF(-1 - 0.000001*i).argument()                                       # needs sage.symbolic
-3.1415916535897934
>>> from sage.all import *
>>> CDF(Integer(6)).argument()
0.0
>>> CDF(i).argument()                                                     # needs sage.symbolic
1.5707963267948966
>>> CDF(-Integer(1)).argument()
3.141592653589793
>>> CDF(-Integer(1) - RealNumber('0.000001')*i).argument()                                       # needs sage.symbolic
-3.1415916535897934
conj()[source]#

This function returns the complex conjugate of the complex number \(z\):

\[\overline{z} = x - i y.\]

EXAMPLES:

sage: z = CDF(2,3); z.conj()
2.0 - 3.0*I
>>> from sage.all import *
>>> z = CDF(Integer(2),Integer(3)); z.conj()
2.0 - 3.0*I
conjugate()[source]#

This function returns the complex conjugate of the complex number \(z\):

\[\overline{z} = x - i y.\]

EXAMPLES:

sage: z = CDF(2,3); z.conjugate()
2.0 - 3.0*I
>>> from sage.all import *
>>> z = CDF(Integer(2),Integer(3)); z.conjugate()
2.0 - 3.0*I
cos()[source]#

This function returns the complex cosine of the complex number \(z\):

\[\cos(z) = \frac{e^{iz} + e^{-iz}}{2}\]

EXAMPLES:

sage: CDF(1,1).cos()   # abs tol 1e-16
0.8337300251311491 - 0.9888977057628651*I
>>> from sage.all import *
>>> CDF(Integer(1),Integer(1)).cos()   # abs tol 1e-16
0.8337300251311491 - 0.9888977057628651*I
cosh()[source]#

This function returns the complex hyperbolic cosine of the complex number \(z\):

\[\cosh(z) = \frac{e^z + e^{-z}}{2}.\]

EXAMPLES:

sage: CDF(1,1).cosh()  # abs tol 1e-16
0.8337300251311491 + 0.9888977057628651*I
>>> from sage.all import *
>>> CDF(Integer(1),Integer(1)).cosh()  # abs tol 1e-16
0.8337300251311491 + 0.9888977057628651*I
cot()[source]#

This function returns the complex cotangent of the complex number \(z\):

\[\cot(z) = \frac{1}{\tan(z)}.\]

EXAMPLES:

sage: CDF(1,1).cot()  # rel tol 1e-15
0.21762156185440268 - 0.8680141428959249*I
>>> from sage.all import *
>>> CDF(Integer(1),Integer(1)).cot()  # rel tol 1e-15
0.21762156185440268 - 0.8680141428959249*I
coth()[source]#

This function returns the complex hyperbolic cotangent of the complex number \(z\):

\[\coth(z) = \frac{1}{\tanh(z)}.\]

EXAMPLES:

sage: CDF(1,1).coth()  # rel tol 1e-15
0.8680141428959249 - 0.21762156185440268*I
>>> from sage.all import *
>>> CDF(Integer(1),Integer(1)).coth()  # rel tol 1e-15
0.8680141428959249 - 0.21762156185440268*I
csc()[source]#

This function returns the complex cosecant of the complex number \(z\):

\[\csc(z) = \frac{1}{\sin(z)}.\]

EXAMPLES:

sage: CDF(1,1).csc()  # rel tol 1e-15
0.6215180171704284 - 0.30393100162842646*I
>>> from sage.all import *
>>> CDF(Integer(1),Integer(1)).csc()  # rel tol 1e-15
0.6215180171704284 - 0.30393100162842646*I
csch()[source]#

This function returns the complex hyperbolic cosecant of the complex number \(z\):

\[{\rm csch}(z) = \frac{1}{{\rm sinh}(z)}.\]

EXAMPLES:

sage: CDF(1,1).csch()  # rel tol 1e-15
0.30393100162842646 - 0.6215180171704284*I
>>> from sage.all import *
>>> CDF(Integer(1),Integer(1)).csch()  # rel tol 1e-15
0.30393100162842646 - 0.6215180171704284*I
dilog()[source]#

Returns the principal branch of the dilogarithm of \(x\), i.e., analytic continuation of the power series

\[\log_2(x) = \sum_{n \ge 1} x^n / n^2.\]

EXAMPLES:

sage: CDF(1,2).dilog()                                                      # needs sage.libs.pari
-0.059474798673809476 + 2.0726479717747566*I
sage: CDF(10000000,10000000).dilog()                                        # needs sage.libs.pari
-134.411774490731 + 38.79396299904504*I
>>> from sage.all import *
>>> CDF(Integer(1),Integer(2)).dilog()                                                      # needs sage.libs.pari
-0.059474798673809476 + 2.0726479717747566*I
>>> CDF(Integer(10000000),Integer(10000000)).dilog()                                        # needs sage.libs.pari
-134.411774490731 + 38.79396299904504*I
eta(omit_frac=0)[source]#

Return the value of the Dedekind \(\eta\) function on self.

INPUT:

  • self – element of the upper half plane (if not, raises a ValueError).

  • omit_frac – (bool, default: False), if True, omit the \(e^{\pi i z / 12}\) factor.

OUTPUT: a complex double number

ALGORITHM: Uses the PARI C library.

The \(\eta\) function is

\[\eta(z) = e^{\pi i z / 12} \prod_{n=1}^{\infty} (1 - e^{2\pi inz})\]

EXAMPLES:

We compute a few values of eta():

sage: CDF(0,1).eta()                                                        # needs sage.libs.pari
0.7682254223260566
sage: CDF(1,1).eta()                                                        # needs sage.libs.pari
0.7420487758365647 + 0.1988313702299107*I
sage: CDF(25,1).eta()                                                       # needs sage.libs.pari
0.7420487758365647 + 0.1988313702299107*I
>>> from sage.all import *
>>> CDF(Integer(0),Integer(1)).eta()                                                        # needs sage.libs.pari
0.7682254223260566
>>> CDF(Integer(1),Integer(1)).eta()                                                        # needs sage.libs.pari
0.7420487758365647 + 0.1988313702299107*I
>>> CDF(Integer(25),Integer(1)).eta()                                                       # needs sage.libs.pari
0.7420487758365647 + 0.1988313702299107*I

eta() works even if the inputs are large:

sage: CDF(0, 10^15).eta()
0.0
sage: CDF(10^15, 0.1).eta()  # abs tol 1e-10                                # needs sage.libs.pari
-0.115342592727 - 0.19977923088*I
>>> from sage.all import *
>>> CDF(Integer(0), Integer(10)**Integer(15)).eta()
0.0
>>> CDF(Integer(10)**Integer(15), RealNumber('0.1')).eta()  # abs tol 1e-10                                # needs sage.libs.pari
-0.115342592727 - 0.19977923088*I

We compute a few values of eta(), but with the fractional power of \(e\) omitted:

sage: CDF(0,1).eta(True)                                                    # needs sage.libs.pari
0.9981290699259585
>>> from sage.all import *
>>> CDF(Integer(0),Integer(1)).eta(True)                                                    # needs sage.libs.pari
0.9981290699259585

We compute eta() to low precision directly from the definition:

sage: z = CDF(1,1); z.eta()                                                 # needs sage.libs.pari
0.7420487758365647 + 0.1988313702299107*I
sage: i = CDF(0,1); pi = CDF(pi)                                            # needs sage.symbolic
sage: exp(pi * i * z / 12) * prod(1 - exp(2*pi*i*n*z)                       # needs sage.libs.pari sage.symbolic
....:                             for n in range(1, 10))
0.7420487758365647 + 0.19883137022991068*I
>>> from sage.all import *
>>> z = CDF(Integer(1),Integer(1)); z.eta()                                                 # needs sage.libs.pari
0.7420487758365647 + 0.1988313702299107*I
>>> i = CDF(Integer(0),Integer(1)); pi = CDF(pi)                                            # needs sage.symbolic
>>> exp(pi * i * z / Integer(12)) * prod(Integer(1) - exp(Integer(2)*pi*i*n*z)                       # needs sage.libs.pari sage.symbolic
...                             for n in range(Integer(1), Integer(10)))
0.7420487758365647 + 0.19883137022991068*I

The optional argument allows us to omit the fractional part:

sage: z.eta(omit_frac=True)                                                 # needs sage.libs.pari
0.9981290699259585
sage: pi = CDF(pi)                                                          # needs sage.symbolic
sage: prod(1 - exp(2*pi*i*n*z) for n in range(1,10))  # abs tol 1e-12       # needs sage.libs.pari sage.symbolic
0.998129069926 + 4.59084695545e-19*I
>>> from sage.all import *
>>> z.eta(omit_frac=True)                                                 # needs sage.libs.pari
0.9981290699259585
>>> pi = CDF(pi)                                                          # needs sage.symbolic
>>> prod(Integer(1) - exp(Integer(2)*pi*i*n*z) for n in range(Integer(1),Integer(10)))  # abs tol 1e-12       # needs sage.libs.pari sage.symbolic
0.998129069926 + 4.59084695545e-19*I

We illustrate what happens when \(z\) is not in the upper half plane:

sage: z = CDF(1)
sage: z.eta()
Traceback (most recent call last):
...
ValueError: value must be in the upper half plane
>>> from sage.all import *
>>> z = CDF(Integer(1))
>>> z.eta()
Traceback (most recent call last):
...
ValueError: value must be in the upper half plane

You can also use functional notation:

sage: z = CDF(1,1)
sage: eta(z)                                                                # needs sage.libs.pari
0.7420487758365647 + 0.1988313702299107*I
>>> from sage.all import *
>>> z = CDF(Integer(1),Integer(1))
>>> eta(z)                                                                # needs sage.libs.pari
0.7420487758365647 + 0.1988313702299107*I
exp()[source]#

This function returns the complex exponential of the complex number \(z\), \(\exp(z)\).

EXAMPLES:

sage: CDF(1,1).exp()  # abs tol 4e-16
1.4686939399158851 + 2.2873552871788423*I
>>> from sage.all import *
>>> CDF(Integer(1),Integer(1)).exp()  # abs tol 4e-16
1.4686939399158851 + 2.2873552871788423*I

We numerically verify a famous identity to the precision of a double:

sage: z = CDF(0, 2*pi); z                                                   # needs sage.symbolic
6.283185307179586*I
sage: exp(z)  # rel tol 1e-4                                                # needs sage.symbolic
1.0 - 2.4492935982947064e-16*I
>>> from sage.all import *
>>> z = CDF(Integer(0), Integer(2)*pi); z                                                   # needs sage.symbolic
6.283185307179586*I
>>> exp(z)  # rel tol 1e-4                                                # needs sage.symbolic
1.0 - 2.4492935982947064e-16*I
gamma()[source]#

Return the gamma function \(\Gamma(z)\) evaluated at self, the complex number \(z\).

EXAMPLES:

sage: # needs sage.libs.pari
sage: CDF(5,0).gamma()
24.0
sage: CDF(1,1).gamma()
0.49801566811835607 - 0.15494982830181067*I
sage: CDF(0).gamma()
Infinity
sage: CDF(-1,0).gamma()
Infinity
>>> from sage.all import *
>>> # needs sage.libs.pari
>>> CDF(Integer(5),Integer(0)).gamma()
24.0
>>> CDF(Integer(1),Integer(1)).gamma()
0.49801566811835607 - 0.15494982830181067*I
>>> CDF(Integer(0)).gamma()
Infinity
>>> CDF(-Integer(1),Integer(0)).gamma()
Infinity
gamma_inc(t)[source]#

Return the incomplete gamma function evaluated at this complex number.

EXAMPLES:

sage: CDF(1,1).gamma_inc(CDF(2,3))                                          # needs sage.libs.pari
0.0020969148636468277 - 0.059981913655449706*I
sage: CDF(1,1).gamma_inc(5)                                                 # needs sage.libs.pari
-0.001378130936215849 + 0.006519820023119819*I
sage: CDF(2,0).gamma_inc(CDF(1,1))                                          # needs sage.libs.pari
0.7070920963459381 - 0.4203536409598115*I
>>> from sage.all import *
>>> CDF(Integer(1),Integer(1)).gamma_inc(CDF(Integer(2),Integer(3)))                                          # needs sage.libs.pari
0.0020969148636468277 - 0.059981913655449706*I
>>> CDF(Integer(1),Integer(1)).gamma_inc(Integer(5))                                                 # needs sage.libs.pari
-0.001378130936215849 + 0.006519820023119819*I
>>> CDF(Integer(2),Integer(0)).gamma_inc(CDF(Integer(1),Integer(1)))                                          # needs sage.libs.pari
0.7070920963459381 - 0.4203536409598115*I
imag()[source]#

Return the imaginary part of this complex double.

EXAMPLES:

sage: a = CDF(3,-2)
sage: a.imag()
-2.0
sage: a.imag_part()
-2.0
>>> from sage.all import *
>>> a = CDF(Integer(3),-Integer(2))
>>> a.imag()
-2.0
>>> a.imag_part()
-2.0
imag_part()[source]#

Return the imaginary part of this complex double.

EXAMPLES:

sage: a = CDF(3,-2)
sage: a.imag()
-2.0
sage: a.imag_part()
-2.0
>>> from sage.all import *
>>> a = CDF(Integer(3),-Integer(2))
>>> a.imag()
-2.0
>>> a.imag_part()
-2.0
is_NaN()[source]#

Check if self is not-a-number.

EXAMPLES:

sage: CDF(1, 2).is_NaN()
False
sage: CDF(NaN).is_NaN()                                                     # needs sage.symbolic
True
sage: (1/CDF(0, 0)).is_NaN()
True
>>> from sage.all import *
>>> CDF(Integer(1), Integer(2)).is_NaN()
False
>>> CDF(NaN).is_NaN()                                                     # needs sage.symbolic
True
>>> (Integer(1)/CDF(Integer(0), Integer(0))).is_NaN()
True
is_infinity()[source]#

Check if self is \(\infty\).

EXAMPLES:

sage: CDF(1, 2).is_infinity()
False
sage: CDF(0, oo).is_infinity()
True
>>> from sage.all import *
>>> CDF(Integer(1), Integer(2)).is_infinity()
False
>>> CDF(Integer(0), oo).is_infinity()
True
is_integer()[source]#

Return True if this number is a integer

EXAMPLES:

sage: CDF(0.5).is_integer()
False
sage: CDF(I).is_integer()                                                   # needs sage.symbolic
False
sage: CDF(2).is_integer()
True
>>> from sage.all import *
>>> CDF(RealNumber('0.5')).is_integer()
False
>>> CDF(I).is_integer()                                                   # needs sage.symbolic
False
>>> CDF(Integer(2)).is_integer()
True
is_negative_infinity()[source]#

Check if self is \(-\infty\).

EXAMPLES:

sage: CDF(1, 2).is_negative_infinity()
False
sage: CDF(-oo, 0).is_negative_infinity()
True
sage: CDF(0, -oo).is_negative_infinity()
False
>>> from sage.all import *
>>> CDF(Integer(1), Integer(2)).is_negative_infinity()
False
>>> CDF(-oo, Integer(0)).is_negative_infinity()
True
>>> CDF(Integer(0), -oo).is_negative_infinity()
False
is_positive_infinity()[source]#

Check if self is \(+\infty\).

EXAMPLES:

sage: CDF(1, 2).is_positive_infinity()
False
sage: CDF(oo, 0).is_positive_infinity()
True
sage: CDF(0, oo).is_positive_infinity()
False
>>> from sage.all import *
>>> CDF(Integer(1), Integer(2)).is_positive_infinity()
False
>>> CDF(oo, Integer(0)).is_positive_infinity()
True
>>> CDF(Integer(0), oo).is_positive_infinity()
False
is_square()[source]#

This function always returns True as \(\CC\) is algebraically closed.

EXAMPLES:

sage: CDF(-1).is_square()
True
>>> from sage.all import *
>>> CDF(-Integer(1)).is_square()
True
log(base=None)[source]#

This function returns the complex natural logarithm to the given base of the complex number \(z\), \(\log(z)\). The branch cut is the negative real axis.

INPUT:

  • base – default: \(e\), the base of the natural logarithm

EXAMPLES:

sage: CDF(1,1).log()
0.34657359027997264 + 0.7853981633974483*I
>>> from sage.all import *
>>> CDF(Integer(1),Integer(1)).log()
0.34657359027997264 + 0.7853981633974483*I

This is the only example different from the GSL:

sage: CDF(0,0).log()
-infinity
>>> from sage.all import *
>>> CDF(Integer(0),Integer(0)).log()
-infinity
log10()[source]#

This function returns the complex base-10 logarithm of the complex number \(z\), \(\log_{10}(z)\).

The branch cut is the negative real axis.

EXAMPLES:

sage: CDF(1,1).log10()
0.15051499783199057 + 0.3410940884604603*I
>>> from sage.all import *
>>> CDF(Integer(1),Integer(1)).log10()
0.15051499783199057 + 0.3410940884604603*I
log_b(b)[source]#

This function returns the complex base-\(b\) logarithm of the complex number \(z\), \(\log_b(z)\). This quantity is computed as the ratio \(\log(z)/\log(b)\).

The branch cut is the negative real axis.

EXAMPLES:

sage: CDF(1,1).log_b(10)  # rel tol 1e-15
0.15051499783199057 + 0.3410940884604603*I
>>> from sage.all import *
>>> CDF(Integer(1),Integer(1)).log_b(Integer(10))  # rel tol 1e-15
0.15051499783199057 + 0.3410940884604603*I
logabs()[source]#

This function returns the natural logarithm of the magnitude of the complex number \(z\), \(\log|z|\).

This allows for an accurate evaluation of \(\log|z|\) when \(|z|\) is close to \(1\). The direct evaluation of log(abs(z)) would lead to a loss of precision in this case.

EXAMPLES:

sage: CDF(1.1,0.1).logabs()
0.09942542937258267
sage: log(abs(CDF(1.1,0.1)))
0.09942542937258259
>>> from sage.all import *
>>> CDF(RealNumber('1.1'),RealNumber('0.1')).logabs()
0.09942542937258267
>>> log(abs(CDF(RealNumber('1.1'),RealNumber('0.1'))))
0.09942542937258259
sage: log(abs(ComplexField(200)(1.1,0.1)))
0.099425429372582595066319157757531449594489450091985182495705
>>> from sage.all import *
>>> log(abs(ComplexField(Integer(200))(RealNumber('1.1'),RealNumber('0.1'))))
0.099425429372582595066319157757531449594489450091985182495705
norm()[source]#

This function returns the squared magnitude \(|z|^2\) of the complex number \(z\), otherwise known as the complex norm. If \(c = a + bi\) is a complex number, then the norm of \(c\) is defined as the product of \(c\) and its complex conjugate:

\[\text{norm}(c) = \text{norm}(a + bi) = c \cdot \overline{c} = a^2 + b^2.\]

The norm of a complex number is different from its absolute value. The absolute value of a complex number is defined to be the square root of its norm. A typical use of the complex norm is in the integral domain \(\ZZ[i]\) of Gaussian integers, where the norm of each Gaussian integer \(c = a + bi\) is defined as its complex norm.

EXAMPLES:

sage: CDF(2,3).norm()
13.0
>>> from sage.all import *
>>> CDF(Integer(2),Integer(3)).norm()
13.0
nth_root(n, all=False)[source]#

The n-th root function.

INPUT:

  • all – bool (default: False); if True, return a list of all n-th roots.

EXAMPLES:

sage: a = CDF(125)
sage: a.nth_root(3)
5.000000000000001
sage: a = CDF(10, 2)
sage: [r^5 for r in a.nth_root(5, all=True)]  # rel tol 1e-14
[9.999999999999998 + 2.0*I, 9.999999999999993 + 2.000000000000002*I, 9.999999999999996 + 1.9999999999999907*I, 9.999999999999993 + 2.0000000000000004*I, 9.999999999999998 + 1.9999999999999802*I]
sage: abs(sum(a.nth_root(111, all=True)))  # rel tol 0.1
1.1057313523818259e-13
>>> from sage.all import *
>>> a = CDF(Integer(125))
>>> a.nth_root(Integer(3))
5.000000000000001
>>> a = CDF(Integer(10), Integer(2))
>>> [r**Integer(5) for r in a.nth_root(Integer(5), all=True)]  # rel tol 1e-14
[9.999999999999998 + 2.0*I, 9.999999999999993 + 2.000000000000002*I, 9.999999999999996 + 1.9999999999999907*I, 9.999999999999993 + 2.0000000000000004*I, 9.999999999999998 + 1.9999999999999802*I]
>>> abs(sum(a.nth_root(Integer(111), all=True)))  # rel tol 0.1
1.1057313523818259e-13
prec()[source]#

Returns the precision of this number (to be more similar to ComplexNumber). Always returns 53.

EXAMPLES:

sage: CDF(0).prec()
53
>>> from sage.all import *
>>> CDF(Integer(0)).prec()
53
real()[source]#

Return the real part of this complex double.

EXAMPLES:

sage: a = CDF(3,-2)
sage: a.real()
3.0
sage: a.real_part()
3.0
>>> from sage.all import *
>>> a = CDF(Integer(3),-Integer(2))
>>> a.real()
3.0
>>> a.real_part()
3.0
real_part()[source]#

Return the real part of this complex double.

EXAMPLES:

sage: a = CDF(3,-2)
sage: a.real()
3.0
sage: a.real_part()
3.0
>>> from sage.all import *
>>> a = CDF(Integer(3),-Integer(2))
>>> a.real()
3.0
>>> a.real_part()
3.0
sec()[source]#

This function returns the complex secant of the complex number \(z\):

\[{\rm sec}(z) = \frac{1}{\cos(z)}.\]

EXAMPLES:

sage: CDF(1,1).sec()  # rel tol 1e-15
0.4983370305551868 + 0.591083841721045*I
>>> from sage.all import *
>>> CDF(Integer(1),Integer(1)).sec()  # rel tol 1e-15
0.4983370305551868 + 0.591083841721045*I
sech()[source]#

This function returns the complex hyperbolic secant of the complex number \(z\):

\[{\rm sech}(z) = \frac{1}{{\rm cosh}(z)}.\]

EXAMPLES:

sage: CDF(1,1).sech()  # rel tol 1e-15
0.4983370305551868 - 0.591083841721045*I
>>> from sage.all import *
>>> CDF(Integer(1),Integer(1)).sech()  # rel tol 1e-15
0.4983370305551868 - 0.591083841721045*I
sin()[source]#

This function returns the complex sine of the complex number \(z\):

\[\sin(z) = \frac{e^{iz} - e^{-iz}}{2i}.\]

EXAMPLES:

sage: CDF(1,1).sin()
1.2984575814159773 + 0.6349639147847361*I
>>> from sage.all import *
>>> CDF(Integer(1),Integer(1)).sin()
1.2984575814159773 + 0.6349639147847361*I
sinh()[source]#

This function returns the complex hyperbolic sine of the complex number \(z\):

\[\sinh(z) = \frac{e^z - e^{-z}}{2}.\]

EXAMPLES:

sage: CDF(1,1).sinh()
0.6349639147847361 + 1.2984575814159773*I
>>> from sage.all import *
>>> CDF(Integer(1),Integer(1)).sinh()
0.6349639147847361 + 1.2984575814159773*I
sqrt(all=False, **kwds)[source]#

The square root function.

INPUT:

  • all – bool (default: False); if True, return a list of all square roots.

If all is False, the branch cut is the negative real axis. The result always lies in the right half of the complex plane.

EXAMPLES:

We compute several square roots:

sage: a = CDF(2,3)
sage: b = a.sqrt(); b  # rel tol 1e-15
1.6741492280355401 + 0.8959774761298381*I
sage: b^2  # rel tol 1e-15
2.0 + 3.0*I
sage: a^(1/2)   # abs tol 1e-16
1.6741492280355401 + 0.895977476129838*I
>>> from sage.all import *
>>> a = CDF(Integer(2),Integer(3))
>>> b = a.sqrt(); b  # rel tol 1e-15
1.6741492280355401 + 0.8959774761298381*I
>>> b**Integer(2)  # rel tol 1e-15
2.0 + 3.0*I
>>> a**(Integer(1)/Integer(2))   # abs tol 1e-16
1.6741492280355401 + 0.895977476129838*I

We compute the square root of -1:

sage: a = CDF(-1)
sage: a.sqrt()
1.0*I
>>> from sage.all import *
>>> a = CDF(-Integer(1))
>>> a.sqrt()
1.0*I

We compute all square roots:

sage: CDF(-2).sqrt(all=True)
[1.4142135623730951*I, -1.4142135623730951*I]
sage: CDF(0).sqrt(all=True)
[0.0]
>>> from sage.all import *
>>> CDF(-Integer(2)).sqrt(all=True)
[1.4142135623730951*I, -1.4142135623730951*I]
>>> CDF(Integer(0)).sqrt(all=True)
[0.0]
tan()[source]#

This function returns the complex tangent of the complex number \(z\):

\[\tan(z) = \frac{\sin(z)}{\cos(z)}.\]

EXAMPLES:

sage: CDF(1,1).tan()
0.27175258531951174 + 1.0839233273386946*I
>>> from sage.all import *
>>> CDF(Integer(1),Integer(1)).tan()
0.27175258531951174 + 1.0839233273386946*I
tanh()[source]#

This function returns the complex hyperbolic tangent of the complex number \(z\):

\[\tanh(z) = \frac{\sinh(z)}{\cosh(z)}.\]

EXAMPLES:

sage: CDF(1,1).tanh()
1.0839233273386946 + 0.27175258531951174*I
>>> from sage.all import *
>>> CDF(Integer(1),Integer(1)).tanh()
1.0839233273386946 + 0.27175258531951174*I
zeta()[source]#

Return the Riemann zeta function evaluated at this complex number.

EXAMPLES:

sage: z = CDF(1, 1)
sage: z.zeta()                                                              # needs sage.libs.pari
0.5821580597520036 - 0.9268485643308071*I
sage: zeta(z)                                                               # needs sage.libs.pari
0.5821580597520036 - 0.9268485643308071*I
sage: zeta(CDF(1))                                                          # needs sage.libs.pari
Infinity
>>> from sage.all import *
>>> z = CDF(Integer(1), Integer(1))
>>> z.zeta()                                                              # needs sage.libs.pari
0.5821580597520036 - 0.9268485643308071*I
>>> zeta(z)                                                               # needs sage.libs.pari
0.5821580597520036 - 0.9268485643308071*I
>>> zeta(CDF(Integer(1)))                                                          # needs sage.libs.pari
Infinity
sage.rings.complex_double.ComplexDoubleField()[source]#

Returns the field of double precision complex numbers.

EXAMPLES:

sage: ComplexDoubleField()
Complex Double Field
sage: ComplexDoubleField() is CDF
True
>>> from sage.all import *
>>> ComplexDoubleField()
Complex Double Field
>>> ComplexDoubleField() is CDF
True
class sage.rings.complex_double.ComplexDoubleField_class[source]#

Bases: ComplexDoubleField

An approximation to the field of complex numbers using double precision floating point numbers. Answers derived from calculations in this approximation may differ from what they would be if those calculations were performed in the true field of complex numbers. This is due to the rounding errors inherent to finite precision calculations.

ALGORITHM:

Arithmetic is done using GSL (the GNU Scientific Library).

algebraic_closure()[source]#

Returns the algebraic closure of self, i.e., the complex double field.

EXAMPLES:

sage: CDF.algebraic_closure()
Complex Double Field
>>> from sage.all import *
>>> CDF.algebraic_closure()
Complex Double Field
characteristic()[source]#

Return the characteristic of the complex double field, which is 0.

EXAMPLES:

sage: CDF.characteristic()
0
>>> from sage.all import *
>>> CDF.characteristic()
0
construction()[source]#

Returns the functorial construction of self, namely, algebraic closure of the real double field.

EXAMPLES:

sage: c, S = CDF.construction(); S
Real Double Field
sage: CDF == c(S)
True
>>> from sage.all import *
>>> c, S = CDF.construction(); S
Real Double Field
>>> CDF == c(S)
True
gen(n=0)[source]#

Return the generator of the complex double field.

EXAMPLES:

sage: CDF.0
1.0*I
sage: CDF.gen(0)
1.0*I
>>> from sage.all import *
>>> CDF.gen(0)
1.0*I
>>> CDF.gen(Integer(0))
1.0*I
is_exact()[source]#

Returns whether or not this field is exact, which is always False.

EXAMPLES:

sage: CDF.is_exact()
False
>>> from sage.all import *
>>> CDF.is_exact()
False
ngens()[source]#

The number of generators of this complex field as an \(\RR\)-algebra.

There is one generator, namely sqrt(-1).

EXAMPLES:

sage: CDF.ngens()
1
>>> from sage.all import *
>>> CDF.ngens()
1
pi()[source]#

Returns \(\pi\) as a double precision complex number.

EXAMPLES:

sage: CDF.pi()
3.141592653589793
>>> from sage.all import *
>>> CDF.pi()
3.141592653589793
prec()[source]#

Return the precision of this complex double field (to be more similar to ComplexField). Always returns 53.

EXAMPLES:

sage: CDF.prec()
53
>>> from sage.all import *
>>> CDF.prec()
53
precision()[source]#

Return the precision of this complex double field (to be more similar to ComplexField). Always returns 53.

EXAMPLES:

sage: CDF.prec()
53
>>> from sage.all import *
>>> CDF.prec()
53
random_element(xmin=-1, xmax=1, ymin=-1, ymax=1)[source]#

Return a random element of this complex double field with real and imaginary part bounded by xmin, xmax, ymin, ymax.

EXAMPLES:

sage: CDF.random_element().parent() is CDF
True
sage: re, im = CDF.random_element()
sage: -1 <= re <= 1, -1 <= im <= 1
(True, True)
sage: re, im = CDF.random_element(-10,10,-10,10)
sage: -10 <= re <= 10, -10 <= im <= 10
(True, True)
sage: re, im = CDF.random_element(-10^20,10^20,-2,2)
sage: -10^20 <= re <= 10^20, -2 <= im <= 2
(True, True)
>>> from sage.all import *
>>> CDF.random_element().parent() is CDF
True
>>> re, im = CDF.random_element()
>>> -Integer(1) <= re <= Integer(1), -Integer(1) <= im <= Integer(1)
(True, True)
>>> re, im = CDF.random_element(-Integer(10),Integer(10),-Integer(10),Integer(10))
>>> -Integer(10) <= re <= Integer(10), -Integer(10) <= im <= Integer(10)
(True, True)
>>> re, im = CDF.random_element(-Integer(10)**Integer(20),Integer(10)**Integer(20),-Integer(2),Integer(2))
>>> -Integer(10)**Integer(20) <= re <= Integer(10)**Integer(20), -Integer(2) <= im <= Integer(2)
(True, True)
real_double_field()[source]#

The real double field, which you may view as a subfield of this complex double field.

EXAMPLES:

sage: CDF.real_double_field()
Real Double Field
>>> from sage.all import *
>>> CDF.real_double_field()
Real Double Field
to_prec(prec)[source]#

Returns the complex field to the specified precision. As doubles have fixed precision, this will only return a complex double field if prec is exactly 53.

EXAMPLES:

sage: CDF.to_prec(53)
Complex Double Field
sage: CDF.to_prec(250)
Complex Field with 250 bits of precision
>>> from sage.all import *
>>> CDF.to_prec(Integer(53))
Complex Double Field
>>> CDF.to_prec(Integer(250))
Complex Field with 250 bits of precision
zeta(n=2)[source]#

Return a primitive \(n\)-th root of unity in this CDF, for \(n \geq 1\).

INPUT:

  • n – a positive integer (default: 2)

OUTPUT: a complex \(n\)-th root of unity.

EXAMPLES:

sage: CDF.zeta(7)  # rel tol 1e-15
0.6234898018587336 + 0.7818314824680298*I
sage: CDF.zeta(1)
1.0
sage: CDF.zeta()
-1.0
sage: CDF.zeta() == CDF.zeta(2)
True
>>> from sage.all import *
>>> CDF.zeta(Integer(7))  # rel tol 1e-15
0.6234898018587336 + 0.7818314824680298*I
>>> CDF.zeta(Integer(1))
1.0
>>> CDF.zeta()
-1.0
>>> CDF.zeta() == CDF.zeta(Integer(2))
True
sage: CDF.zeta(0.5)
Traceback (most recent call last):
...
ValueError: n must be a positive integer
sage: CDF.zeta(0)
Traceback (most recent call last):
...
ValueError: n must be a positive integer
sage: CDF.zeta(-1)
Traceback (most recent call last):
...
ValueError: n must be a positive integer
>>> from sage.all import *
>>> CDF.zeta(RealNumber('0.5'))
Traceback (most recent call last):
...
ValueError: n must be a positive integer
>>> CDF.zeta(Integer(0))
Traceback (most recent call last):
...
ValueError: n must be a positive integer
>>> CDF.zeta(-Integer(1))
Traceback (most recent call last):
...
ValueError: n must be a positive integer
class sage.rings.complex_double.ComplexToCDF[source]#

Bases: Morphism

Fast morphism for anything such that the elements have attributes .real and .imag (e.g. numpy complex types).

EXAMPLES:

sage: # needs numpy
sage: import numpy
sage: f = CDF.coerce_map_from(numpy.complex_)
sage: f(numpy.complex_(I))
1.0*I
sage: f(numpy.complex_(I)).parent()
Complex Double Field
>>> from sage.all import *
>>> # needs numpy
>>> import numpy
>>> f = CDF.coerce_map_from(numpy.complex_)
>>> f(numpy.complex_(I))
1.0*I
>>> f(numpy.complex_(I)).parent()
Complex Double Field
class sage.rings.complex_double.FloatToCDF[source]#

Bases: Morphism

Fast morphism from anything with a __float__ method to a CDF element.

EXAMPLES:

sage: f = CDF.coerce_map_from(ZZ); f
Native morphism:
  From: Integer Ring
  To:   Complex Double Field
sage: f(4)
4.0
sage: f = CDF.coerce_map_from(QQ); f
Native morphism:
  From: Rational Field
  To:   Complex Double Field
sage: f(1/2)
0.5
sage: f = CDF.coerce_map_from(int); f
Native morphism:
  From: Set of Python objects of class 'int'
  To:   Complex Double Field
sage: f(3r)
3.0
sage: f = CDF.coerce_map_from(float); f
Native morphism:
  From: Set of Python objects of class 'float'
  To:   Complex Double Field
sage: f(3.5)
3.5
>>> from sage.all import *
>>> f = CDF.coerce_map_from(ZZ); f
Native morphism:
  From: Integer Ring
  To:   Complex Double Field
>>> f(Integer(4))
4.0
>>> f = CDF.coerce_map_from(QQ); f
Native morphism:
  From: Rational Field
  To:   Complex Double Field
>>> f(Integer(1)/Integer(2))
0.5
>>> f = CDF.coerce_map_from(int); f
Native morphism:
  From: Set of Python objects of class 'int'
  To:   Complex Double Field
>>> f(3)
3.0
>>> f = CDF.coerce_map_from(float); f
Native morphism:
  From: Set of Python objects of class 'float'
  To:   Complex Double Field
>>> f(RealNumber('3.5'))
3.5
sage.rings.complex_double.is_ComplexDoubleElement(x)[source]#

Return True if x is a ComplexDoubleElement.

EXAMPLES:

sage: from sage.rings.complex_double import is_ComplexDoubleElement
sage: is_ComplexDoubleElement(0)
doctest:warning...
DeprecationWarning: The function is_ComplexDoubleElement is deprecated;
use 'isinstance(..., ComplexDoubleElement)' instead.
See https://github.com/sagemath/sage/issues/38128 for details.
False
sage: is_ComplexDoubleElement(CDF(0))
True
>>> from sage.all import *
>>> from sage.rings.complex_double import is_ComplexDoubleElement
>>> is_ComplexDoubleElement(Integer(0))
doctest:warning...
DeprecationWarning: The function is_ComplexDoubleElement is deprecated;
use 'isinstance(..., ComplexDoubleElement)' instead.
See https://github.com/sagemath/sage/issues/38128 for details.
False
>>> is_ComplexDoubleElement(CDF(Integer(0)))
True