Arbitrary precision complex balls using Arb

This is a binding to the Arb library; it may be useful to refer to its documentation for more details.

Parts of the documentation for this module are copied or adapted from Arb’s own documentation, licenced under the GNU General Public License version 2, or later.

Data Structure

A ComplexBall represents a complex number with error bounds. It wraps an Arb object of type acb_t, which consists of a pair of real number balls representing the real and imaginary part with separate error bounds. (See the documentation of sage.rings.real_arb for more information.)

A ComplexBall thus represents a rectangle \([m_1-r_1, m_1+r_1] + [m_2-r_2, m_2+r_2] i\) in the complex plane. This is used in Arb instead of a disk or square representation (consisting of a complex floating-point midpoint with a single radius), since it allows implementing many operations more conveniently by splitting into ball operations on the real and imaginary parts. It also allows tracking when complex numbers have an exact (for example exactly zero) real part and an inexact imaginary part, or vice versa.

The parents of complex balls are instances of ComplexBallField. The name CBF is bound to the complex ball field with the default precision of 53 bits:

sage: CBF is ComplexBallField() is ComplexBallField(53)
True

Comparison

Warning

In accordance with the semantics of Arb, identical ComplexBall objects are understood to give permission for algebraic simplification. This assumption is made to improve performance. For example, setting z = x*x sets \(z\) to a ball enclosing the set \(\{t^2 : t \in x\}\) and not the (generally larger) set \(\{tu : t \in x, u \in x\}\).

Two elements are equal if and only if they are the same object or if both are exact and equal:

sage: a = CBF(1, 2)
sage: b = CBF(1, 2)
sage: a is b
False
sage: a == b
True
sage: a = CBF(1/3, 1/5)
sage: b = CBF(1/3, 1/5)
sage: a.is_exact()
False
sage: b.is_exact()
False
sage: a is b
False
sage: a == b
False

A ball is non-zero in the sense of usual comparison if and only if it does not contain zero:

sage: a = CBF(RIF(-0.5, 0.5))
sage: a != 0
False
sage: b = CBF(1/3, 1/5)
sage: b != 0
True

However, bool(b) returns False for a ball b only if b is exactly zero:

sage: bool(a)
True
sage: bool(b)
True
sage: bool(CBF.zero())
False

Coercion

Automatic coercions work as expected:

sage: bpol = 1/3*CBF(i) + AA(sqrt(2)) + (polygen(RealBallField(20), 'x') + QQbar(i))
sage: bpol
x + [1.41421 +/- 5.09e-6] + [1.33333 +/- 3.97e-6]*I
sage: bpol.parent()
Univariate Polynomial Ring in x over Complex ball field with 20 bits precision
sage: bpol/3
([0.333333 +/- 4.93e-7])*x + [0.47140 +/- 5.39e-6] + [0.44444 +/- 4.98e-6]*I
sage: SR.coerce(CBF(0.42 + 3.33*I))
[0.4200000000000000 +/- 1.56e-17] + [3.330000000000000 +/- 7.11e-17]*I

Check that trac ticket #19839 is fixed:

sage: log(SR(CBF(0.42))).pyobject().parent()
Complex ball field with 53 bits precision

Classes and Methods

class sage.rings.complex_arb.ComplexBall

Bases: sage.structure.element.RingElement

Hold one acb_t of the Arb library

EXAMPLES:

sage: a = ComplexBallField()(1, 1)
sage: a
1.000000000000000 + 1.000000000000000*I
above_abs()

Return an upper bound for the absolute value of this complex ball.

OUTPUT:

A ball with zero radius

EXAMPLES:

sage: b = ComplexBallField(8)(1+i).above_abs()
sage: b
[1.4 +/- 0.0219]
sage: b.is_exact()
True
sage: QQ(b)*128
182

See also

below_abs()

accuracy()

Return the effective relative accuracy of this ball measured in bits.

This is computed as if calling accuracy() on the real ball whose midpoint is the larger out of the real and imaginary midpoints of this complex ball, and whose radius is the larger out of the real and imaginary radii of this complex ball.

EXAMPLES:

sage: CBF(exp(I*pi/3)).accuracy()
52
sage: CBF(I/2).accuracy() == CBF.base().maximal_accuracy()
True
sage: CBF('nan', 'inf').accuracy() == -CBF.base().maximal_accuracy()
True
add_error(ampl)

Increase the radii of the real and imaginary parts by (an upper bound on) ampl.

If ampl is negative, the radii remain unchanged.

INPUT:

  • ampl - A real ball (or an object that can be coerced to a real ball).

OUTPUT:

A new complex ball.

EXAMPLES:

sage: CBF(1+i).add_error(10^-16)
[1.000000000000000 +/- 1.01e-16] + [1.000000000000000 +/- 1.01e-16]*I
agm1()

Return the arithmetic-geometric mean of 1 and self.

The arithmetic-geometric mean is defined such that the function is continuous in the complex plane except for a branch cut along the negative half axis (where it is continuous from above). This corresponds to always choosing an “optimal” branch for the square root in the arithmetic-geometric mean iteration.

EXAMPLES:

sage: CBF(0, -1).agm1()
[0.5990701173678 +/- 1.15e-14] + [-0.5990701173678 +/- 1.19e-14]*I
airy()

Return the Airy functions Ai, Ai’, Bi, Bi’ with argument self, evaluated simultaneously.

EXAMPLES:

sage: CBF(10*pi).airy()
([1.2408955946101e-52 +/- 5.50e-66],
 [-6.965048886977e-52 +/- 5.23e-65],
 [2.288295683344e+50 +/- 5.10e+37],
 [1.2807602335816e+51 +/- 4.97e+37])
sage: ai, aip, bi, bip = CBF(1,2).airy()
sage: (ai * bip - bi * aip) * CBF(pi)
[1.0000000000000 +/- 1.25e-15] + [+/- 3.27e-16]*I
airy_ai()

Return the Airy function Ai with argument self.

EXAMPLES:

sage: CBF(1,2).airy_ai()
[-0.2193862549814276 +/- 7.47e-17] + [-0.1753859114081094 +/- 4.06e-17]*I
airy_ai_prime()

Return the Airy function derivative Ai’ with argument self.

EXAMPLES:

sage: CBF(1,2).airy_ai_prime()
[0.1704449781789148 +/- 3.12e-17] + [0.387622439413295 +/- 1.06e-16]*I
airy_bi()

Return the Airy function Bi with argument self.

EXAMPLES:

sage: CBF(1,2).airy_bi()
[0.0488220324530612 +/- 1.30e-17] + [0.1332740579917484 +/- 6.25e-17]*I
airy_bi_prime()

Return the Airy function derivative Bi’ with argument self.

EXAMPLES:

sage: CBF(1,2).airy_bi_prime()
[-0.857239258605362 +/- 3.47e-16] + [0.4955063363095674 +/- 9.22e-17]*I
arccos()

Return the arccosine of this ball.

EXAMPLES:

sage: CBF(1+i).arccos()
[0.90455689430238 +/- 2.18e-15] + [-1.06127506190504 +/- 5.04e-15]*I
arccosh()

Return the hyperbolic arccosine of this ball.

EXAMPLES:

sage: CBF(1+i).arccosh()
[1.061275061905035 +/- 8.44e-16] + [0.904556894302381 +/- 8.22e-16]*I
arcsin()

Return the arcsine of this ball.

EXAMPLES:

sage: CBF(1+i).arcsin()
[0.66623943249252 +/- 5.40e-15] + [1.06127506190504 +/- 5.04e-15]*I
arcsinh()

Return the hyperbolic arcsine of this ball.

EXAMPLES:

sage: CBF(1+i).arcsinh()
[1.06127506190504 +/- 5.04e-15] + [0.66623943249252 +/- 5.40e-15]*I
arctan()

Return the arctangent of this ball.

EXAMPLES:

sage: CBF(1+i).arctan()
[1.017221967897851 +/- 4.93e-16] + [0.4023594781085251 +/- 8.52e-17]*I
sage: CBF(i).arctan()
nan + nan*I
arctanh()

Return the hyperbolic arctangent of this ball.

EXAMPLES:

sage: CBF(1+i).arctanh()
[0.4023594781085251 +/- 8.52e-17] + [1.017221967897851 +/- 4.93e-16]*I
arg()

Return the argument of this complex ball.

EXAMPLES:

sage: CBF(1 + i).arg()
[0.785398163397448 +/- 3.91e-16]
sage: CBF(-1).arg()
[3.141592653589793 +/- 5.61e-16]
sage: CBF(-1).arg().parent()
Real ball field with 53 bits precision
barnes_g()

Return the Barnes G-function of self.

EXAMPLES:

sage: CBF(-4).barnes_g()
0
sage: CBF(8).barnes_g()
24883200.00000000
sage: CBF(500,10).barnes_g()
[4.54078781e+254873 +/- 5.43e+254864] + [8.65835455e+254873 +/- 7.28e+254864]*I
below_abs(test_zero=False)

Return a lower bound for the absolute value of this complex ball.

INPUT:

  • test_zero (boolean, default False) – if True, make sure that the returned lower bound is positive, raising an error if the ball contains zero.

OUTPUT:

A ball with zero radius

EXAMPLES:

sage: b = ComplexBallField(8)(1+i).below_abs()
sage: b
[1.4 +/- 0.0141]
sage: b.is_exact()
True
sage: QQ(b)*128
181
sage: (CBF(1/3) - 1/3).below_abs()
0
sage: (CBF(1/3) - 1/3).below_abs(test_zero=True)
Traceback (most recent call last):
...
ValueError: ball contains zero

See also

above_abs()

bessel_I(nu)

Return the modified Bessel function of the first kind with argument self and index nu.

EXAMPLES:

sage: CBF(1, 1).bessel_I(1)
[0.365028028827088 +/- 3.96e-16] + [0.614160334922903 +/- 6.38e-16]*I
sage: CBF(100, -100).bessel_I(1/3)
[5.4362189595644e+41 +/- 6.40e+27] + [7.1989436985321e+41 +/- 2.92e+27]*I
bessel_J(nu)

Return the Bessel function of the first kind with argument self and index nu.

EXAMPLES:

sage: CBF(1, 1).bessel_J(1)
[0.614160334922903 +/- 6.38e-16] + [0.365028028827088 +/- 3.96e-16]*I
sage: CBF(100, -100).bessel_J(1/3)
[1.108431870251e+41 +/- 5.53e+28] + [-8.952577603125e+41 +/- 2.93e+28]*I
bessel_J_Y(nu)

Return the Bessel function of the first and second kind with argument self and index nu, computed simultaneously.

EXAMPLES:

sage: J, Y = CBF(1, 1).bessel_J_Y(1)
sage: J - CBF(1, 1).bessel_J(1)
[+/- 3.75e-16] + [+/- 2.64e-16]*I
sage: Y - CBF(1, 1).bessel_Y(1)
[+/- 1.64e-14] + [+/- 1.62e-14]*I
bessel_K(nu)

Return the modified Bessel function of the second kind with argument self and index nu.

EXAMPLES:

sage: CBF(1, 1).bessel_K(0)
[0.08019772694652 +/- 3.19e-15] + [-0.35727745928533 +/- 1.08e-15]*I
sage: CBF(1, 1).bessel_K(1)
[0.02456830552374 +/- 4.84e-15] + [-0.45971947380119 +/- 5.35e-15]*I
sage: CBF(100, 100).bessel_K(QQbar(i))
[3.8693896656383e-45 +/- 2.76e-59] + [5.507100423418e-46 +/- 4.01e-59]*I
bessel_Y(nu)

Return the Bessel function of the second kind with argument self and index nu.

EXAMPLES:

sage: CBF(1, 1).bessel_Y(1)
[-0.6576945355913 +/- 5.29e-14] + [0.6298010039929 +/- 2.45e-14]*I
sage: CBF(100, -100).bessel_Y(1/3)
[-8.952577603125e+41 +/- 4.66e+28] + [-1.108431870251e+41 +/- 6.31e+28]*I
chebyshev_T(n)

Return the Chebyshev function of the first kind of order n evaluated at self.

EXAMPLES:

sage: CBF(1/3).chebyshev_T(20)
[0.8710045668809 +/- 6.15e-14]
sage: CBF(1/3).chebyshev_T(CBF(5,1))
[1.8429685451876 +/- 3.57e-14] + [0.20053614301799 +/- 7.05e-15]*I
chebyshev_U(n)

Return the Chebyshev function of the second kind of order n evaluated at self.

EXAMPLES:

sage: CBF(1/3).chebyshev_U(20)
[0.6973126541184 +/- 2.83e-14]
sage: CBF(1/3).chebyshev_U(CBF(5,1))
[1.7588496489342 +/- 5.99e-14] + [0.7497317165104 +/- 4.35e-14]*I
chi()

Return the hyperbolic cosine integral with argument self.

EXAMPLES:

sage: CBF(1, 1).chi()
[0.882172180555936 +/- 4.85e-16] + [1.28354719327494 +/- 1.07e-15]*I
sage: CBF(0).chi()
nan + nan*I
ci()

Return the cosine integral with argument self.

EXAMPLES:

sage: CBF(1, 1).ci()
[0.882172180555936 +/- 4.85e-16] + [0.287249133519956 +/- 3.92e-16]*I
sage: CBF(0).ci()
nan + nan*I
conjugate()

Return the complex conjugate of this ball.

EXAMPLES:

sage: CBF(-2 + I/3).conjugate()
-2.000000000000000 + [-0.3333333333333333 +/- 7.04e-17]*I
contains_exact(other)

Return True iff other is contained in self.

Use other in self for a test that works for a wider range of inputs but may return false negatives.

INPUT:

EXAMPLES:

sage: CBF(RealBallField(100)(1/3), 0).contains_exact(1/3)
True
sage: CBF(1).contains_exact(1)
True
sage: CBF(1).contains_exact(CBF(1))
True

sage: CBF(sqrt(2)).contains_exact(sqrt(2))
Traceback (most recent call last):
...
TypeError: unsupported type: <type 'sage.symbolic.expression.Expression'>
contains_integer()

Return True iff this ball contains any integer.

EXAMPLES:

sage: CBF(3, RBF(0.1)).contains_integer()
False
sage: CBF(3, RBF(0.1,0.1)).contains_integer()
True
contains_zero()

Return True iff this ball contains zero.

EXAMPLES:

sage: CBF(0).contains_zero()
True
sage: CBF(RIF(-1,1)).contains_zero()
True
sage: CBF(i).contains_zero()
False
cos()

Return the cosine of this ball.

EXAMPLES:

sage: CBF(i*pi).cos()
[11.59195327552152 +/- 8.38e-15]
cot()

Return the cotangent of this ball.

EXAMPLES:

sage: CBF(pi, 1/10).cot()
[+/- 5.74e-14] + [-10.0333111322540 +/- 2.81e-14]*I
sage: CBF(pi).cot()
[+/- inf]
cube()

Return the cube of this ball.

The result is computed efficiently using two real squarings, two real multiplications, and scalar operations.

EXAMPLES:

sage: CBF(1, 1).cube()
-2.000000000000000 + 2.000000000000000*I
diameter()

Return the diameter of this ball.

EXAMPLES:

sage: CBF(1 + i).diameter()
0.00000000
sage: CBF(i/3).diameter()
2.2204460e-16
sage: CBF(i/3).diameter().parent()
Real Field with 30 bits of precision
sage: CBF(CIF(RIF(1.02, 1.04), RIF(2.1, 2.2))).diameter()
0.20000000

See also

rad(), mid()

ei()

Return the exponential integral with argument self.

EXAMPLES:

sage: CBF(1, 1).ei()
[1.76462598556385 +/- 5.82e-15] + [2.38776985151052 +/- 4.29e-15]*I
sage: CBF(0).ei()
nan
eisenstein(n)

Return the first n entries in the sequence of Eisenstein series \(G_4(\tau), G_6(\tau), G_8(\tau), \ldots\) where tau is given by self. The output is a list.

EXAMPLES:

sage: a, b, c, d = 2, 5, 1, 3
sage: tau = CBF(1,3)
sage: tau.eisenstein(4)
[[2.1646498507193 +/- 6.30e-14],
 [2.0346794456073 +/- 2.44e-14],
 [2.0081609898081 +/- 3.67e-14],
 [2.0019857082706 +/- 4.60e-14]]
sage: ((a*tau+b)/(c*tau+d)).eisenstein(3)[2]
[331011.2004330 +/- 9.33e-8] + [-711178.1655746 +/- 7.51e-8]*I
sage: (c*tau+d)^8 * tau.eisenstein(3)[2]
[331011.20043304 +/- 7.62e-9] + [-711178.1655746 +/- 1.34e-8]*I
elliptic_e()

Return the complete elliptic integral of the second kind evaluated at m given by self.

EXAMPLES:

sage: CBF(2,3).elliptic_e()
[1.472797144959 +/- 4.82e-13] + [-1.231604783936 +/- 1.25e-13]*I
elliptic_k()

Return the complete elliptic integral of the first kind evaluated at m given by self.

EXAMPLES:

sage: CBF(2,3).elliptic_k()
[1.0429132919285 +/- 3.65e-14] + [0.6296824723086 +/- 6.15e-14]*I
elliptic_p(tau, n=None)

Return the Weierstrass elliptic function with lattice parameter tau, evaluated at self. The function is doubly periodic in self with periods 1 and tau, which should lie in the upper half plane.

If n is given, return a list containing the first n terms in the Taylor expansion at self. In particular, with n = 2, compute the Weierstrass elliptic function together with its derivative, which generate the field of elliptic functions with periods 1 and tau.

EXAMPLES:

sage: tau = CBF(1,4)
sage: z = CBF(sqrt(2), sqrt(3))
sage: z.elliptic_p(tau)
[-3.28920996772709 +/- 7.63e-15] + [-0.0003673767302933 +/- 6.04e-17]*I
sage: (z + tau).elliptic_p(tau)
[-3.28920996772709 +/- 7.97e-15] + [-0.000367376730293 +/- 6.51e-16]*I
sage: (z + 1).elliptic_p(tau)
[-3.28920996772709 +/- 7.63e-15] + [-0.0003673767302933 +/- 6.04e-17]*I

sage: z.elliptic_p(tau, 3)
[[-3.28920996772709 +/- 7.62e-15] + [-0.0003673767302933 +/- 5.40e-17]*I,
 [0.002473055794309 +/- 5.01e-16] + [0.003859554040267 +/- 4.45e-16]*I,
 [-0.01299087561709 +/- 4.72e-15] + [0.00725027521915 +/- 4.32e-15]*I]
sage: (z + 3 + 4*tau).elliptic_p(tau, 3)
[[-3.2892099677271 +/- 2.29e-14] + [-0.00036737673029 +/- 8.58e-15]*I,
 [0.002473055794 +/- 6.59e-13] + [0.003859554040 +/- 6.17e-13]*I,
 [-0.0129908756 +/- 3.39e-11] + [0.0072502752 +/- 3.60e-11]*I]
erf()

Return the error function with argument self.

EXAMPLES:

sage: CBF(1, 1).erf()
[1.316151281697947 +/- 7.26e-16] + [0.1904534692378347 +/- 6.03e-17]*I
erfc()

Compute the complementary error function with argument self.

EXAMPLES:

sage: CBF(20).erfc()
[5.39586561160790e-176 +/- 6.73e-191]
sage: CBF(100, 100).erfc()
[0.00065234366376858 +/- 8.37e-18] + [-0.00393572636292141 +/- 7.21e-18]*I
exp()

Return the exponential of this ball.

See also

exppii()

EXAMPLES:

sage: CBF(i*pi).exp()
[-1.00000000000000 +/- 6.67e-16] + [+/- 5.68e-16]*I
exp_integral_e(s)

Return the image of this ball by the generalized exponential integral with index s.

EXAMPLES:

sage: CBF(1+i).exp_integral_e(1)
[0.00028162445198 +/- 2.79e-15] + [-0.17932453503936 +/- 2.12e-15]*I
sage: CBF(1+i).exp_integral_e(QQbar(i))
[-0.10396361883964 +/- 3.78e-15] + [-0.16268401277783 +/- 3.69e-15]*I
exppii()

Return exp(pi*i*self).

EXAMPLES:

sage: CBF(1/2).exppii()
1.000000000000000*I
sage: CBF(0, -1/pi).exppii()
[2.71828182845904 +/- 6.20e-15]
gamma(z=None)

Return the image of this ball by the Euler Gamma function (if z = None) or the incomplete Gamma function (otherwise).

EXAMPLES:

sage: CBF(1, 1).gamma()
[0.498015668118356 +/- 9.16e-16] + [-0.154949828301811 +/- 7.08e-16]*I
sage: CBF(-1).gamma()
nan
sage: CBF(1, 1).gamma(0)
[0.498015668118356 +/- 9.16e-16] + [-0.154949828301811 +/- 7.08e-16]*I
sage: CBF(1, 1).gamma(100)
[-3.6143867454139e-45 +/- 6.88e-59] + [-3.7022961377791e-44 +/- 4.41e-58]*I
sage: CBF(1, 1).gamma(CLF(i))
[0.32886684193500 +/- 5.04e-15] + [-0.18974945045621 +/- 1.26e-15]*I
gegenbauer_C(n, m)

Return the Gegenbauer polynomial (or function) \(C_n^m(z)\) evaluated at self.

EXAMPLES:

sage: CBF(-10).gegenbauer_C(7, 1/2)
[-263813415.6250000 +/- 9.57e-8]
hermite_H(n)

Return the Hermite function (or polynomial) of order n evaluated at self.

EXAMPLES:

sage: CBF(10).hermite_H(1)
20.00000000000000
sage: CBF(10).hermite_H(30)
[8.0574670961707e+37 +/- 3.28e+23]
hypergeometric(a, b, regularized=False)

Return the generalized hypergeometric function of self.

INPUT:

  • a – upper parameters, list of complex numbers that coerce into this ball’s parent;
  • b – lower parameters, list of complex numbers that coerce into this ball’s parent.
  • regularized – if True, the regularized generalized hypergeometric function is computed.

OUTPUT:

The generalized hypergeometric function defined by

\[{}_pF_q(a_1,\ldots,a_p;b_1,\ldots,b_q;z) = \sum_{k=0}^\infty \frac{(a_1)_k\dots(a_p)_k}{(b_1)_k\dots(b_q)_k} \frac {z^k} {k!}\]

extended using analytic continuation or regularization when the sum does not converge.

The regularized generalized hypergeometric function

\[{}_pF_q(a_1,\ldots,a_p;b_1,\ldots,b_q;z) = \sum_{k=0}^\infty \frac{(a_1)_k\dots(a_p)_k}{\Gamma(b_1+k)\dots\Gamma(b_q+k)} \frac {z^k} {k!}\]

is well-defined even when the lower parameters are nonpositive integers. Currently, this is only supported for some \(p\) and \(q\).

EXAMPLES:

sage: CBF(1, pi/2).hypergeometric([], [])
[+/- 7.72e-16] + [2.71828182845904 +/- 6.45e-15]*I

sage: CBF(1, pi).hypergeometric([1/4], [1/4])
[-2.7182818284590 +/- 7.11e-14] + [+/- 2.25e-14]*I

sage: CBF(1000, 1000).hypergeometric([10], [AA(sqrt(2))])
[9.79300951360e+454 +/- 5.01e+442] + [5.522579106816e+455 +/- 3.56e+442]*I
sage: CBF(1000, 1000).hypergeometric([100], [AA(sqrt(2))])
[1.27967355557e+590 +/- 8.60e+578] + [-9.32333491987e+590 +/- 8.18e+578]*I

sage: CBF(0, 1).hypergeometric([], [1/2, 1/3, 1/4])
[-3.7991962344383 +/- 8.78e-14] + [23.878097177805 +/- 3.87e-13]*I

sage: CBF(0).hypergeometric([1], [])
1.000000000000000
sage: CBF(1, 1).hypergeometric([1], [])
1.000000000000000*I

sage: CBF(2+3*I).hypergeometric([1/4,1/3],[1/2])
[0.7871684267473 +/- 7.34e-14] + [0.2749254173721 +/- 9.23e-14]*I
sage: CBF(2+3*I).hypergeometric([1/4,1/3],[1/2],regularized=True)
[0.4441122268685 +/- 3.96e-14] + [0.1551100567338 +/- 5.75e-14]*I

sage: CBF(5).hypergeometric([2,3], [-5])
nan + nan*I
sage: CBF(5).hypergeometric([2,3], [-5], regularized=True)
[5106.925964355 +/- 5.41e-10]

sage: CBF(2016).hypergeometric([], [2/3])
[2.025642692328e+38 +/- 3.00e+25]
sage: CBF(-2016).hypergeometric([], [2/3], regularized=True)
[-0.0005428550847 +/- 5.00e-14]

sage: CBF(-7).hypergeometric([4], [])
0.0002441406250000000

sage: CBF(0, 3).hypergeometric([CBF(1,1)], [-4], regularized=True)
[239.514000752841 +/- 8.03e-13] + [105.175157349015 +/- 6.28e-13]*I
hypergeometric_U(a, b)

Return the Tricomi confluent hypergeometric function U(a, b, self) of this ball.

EXAMPLES:

sage: CBF(1000, 1000).hypergeometric_U(RLF(pi), -100)
[-7.261605907166e-11 +/- 5.04e-24] + [-7.928136216391e-11 +/- 5.52e-24]*I
sage: CBF(1000, 1000).hypergeometric_U(0, -100)
1.000000000000000
identical(other)

Return whether self and other represent the same ball.

INPUT:

OUTPUT:

Return True iff self and other are equal as sets, i.e. if their real and imaginary parts each have the same midpoint and radius.

Note that this is not the same thing as testing whether both self and other certainly represent the complex real number, unless either self or other is exact (and neither contains NaN). To test whether both operands might represent the same mathematical quantity, use overlaps() or in, depending on the circumstance.

EXAMPLES:

sage: CBF(1, 1/3).identical(1 + CBF(0, 1)/3)
True
sage: CBF(1, 1).identical(1 + CBF(0, 1/3)*3)
False
imag()

Return the imaginary part of this ball.

OUTPUT:

A RealBall.

EXAMPLES:

sage: a = CBF(1/3, 1/5)
sage: a.imag()
[0.2000000000000000 +/- 4.45e-17]
is_NaN()

Return True iff either the real or the imaginary part is not-a-number.

EXAMPLES:

sage: CBF(NaN).is_NaN()
True
sage: CBF(-5).gamma().is_NaN()
True
sage: CBF(oo).is_NaN()
False
sage: CBF(42+I).is_NaN()
False
is_exact()

Return True iff the radius of this ball is zero.

EXAMPLES:

sage: CBF(1).is_exact()
True
sage: CBF(1/3, 1/3).is_exact()
False
is_nonzero()

Return True iff zero is not contained in the interval represented by this ball.

Note

This method is not the negation of is_zero(): it only returns True if zero is known not to be contained in the ball.

Use bool(b) (or, equivalently, not b.is_zero()) to check if a ball b may represent a nonzero number (for instance, to determine the “degree” of a polynomial with ball coefficients).

EXAMPLES:

sage: CBF(pi, 1/3).is_nonzero()
True
sage: CBF(RIF(-0.5, 0.5), 1/3).is_nonzero()
True
sage: CBF(1/3, RIF(-0.5, 0.5)).is_nonzero()
True
sage: CBF(RIF(-0.5, 0.5), RIF(-0.5, 0.5)).is_nonzero()
False

See also

is_zero()

is_real()

Return True iff the imaginary part of this ball is exactly zero.

EXAMPLES:

sage: CBF(1/3, 0).is_real()
True
sage: (CBF(i/3) - CBF(1, 1/3)).is_real()
False
sage: CBF('inf').is_real()
True
is_zero()

Return True iff the midpoint and radius of this ball are both zero.

EXAMPLES:

sage: CBF(0).is_zero()
True
sage: CBF(RIF(-0.5, 0.5)).is_zero()
False

See also

is_nonzero()

jacobi_P(n, a, b)

Return the Jacobi polynomial (or function) \(P_n^{(a,b)}(z)\) evaluated at self.

EXAMPLES:

sage: CBF(5,-6).jacobi_P(8, CBF(1,2), CBF(2,3))
[-920983000.45982 +/- 2.22e-6] + [6069919969.92857 +/- 4.77e-6]*I
jacobi_theta(tau)

Return the four Jacobi theta functions evaluated at the argument self (representing \(z\)) and the parameter tau which should lie in the upper half plane.

The following definitions are used:

\[ \begin{align}\begin{aligned}\theta_1(z,\tau) = 2 q_{1/4} \sum_{n=0}^{\infty} (-1)^n q^{n(n+1)} \sin((2n+1) \pi z)\\\theta_2(z,\tau) = 2 q_{1/4} \sum_{n=0}^{\infty} q^{n(n+1)} \cos((2n+1) \pi z)\\\theta_3(z,\tau) = 1 + 2 \sum_{n=1}^{\infty} q^{n^2} \cos(2n \pi z)\\\theta_4(z,\tau) = 1 + 2 \sum_{n=1}^{\infty} (-1)^n q^{n^2} \cos(2n \pi z)\end{aligned}\end{align} \]

where \(q = \exp(\pi i \tau)\) and \(q_{1/4} = \exp(\pi i \tau / 4)\). Note that \(z\) is multiplied by \(\pi\); some authors omit this factor.

EXAMPLES:

sage: CBF(3,-1/2).jacobi_theta(CBF(1/4,2))
([-0.186580562274757 +/- 5.52e-16] + [0.93841744788594 +/- 2.48e-15]*I,
 [-1.02315311037951 +/- 4.10e-15] + [-0.203600094532010 +/- 7.33e-16]*I,
 [1.030613911309632 +/- 4.25e-16] + [0.030613917822067 +/- 1.89e-16]*I,
 [0.969386075665498 +/- 4.65e-16] + [-0.030613917822067 +/- 1.89e-16]*I)

sage: CBF(3,-1/2).jacobi_theta(CBF(1/4,-2))
(nan + nan*I, nan + nan*I, nan + nan*I, nan + nan*I)

sage: CBF(0).jacobi_theta(CBF(0,1))
(0,
 [0.91357913815612 +/- 3.96e-15],
 [1.086434811213308 +/- 8.16e-16],
 [0.913579138156117 +/- 8.89e-16])
laguerre_L(n, m=0)

Return the Laguerre polynomial (or function) \(L_n^m(z)\) evaluated at self.

EXAMPLES:

sage: CBF(10).laguerre_L(3)
[-45.6666666666666 +/- 9.28e-14]
sage: CBF(10).laguerre_L(3, 2)
[-6.666666666667 +/- 4.15e-13]
sage: CBF(5,7).laguerre_L(CBF(2,3), CBF(1,-2))
[5515.315030271 +/- 4.37e-10] + [-12386.942845271 +/- 5.47e-10]*I
legendre_P(n, m=0, type=2)

Return the Legendre function of the first kind \(P_n^m(z)\) evaluated at self.

The type parameter can be either 2 or 3. This selects between different branch cut conventions. The definitions of the “type 2” and “type 3” functions are the same as those used by Mathematica and mpmath.

EXAMPLES:

sage: CBF(1/2).legendre_P(5)
0.08984375000000000
sage: CBF(1,2).legendre_P(CBF(2,3), CBF(0,1))
[0.10996180744364 +/- 7.45e-15] + [0.14312767804055 +/- 8.38e-15]*I
sage: CBF(-10).legendre_P(5, 325/100)
[-22104403.487377 +/- 6.81e-7] + [53364750.687392 +/- 7.25e-7]*I
sage: CBF(-10).legendre_P(5, 325/100, type=3)
[-57761589.914581 +/- 6.99e-7] + [+/- 5.14e-7]*I
legendre_Q(n, m=0, type=2)

Return the Legendre function of the second kind \(Q_n^m(z)\) evaluated at self.

The type parameter can be either 2 or 3. This selects between different branch cut conventions. The definitions of the “type 2” and “type 3” functions are the same as those used by Mathematica and mpmath.

EXAMPLES:

sage: CBF(1/2).legendre_Q(5)
[0.55508089057168 +/- 2.79e-15]
sage: CBF(1,2).legendre_Q(CBF(2,3), CBF(0,1))
[0.167678710 +/- 4.60e-10] + [-0.161558598 +/- 7.47e-10]*I
sage: CBF(-10).legendre_Q(5, 325/100)
[-83825154.36008 +/- 4.94e-6] + [-34721515.80396 +/- 5.40e-6]*I
sage: CBF(-10).legendre_Q(5, 325/100, type=3)
[-4.797306921692e-6 +/- 6.82e-19] + [-4.797306921692e-6 +/- 6.57e-19]*I
li(offset=False)

Return the logarithmic integral with argument self.

If offset is True, return the offset logarithmic integral.

EXAMPLES:

sage: CBF(1, 1).li()
[0.61391166922120 +/- 6.40e-15] + [2.05958421419258 +/- 5.61e-15]*I
sage: CBF(0).li()
0
sage: CBF(0).li(offset=True)
[-1.045163780117493 +/- 5.54e-16]
sage: li(0).n()
0.000000000000000
sage: Li(0).n()
-1.04516378011749
log(base=None)

General logarithm (principal branch).

INPUT:

  • base (optional, complex ball or number) – if None, return the principal branch of the natural logarithm ln(self), otherwise, return the general logarithm ln(self)/ln(base)

EXAMPLES:

sage: CBF(2*i).log()
[0.6931471805599453 +/- 4.16e-17] + [1.570796326794897 +/- 6.65e-16]*I
sage: CBF(-1).log()
[3.141592653589793 +/- 5.61e-16]*I

sage: CBF(2*i).log(2)
[1.000000000000000 +/- 8.01e-17] + [2.26618007091360 +/- 4.23e-15]*I
sage: CBF(2*i).log(CBF(i))
[1.000000000000000 +/- 2.83e-16] + [-0.441271200305303 +/- 2.82e-16]*I

sage: CBF('inf').log()
nan + nan*I
sage: CBF(2).log(0)
nan + nan*I
log1p()

Return log(1 + self), computed accurately when self is close to zero.

EXAMPLES:

sage: eps = RBF(1e-50)
sage: CBF(1+eps, eps).log()
[+/- 2.23e-16] + [1.000000000000000e-50 +/- 2.30e-66]*I
sage: CBF(eps, eps).log1p()
[1.000000000000000e-50 +/- 7.63e-68] + [1.00000000000000e-50 +/- 2.30e-66]*I
log_barnes_g()

Return the logarithmic Barnes G-function of self.

EXAMPLES:

sage: CBF(10^100).log_barnes_g()
[1.14379254649702e+202 +/- 4.09e+187]
sage: CBF(0,1000).log_barnes_g()
[-2702305.04929258 +/- 2.60e-9] + [-790386.325561423 +/- 9.72e-10]*I
log_gamma()

Return the image of this ball by the logarithmic Gamma function.

The branch cut of the logarithmic gamma function is placed on the negative half-axis, which means that log_gamma(z) + log z = log_gamma(z+1) holds for all \(z\), whereas log_gamma(z) != log(gamma(z)) in general.

EXAMPLES:

sage: CBF(1000, 1000).log_gamma()
[5466.22252162990 +/- 3.05e-12] + [7039.33429191119 +/- 3.81e-12]*I
sage: CBF(-1/2).log_gamma()
[1.265512123484645 +/- 8.82e-16] + [-3.141592653589793 +/- 5.68e-16]*I
sage: CBF(-1).log_gamma()
nan + [-3.141592653589793 +/- 5.68e-16]*I
mid()

Return the midpoint of this ball.

OUTPUT:

ComplexNumber, floating-point complex number formed by the centers of the real and imaginary parts of this ball.

EXAMPLES:

sage: CBF(1/3, 1).mid()
0.333333333333333 + 1.00000000000000*I
sage: CBF(1/3, 1).mid().parent()
Complex Field with 53 bits of precision
sage: CBF('inf', 'nan').mid()
+infinity + NaN*I
sage: CBF('nan', 'inf').mid()
NaN + +infinity*I
sage: CBF('nan').mid()
NaN
sage: CBF('inf').mid()
+infinity
sage: CBF(0, 'inf').mid()
+infinity*I

See also

squash()

modular_delta()

Return the modular discriminant with tau given by self.

EXAMPLES:

sage: CBF(0,1).modular_delta()
[0.0017853698506421 +/- 6.15e-17]
sage: a, b, c, d = 2, 5, 1, 3
sage: tau = CBF(1,3)
sage: ((a*tau+b)/(c*tau+d)).modular_delta()
[0.20921376655 +/- 6.94e-12] + [1.5761192552 +/- 3.47e-11]*I
sage: (c*tau+d)^12 * tau.modular_delta()
[0.20921376654986 +/- 4.89e-15] + [1.5761192552253 +/- 4.45e-14]*I
modular_eta()

Return the Dedekind eta function with tau given by self.

EXAMPLES:

sage: CBF(0,1).modular_eta()
[0.768225422326057 +/- 9.18e-16]
sage: CBF(12,1).modular_eta()
[-0.768225422326057 +/- 9.18e-16]
modular_j()

Return the modular j-invariant with tau given by self.

EXAMPLES:

sage: CBF(0,1).modular_j()
[1728.0000000000 +/- 5.33e-11]
modular_lambda()

Return the modular lambda function with tau given by self.

EXAMPLES:

sage: tau = CBF(sqrt(2),pi)
sage: tau.modular_lambda()
[-0.00022005123884157 +/- 6.41e-18] + [-0.0007978734645994 +/- 5.15e-17]*I
sage: (tau + 2).modular_lambda()
[-0.00022005123884157 +/- 6.41e-18] + [-0.0007978734645994 +/- 5.15e-17]*I
sage: (tau / (1 - 2*tau)).modular_lambda()
[-0.00022005123884 +/- 2.53e-15] + [-0.00079787346460 +/- 2.85e-15]*I
overlaps(other)

Return True iff self and other have some point in common.

INPUT:

EXAMPLES:

sage: CBF(1, 1).overlaps(1 + CBF(0, 1/3)*3)
True
sage: CBF(1, 1).overlaps(CBF(1, 'nan'))
True
sage: CBF(1, 1).overlaps(CBF(0, 'nan'))
False
polylog(s)

Return the polylogarithm \(\operatorname{Li}_s(\mathrm{self})\).

EXAMPLES:

sage: CBF(2).polylog(1)
[+/- 4.65e-15] + [-3.14159265358979 +/- 8.15e-15]*I
sage: CBF(1, 1).polylog(CBF(1, 1))
[0.3708160030469 +/- 2.38e-14] + [2.7238016577979 +/- 4.20e-14]*I
psi(n=None)

Compute the digamma function with argument self.

If n is provided, compute the polygamma function of order n and argument self.

EXAMPLES:

sage: CBF(1, 1).psi()
[0.0946503206224770 +/- 7.74e-17] + [1.076674047468581 +/- 2.58e-16]*I
sage: CBF(-1).psi()
nan
sage: CBF(1,1).psi(10)
[56514.8269344249 +/- 4.70e-11] + [56215.1218005823 +/- 5.70e-11]*I
rad()

Return an upper bound for the error radius of this ball.

OUTPUT:

A RealNumber of the same precision as the radii of real balls.

Warning

Unlike a RealBall, a ComplexBall is not defined by its midpoint and radius. (Instances of ComplexBall are actually rectangles, not balls.)

EXAMPLES:

sage: CBF(1 + i).rad()
0.00000000
sage: CBF(i/3).rad()
1.1102230e-16
sage: CBF(i/3).rad().parent()
Real Field with 30 bits of precision

See also

diameter(), mid()

real()

Return the real part of this ball.

OUTPUT:

A RealBall.

EXAMPLES:

sage: a = CBF(1/3, 1/5)
sage: a.real()
[0.3333333333333333 +/- 7.04e-17]
rgamma()

Compute the reciprocal gamma function with argument self.

EXAMPLES:

sage: CBF(6).rgamma()
[0.00833333333333333 +/- 4.96e-18]
sage: CBF(-1).rgamma()
0
rising_factorial(n)

Return the n-th rising factorial of this ball.

The \(n\)-th rising factorial of \(x\) is equal to \(x (x+1) \cdots (x+n-1)\).

For complex \(n\), it is a quotient of gamma functions.

EXAMPLES:

sage: CBF(1).rising_factorial(5)
120.0000000000000
sage: CBF(1/3, 1/2).rising_factorial(300)
[-3.87949484514e+612 +/- 5.23e+600] + [-3.52042209763e+612 +/- 5.55e+600]*I

sage: CBF(1).rising_factorial(-1)
nan
sage: CBF(1).rising_factorial(2**64)
[+/- 2.30e+347382171305201370464]
sage: ComplexBallField(128)(1).rising_factorial(2**64)
[2.343691126796861348e+347382171305201285713 +/- 4.71e+347382171305201285694]
sage: CBF(1/2).rising_factorial(CBF(2,3))
[-0.123060451458124 +/- 4.43e-16] + [0.040641263167655 +/- 3.72e-16]*I
round()

Return a copy of this ball rounded to the precision of the parent.

EXAMPLES:

It is possible to create balls whose midpoint is more precise that their parent’s nominal precision (see real_arb for more information):

sage: b = CBF(exp(I*pi/3).n(100))
sage: b.mid()
0.50000000000000000000000000000 + 0.86602540378443864676372317075*I

The round() method rounds such a ball to its parent’s precision:

sage: b.round().mid()
0.500000000000000 + 0.866025403784439*I

See also

trim()

rsqrt()

Return the reciprocal square root of self.

If either the real or imaginary part is exactly zero, only a single real reciprocal square root is needed.

EXAMPLES:

sage: CBF(-2).rsqrt()
[-0.707106781186547 +/- 5.73e-16]*I
sage: CBF(0, 1/2).rsqrt()
1.000000000000000 - 1.000000000000000*I
sage: CBF(0).rsqrt()
nan
shi()

Return the hyperbolic sine integral with argument self.

EXAMPLES:

sage: CBF(1, 1).shi()
[0.88245380500792 +/- 3.36e-15] + [1.10422265823558 +/- 2.48e-15]*I
sage: CBF(0).shi()
0
si()

Return the sine integral with argument self.

EXAMPLES:

sage: CBF(1, 1).si()
[1.10422265823558 +/- 2.48e-15] + [0.88245380500792 +/- 3.36e-15]*I
sage: CBF(0).si()
0
sin()

Return the sine of this ball.

EXAMPLES:

sage: CBF(i*pi).sin()
[11.5487393572577 +/- 5.34e-14]*I
spherical_harmonic(phi, n, m)

Return the spherical harmonic \(Y_n^m(\theta,\phi)\) evaluated at \(\theta\) given by self. In the current implementation, n and m must be small integers.

EXAMPLES:

sage: CBF(1+I).spherical_harmonic(1/2, -3, -2)
[0.80370071745224 +/- 4.02e-15] + [-0.07282031864711 +/- 4.69e-15]*I
sqrt()

Return the square root of this ball.

If either the real or imaginary part is exactly zero, only a single real square root is needed.

EXAMPLES:

sage: CBF(-2).sqrt()
[1.414213562373095 +/- 2.99e-16]*I
squash()

Return an exact ball with the same midpoint as this ball.

OUTPUT:

A ComplexBall.

EXAMPLES:

sage: mid = CBF(1/3, 1/10).squash()
sage: mid
[0.3333333333333333 +/- 1.49e-17] + [0.09999999999999999 +/- 1.68e-18]*I
sage: mid.parent()
Complex ball field with 53 bits precision
sage: mid.is_exact()
True

See also

mid()

tan()

Return the tangent of this ball.

EXAMPLES:

sage: CBF(pi/2, 1/10).tan()
[+/- 2.87e-14] + [10.0333111322540 +/- 2.36e-14]*I
sage: CBF(pi/2).tan()
[+/- inf]
trim()

Return a trimmed copy of this ball.

Return a copy of this ball with both the real and imaginary parts trimmed (see trim()).

EXAMPLES:

sage: b = CBF(1/3, RBF(1/3, rad=.01))
sage: b.mid()
0.333333333333333 + 0.333333333333333*I
sage: b.trim().mid()
0.333333333333333 + 0.333333015441895*I

See also

round()

zeta(a=None)

Return the image of this ball by the Hurwitz zeta function.

For a = None, this computes the Riemann zeta function.

EXAMPLES:

sage: CBF(1, 1).zeta()
[0.5821580597520036 +/- 5.26e-17] + [-0.9268485643308071 +/- 2.80e-17]*I
sage: CBF(1, 1).zeta(1)
[0.5821580597520036 +/- 5.26e-17] + [-0.9268485643308071 +/- 2.80e-17]*I
sage: CBF(1, 1).zeta(1/2)
[1.497919876084167 +/- 2.91e-16] + [0.2448655353684164 +/- 4.22e-17]*I
sage: CBF(1, 1).zeta(CBF(1, 1))
[-0.3593983122202835 +/- 3.01e-17] + [-2.875283329756940 +/- 4.52e-16]*I
sage: CBF(1, 1).zeta(-1)
nan + nan*I
class sage.rings.complex_arb.ComplexBallField(precision, category)

Bases: sage.structure.unique_representation.UniqueRepresentation, sage.rings.ring.Field

An approximation of the field of complex numbers using pairs of mid-rad intervals.

INPUT:

  • precision – an integer \(\ge 2\).

EXAMPLES:

sage: CBF(1)
1.000000000000000
Element

alias of ComplexBall

characteristic()

Complex ball fields have characteristic zero.

EXAMPLES:

sage: ComplexBallField().characteristic()
0
complex_field()

Return the complex ball field with the same precision, i.e. self

EXAMPLES:

sage: CBF.complex_field() is CBF
True
construction()

Return the construction of a complex ball field as the algebraic closure of the real ball field with the same precision.

EXAMPLES:

sage: functor, base = CBF.construction()
sage: functor, base
(AlgebraicClosureFunctor, Real ball field with 53 bits precision)
sage: functor(base) is CBF
True
gen(i)

For i = 0, return the imaginary unit in this complex ball field.

EXAMPLES:

sage: CBF.0
1.000000000000000*I
sage: CBF.gen(1)
Traceback (most recent call last):
...
ValueError: only one generator
gens()

Return the tuple of generators of this complex ball field, i.e. (i,).

EXAMPLES:

sage: CBF.gens()
(1.000000000000000*I,)
sage: CBF.gens_dict()
{'1.000000000000000*I': 1.000000000000000*I}
is_exact()

Complex ball fields are not exact.

EXAMPLES:

sage: ComplexBallField().is_exact()
False
is_finite()

Complex ball fields are infinite.

They already specify it via their category, but we currently need to re-implement this method due to the legacy implementation in sage.rings.ring.Ring.

EXAMPLES:

sage: ComplexBallField().is_finite()
False
ngens()

Return 1 as the only generator is the imaginary unit.

EXAMPLES:

sage: CBF.ngens()
1
precision()

Return the bit precision used for operations on elements of this field.

EXAMPLES:

sage: ComplexBallField().precision()
53
some_elements()

Complex ball fields contain elements with exact, inexact, infinite, or undefined real and imaginary parts.

EXAMPLES:

sage: CBF.some_elements()
[1.000000000000000,
 -0.5000000000000000*I,
 1.000000000000000 + [0.3333333333333333 +/- 1.49e-17]*I,
 [-0.3333333333333333 +/- 1.49e-17] + 0.2500000000000000*I,
 [-2.175556475109056e+181961467118333366510562 +/- 1.29e+181961467118333366510545],
 [+/- inf],
 [0.3333333333333333 +/- 1.49e-17] + [+/- inf]*I,
 [+/- inf] + [+/- inf]*I,
 nan,
 nan + nan*I,
 [+/- inf] + nan*I]