Arbitrary precision real balls#

This is a binding to the arb module of FLINT. 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 (now FLINT’s) own documentation, licenced (at the time) under the GNU General Public License version 2, or later.

Data Structure#

Ball arithmetic, also known as mid-rad interval arithmetic, is an extension of floating-point arithmetic in which an error bound is attached to each variable. This allows doing rigorous computations over the real numbers, while avoiding the overhead of traditional (inf-sup) interval arithmetic at high precision, and eliminating much of the need for time-consuming and bug-prone manual error analysis associated with standard floating-point arithmetic.

Sage RealBall objects wrap FLINT objects of type arb_t. A real ball represents a ball over the real numbers, that is, an interval \([m-r,m+r]\) where the midpoint \(m\) and the radius \(r\) are (extended) real numbers:

sage: RBF(pi)                                                                       # needs sage.symbolic
[3.141592653589793 +/- ...e-16]
sage: RBF(pi).mid(), RBF(pi).rad()                                                  # needs sage.symbolic
(3.14159265358979, ...e-16)
>>> from sage.all import *
>>> RBF(pi)                                                                       # needs sage.symbolic
[3.141592653589793 +/- ...e-16]
>>> RBF(pi).mid(), RBF(pi).rad()                                                  # needs sage.symbolic
(3.14159265358979, ...e-16)

The midpoint is represented as an arbitrary-precision floating-point number with arbitrary-precision exponent. The radius is a floating-point number with fixed-precision mantissa and arbitrary-precision exponent.

sage: RBF(2)^(2^100)
[2.285367694229514e+381600854690147056244358827360 +/- ...e+381600854690147056244358827344]
>>> from sage.all import *
>>> RBF(Integer(2))**(Integer(2)**Integer(100))
[2.285367694229514e+381600854690147056244358827360 +/- ...e+381600854690147056244358827344]

RealBallField objects (the parents of real balls) model the field of real numbers represented by balls on which computations are carried out with a certain precision:

sage: RBF
Real ball field with 53 bits of precision
>>> from sage.all import *
>>> RBF
Real ball field with 53 bits of precision

It is possible to construct a ball whose parent is the real ball field with precision \(p\) but whose midpoint does not fit on \(p\) bits. However, the results of operations involving such a ball will (usually) be rounded to its parent’s precision:

sage: RBF(factorial(50)).mid(), RBF(factorial(50)).rad()
(3.0414093201713378043612608166064768844377641568961e64, 0.00000000)
sage: (RBF(factorial(50)) + 0).mid()
3.04140932017134e64
>>> from sage.all import *
>>> RBF(factorial(Integer(50))).mid(), RBF(factorial(Integer(50))).rad()
(3.0414093201713378043612608166064768844377641568961e64, 0.00000000)
>>> (RBF(factorial(Integer(50))) + Integer(0)).mid()
3.04140932017134e64

Comparison#

Warning

In accordance with the semantics of FLINT/Arb, identical RealBall objects are understood to give permission for algebraic simplification. This assumption is made to improve performance. For example, setting z = x*x may set \(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 exact and equal (in spite of the above warning, inexact balls are not considered equal to themselves):

sage: a = RBF(1)
sage: b = RBF(1)
sage: a is b
False
sage: a == a
True
sage: a == b
True
>>> from sage.all import *
>>> a = RBF(Integer(1))
>>> b = RBF(Integer(1))
>>> a is b
False
>>> a == a
True
>>> a == b
True
sage: a = RBF(1/3)
sage: b = RBF(1/3)
sage: a.is_exact()
False
sage: b.is_exact()
False
sage: a is b
False
sage: a == a
False
sage: a == b
False
>>> from sage.all import *
>>> a = RBF(Integer(1)/Integer(3))
>>> b = RBF(Integer(1)/Integer(3))
>>> a.is_exact()
False
>>> b.is_exact()
False
>>> a is b
False
>>> a == a
False
>>> a == b
False

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

sage: a = RBF(RIF(-0.5, 0.5))
sage: a != 0
False
sage: b = RBF(1/3)
sage: b != 0
True
>>> from sage.all import *
>>> a = RBF(RIF(-RealNumber('0.5'), RealNumber('0.5')))
>>> a != Integer(0)
False
>>> b = RBF(Integer(1)/Integer(3))
>>> b != Integer(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(RBF.zero())
False
>>> from sage.all import *
>>> bool(a)
True
>>> bool(b)
True
>>> bool(RBF.zero())
False

A ball left is less than a ball right if all elements of left are less than all elements of right.

sage: a = RBF(RIF(1, 2))
sage: b = RBF(RIF(3, 4))
sage: a < b
True
sage: a <= b
True
sage: a > b
False
sage: a >= b
False
sage: a = RBF(RIF(1, 3))
sage: b = RBF(RIF(2, 4))
sage: a < b
False
sage: a <= b
False
sage: a > b
False
sage: a >= b
False
>>> from sage.all import *
>>> a = RBF(RIF(Integer(1), Integer(2)))
>>> b = RBF(RIF(Integer(3), Integer(4)))
>>> a < b
True
>>> a <= b
True
>>> a > b
False
>>> a >= b
False
>>> a = RBF(RIF(Integer(1), Integer(3)))
>>> b = RBF(RIF(Integer(2), Integer(4)))
>>> a < b
False
>>> a <= b
False
>>> a > b
False
>>> a >= b
False

Comparisons with Sage symbolic infinities work with some limitations:

sage: -infinity < RBF(1) < +infinity
True
sage: -infinity < RBF(infinity)
True
sage: RBF(infinity) < infinity
False
sage: RBF(NaN) < infinity                                                           # needs sage.symbolic
Traceback (most recent call last):
...
ValueError: infinite but not with +/- phase
sage: 1/RBF(0) <= infinity
Traceback (most recent call last):
...
ValueError: infinite but not with +/- phase
>>> from sage.all import *
>>> -infinity < RBF(Integer(1)) < +infinity
True
>>> -infinity < RBF(infinity)
True
>>> RBF(infinity) < infinity
False
>>> RBF(NaN) < infinity                                                           # needs sage.symbolic
Traceback (most recent call last):
...
ValueError: infinite but not with +/- phase
>>> Integer(1)/RBF(Integer(0)) <= infinity
Traceback (most recent call last):
...
ValueError: infinite but not with +/- phase

Comparisons between elements of real ball fields, however, support special values and should be preferred:

sage: RBF(NaN) < RBF(infinity)                                                      # needs sage.symbolic
False
sage: RBF(0).add_error(infinity) <= RBF(infinity)
True
>>> from sage.all import *
>>> RBF(NaN) < RBF(infinity)                                                      # needs sage.symbolic
False
>>> RBF(Integer(0)).add_error(infinity) <= RBF(infinity)
True

Classes and Methods#

class sage.rings.real_arb.RealBall[source]#

Bases: RingElement

Hold one arb_t

EXAMPLES:

sage: a = RealBallField()(RIF(1))                     # indirect doctest
sage: b = a.psi()
sage: b # abs tol 1e-15
[-0.5772156649015329 +/- 4.84e-17]
sage: RIF(b)
-0.577215664901533?
>>> from sage.all import *
>>> a = RealBallField()(RIF(Integer(1)))                     # indirect doctest
>>> b = a.psi()
>>> b # abs tol 1e-15
[-0.5772156649015329 +/- 4.84e-17]
>>> RIF(b)
-0.577215664901533?
Chi()[source]#

Hyperbolic cosine integral

EXAMPLES:

sage: RBF(1).Chi()  # abs tol 1e-17
[0.837866940980208 +/- 4.72e-16]
>>> from sage.all import *
>>> RBF(Integer(1)).Chi()  # abs tol 1e-17
[0.837866940980208 +/- 4.72e-16]
Ci()[source]#

Cosine integral

EXAMPLES:

sage: RBF(1).Ci()  # abs tol 5e-16
[0.337403922900968 +/- 3.25e-16]
>>> from sage.all import *
>>> RBF(Integer(1)).Ci()  # abs tol 5e-16
[0.337403922900968 +/- 3.25e-16]
Ei()[source]#

Exponential integral

EXAMPLES:

sage: RBF(1).Ei()  # abs tol 5e-16
[1.89511781635594 +/- 4.94e-15]
>>> from sage.all import *
>>> RBF(Integer(1)).Ei()  # abs tol 5e-16
[1.89511781635594 +/- 4.94e-15]
Li()[source]#

Offset logarithmic integral

EXAMPLES:

sage: RBF(3).Li()  # abs tol 1e-15
[1.11842481454970 +/- 7.61e-15]
>>> from sage.all import *
>>> RBF(Integer(3)).Li()  # abs tol 1e-15
[1.11842481454970 +/- 7.61e-15]
Shi()[source]#

Hyperbolic sine integral

EXAMPLES:

sage: RBF(1).Shi()
[1.05725087537573 +/- 2.77e-15]
>>> from sage.all import *
>>> RBF(Integer(1)).Shi()
[1.05725087537573 +/- 2.77e-15]
Si()[source]#

Sine integral

EXAMPLES:

sage: RBF(1).Si() # abs tol 1e-15
[0.946083070367183 +/- 9.22e-16]
>>> from sage.all import *
>>> RBF(Integer(1)).Si() # abs tol 1e-15
[0.946083070367183 +/- 9.22e-16]
above_abs()[source]#

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

OUTPUT:

A ball with zero radius

EXAMPLES:

sage: b = RealBallField(8)(1/3).above_abs()
sage: b
[0.33 +/- ...e-3]
sage: b.is_exact()
True
sage: QQ(b)
171/512
>>> from sage.all import *
>>> b = RealBallField(Integer(8))(Integer(1)/Integer(3)).above_abs()
>>> b
[0.33 +/- ...e-3]
>>> b.is_exact()
True
>>> QQ(b)
171/512

See also

below_abs()

accuracy()[source]#

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

The accuracy is defined as the difference between the position of the top bit in the midpoint and the top bit in the radius, minus one. The result is clamped between plus/minus maximal_accuracy().

EXAMPLES:

sage: RBF(pi).accuracy()                                                    # needs sage.symbolic
52
sage: RBF(1).accuracy() == RBF.maximal_accuracy()
True
sage: RBF(NaN).accuracy() == -RBF.maximal_accuracy()                        # needs sage.symbolic
True
>>> from sage.all import *
>>> RBF(pi).accuracy()                                                    # needs sage.symbolic
52
>>> RBF(Integer(1)).accuracy() == RBF.maximal_accuracy()
True
>>> RBF(NaN).accuracy() == -RBF.maximal_accuracy()                        # needs sage.symbolic
True
add_error(ampl)[source]#

Increase the radius of this ball by (an upper bound on) ampl.

If ampl is negative, the radius is unchanged.

INPUT:

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

OUTPUT:

A new real ball.

EXAMPLES:

sage: err = RBF(10^-16)
sage: RBF(1).add_error(err)
[1.000000000000000 +/- ...e-16]
>>> from sage.all import *
>>> err = RBF(Integer(10)**-Integer(16))
>>> RBF(Integer(1)).add_error(err)
[1.000000000000000 +/- ...e-16]
agm(other)[source]#

Return the arithmetic-geometric mean of self and other.

EXAMPLES:

sage: RBF(1).agm(1)
1.000000000000000
sage: RBF(sqrt(2)).agm(1)^(-1)                                              # needs sage.symbolic
[0.8346268416740...]
>>> from sage.all import *
>>> RBF(Integer(1)).agm(Integer(1))
1.000000000000000
>>> RBF(sqrt(Integer(2))).agm(Integer(1))**(-Integer(1))                                              # needs sage.symbolic
[0.8346268416740...]
arccos()[source]#

Return the arccosine of this ball.

EXAMPLES:

sage: RBF(1).arccos()
0
sage: RBF(1, rad=.125r).arccos()
nan
>>> from sage.all import *
>>> RBF(Integer(1)).arccos()
0
>>> RBF(Integer(1), rad=.125).arccos()
nan
arccosh()[source]#

Return the inverse hyperbolic cosine of this ball.

EXAMPLES:

sage: RBF(2).arccosh()
[1.316957896924817 +/- ...e-16]
sage: RBF(1).arccosh()
0
sage: RBF(0).arccosh()
nan
>>> from sage.all import *
>>> RBF(Integer(2)).arccosh()
[1.316957896924817 +/- ...e-16]
>>> RBF(Integer(1)).arccosh()
0
>>> RBF(Integer(0)).arccosh()
nan
arcsin()[source]#

Return the arcsine of this ball.

EXAMPLES:

sage: RBF(1).arcsin()
[1.570796326794897 +/- ...e-16]
sage: RBF(1, rad=.125r).arcsin()
nan
>>> from sage.all import *
>>> RBF(Integer(1)).arcsin()
[1.570796326794897 +/- ...e-16]
>>> RBF(Integer(1), rad=.125).arcsin()
nan
arcsinh()[source]#

Return the inverse hyperbolic sine of this ball.

EXAMPLES:

sage: RBF(1).arcsinh()
[0.881373587019543 +/- ...e-16]
sage: RBF(0).arcsinh()
0
>>> from sage.all import *
>>> RBF(Integer(1)).arcsinh()
[0.881373587019543 +/- ...e-16]
>>> RBF(Integer(0)).arcsinh()
0
arctan()[source]#

Return the arctangent of this ball.

EXAMPLES:

sage: RBF(1).arctan()
[0.7853981633974483 +/- ...e-17]
>>> from sage.all import *
>>> RBF(Integer(1)).arctan()
[0.7853981633974483 +/- ...e-17]
arctanh()[source]#

Return the inverse hyperbolic tangent of this ball.

EXAMPLES:

sage: RBF(0).arctanh()
0
sage: RBF(1/2).arctanh()
[0.549306144334055 +/- ...e-16]
sage: RBF(1).arctanh()
nan
>>> from sage.all import *
>>> RBF(Integer(0)).arctanh()
0
>>> RBF(Integer(1)/Integer(2)).arctanh()
[0.549306144334055 +/- ...e-16]
>>> RBF(Integer(1)).arctanh()
nan
below_abs(test_zero=False)[source]#

Return a lower bound for the absolute value of this 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: RealBallField(8)(1/3).below_abs()
[0.33 +/- ...e-5]
sage: b = RealBallField(8)(1/3).below_abs()
sage: b
[0.33 +/- ...e-5]
sage: b.is_exact()
True
sage: QQ(b)
169/512

sage: RBF(0).below_abs()
0
sage: RBF(0).below_abs(test_zero=True)
Traceback (most recent call last):
...
ValueError: ball contains zero
>>> from sage.all import *
>>> RealBallField(Integer(8))(Integer(1)/Integer(3)).below_abs()
[0.33 +/- ...e-5]
>>> b = RealBallField(Integer(8))(Integer(1)/Integer(3)).below_abs()
>>> b
[0.33 +/- ...e-5]
>>> b.is_exact()
True
>>> QQ(b)
169/512

>>> RBF(Integer(0)).below_abs()
0
>>> RBF(Integer(0)).below_abs(test_zero=True)
Traceback (most recent call last):
...
ValueError: ball contains zero

See also

above_abs()

beta(a, z=1)[source]#

(Incomplete) beta function

INPUT:

  • a, z (optional) – real balls

OUTPUT:

The lower incomplete beta function \(B(self, a, z)\).

With the default value of z, the complete beta function \(B(self, a)\).

EXAMPLES:

sage: RBF(sin(3)).beta(RBF(2/3).sqrt())  # abs tol 1e-13                    # needs sage.symbolic
[7.407661629415 +/- 1.07e-13]
sage: RealBallField(100)(7/2).beta(1)  # abs tol 1e-30
[0.28571428571428571428571428571 +/- 5.23e-30]
sage: RealBallField(100)(7/2).beta(1, 1/2)
[0.025253813613805268728601584361 +/- 2.53e-31]
>>> from sage.all import *
>>> RBF(sin(Integer(3))).beta(RBF(Integer(2)/Integer(3)).sqrt())  # abs tol 1e-13                    # needs sage.symbolic
[7.407661629415 +/- 1.07e-13]
>>> RealBallField(Integer(100))(Integer(7)/Integer(2)).beta(Integer(1))  # abs tol 1e-30
[0.28571428571428571428571428571 +/- 5.23e-30]
>>> RealBallField(Integer(100))(Integer(7)/Integer(2)).beta(Integer(1), Integer(1)/Integer(2))
[0.025253813613805268728601584361 +/- 2.53e-31]

Todo

At the moment RBF(beta(a,b)) does not work, one needs RBF(a).beta(b) for this to work. See Issue #32851 and Issue #24641.

ceil()[source]#

Return the ceil of this ball.

EXAMPLES:

sage: RBF(1000+1/3, rad=1.r).ceil()
[1.00e+3 +/- 2.01]
>>> from sage.all import *
>>> RBF(Integer(1000)+Integer(1)/Integer(3), rad=1.).ceil()
[1.00e+3 +/- 2.01]
center()[source]#

Return the center of this ball.

EXAMPLES:

sage: RealBallField(16)(1/3).mid()
0.3333
sage: RealBallField(16)(1/3).mid().parent()
Real Field with 16 bits of precision
sage: RealBallField(16)(RBF(1/3)).mid().parent()
Real Field with 53 bits of precision
sage: RBF('inf').mid()
+infinity
>>> from sage.all import *
>>> RealBallField(Integer(16))(Integer(1)/Integer(3)).mid()
0.3333
>>> RealBallField(Integer(16))(Integer(1)/Integer(3)).mid().parent()
Real Field with 16 bits of precision
>>> RealBallField(Integer(16))(RBF(Integer(1)/Integer(3))).mid().parent()
Real Field with 53 bits of precision
>>> RBF('inf').mid()
+infinity
sage: b = RBF(2)^(2^1000)
sage: b.mid()
+infinity
>>> from sage.all import *
>>> b = RBF(Integer(2))**(Integer(2)**Integer(1000))
>>> b.mid()
+infinity

See also

rad(), squash()

chebyshev_T(n)[source]#

Evaluate the Chebyshev polynomial of the first kind T_n at this ball.

EXAMPLES:

sage: # needs sage.symbolic
sage: RBF(pi).chebyshev_T(0)
1.000000000000000
sage: RBF(pi).chebyshev_T(1)
[3.141592653589793 +/- ...e-16]
sage: RBF(pi).chebyshev_T(10**20)
Traceback (most recent call last):
...
ValueError: index too large
sage: RBF(pi).chebyshev_T(-1)
Traceback (most recent call last):
...
ValueError: expected a nonnegative index
>>> from sage.all import *
>>> # needs sage.symbolic
>>> RBF(pi).chebyshev_T(Integer(0))
1.000000000000000
>>> RBF(pi).chebyshev_T(Integer(1))
[3.141592653589793 +/- ...e-16]
>>> RBF(pi).chebyshev_T(Integer(10)**Integer(20))
Traceback (most recent call last):
...
ValueError: index too large
>>> RBF(pi).chebyshev_T(-Integer(1))
Traceback (most recent call last):
...
ValueError: expected a nonnegative index
chebyshev_U(n)[source]#

Evaluate the Chebyshev polynomial of the second kind U_n at this ball.

EXAMPLES:

sage: # needs sage.symbolic
sage: RBF(pi).chebyshev_U(0)
1.000000000000000
sage: RBF(pi).chebyshev_U(1)
[6.283185307179586 +/- ...e-16]
sage: RBF(pi).chebyshev_U(10**20)
Traceback (most recent call last):
...
ValueError: index too large
sage: RBF(pi).chebyshev_U(-1)
Traceback (most recent call last):
...
ValueError: expected a nonnegative index
>>> from sage.all import *
>>> # needs sage.symbolic
>>> RBF(pi).chebyshev_U(Integer(0))
1.000000000000000
>>> RBF(pi).chebyshev_U(Integer(1))
[6.283185307179586 +/- ...e-16]
>>> RBF(pi).chebyshev_U(Integer(10)**Integer(20))
Traceback (most recent call last):
...
ValueError: index too large
>>> RBF(pi).chebyshev_U(-Integer(1))
Traceback (most recent call last):
...
ValueError: expected a nonnegative index
contains_exact(other)[source]#

Return True iff the given number (or ball) other is contained in the interval represented by self.

If self contains NaN, this function always returns True (as it could represent anything, and in particular could represent all the points included in other). If other contains NaN and self does not, it always returns False.

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

EXAMPLES:

sage: b = RBF(1)
sage: b.contains_exact(1)
True
sage: b.contains_exact(QQ(1))
True
sage: b.contains_exact(1.)
True
sage: b.contains_exact(b)
True
>>> from sage.all import *
>>> b = RBF(Integer(1))
>>> b.contains_exact(Integer(1))
True
>>> b.contains_exact(QQ(Integer(1)))
True
>>> b.contains_exact(RealNumber('1.'))
True
>>> b.contains_exact(b)
True
sage: RBF(1/3).contains_exact(1/3)
True
sage: RBF(sqrt(2)).contains_exact(sqrt(2))                                  # needs sage.symbolic
Traceback (most recent call last):
...
TypeError: unsupported type: <class 'sage.symbolic.expression.Expression'>
>>> from sage.all import *
>>> RBF(Integer(1)/Integer(3)).contains_exact(Integer(1)/Integer(3))
True
>>> RBF(sqrt(Integer(2))).contains_exact(sqrt(Integer(2)))                                  # needs sage.symbolic
Traceback (most recent call last):
...
TypeError: unsupported type: <class 'sage.symbolic.expression.Expression'>
contains_integer()[source]#

Return True iff this ball contains any integer.

EXAMPLES:

sage: RBF(3.1, 0.1).contains_integer()
True
sage: RBF(3.1, 0.05).contains_integer()
False
>>> from sage.all import *
>>> RBF(RealNumber('3.1'), RealNumber('0.1')).contains_integer()
True
>>> RBF(RealNumber('3.1'), RealNumber('0.05')).contains_integer()
False
contains_zero()[source]#

Return True iff this ball contains zero.

EXAMPLES:

sage: RBF(0).contains_zero()
True
sage: RBF(RIF(-1, 1)).contains_zero()
True
sage: RBF(1/3).contains_zero()
False
>>> from sage.all import *
>>> RBF(Integer(0)).contains_zero()
True
>>> RBF(RIF(-Integer(1), Integer(1))).contains_zero()
True
>>> RBF(Integer(1)/Integer(3)).contains_zero()
False
cos()[source]#

Return the cosine of this ball.

EXAMPLES:

sage: RBF(pi).cos()                                                         # needs sage.symbolic
[-1.00000000000000 +/- ...e-16]
>>> from sage.all import *
>>> RBF(pi).cos()                                                         # needs sage.symbolic
[-1.00000000000000 +/- ...e-16]

See also

cospi()

cos_integral()[source]#

Cosine integral

EXAMPLES:

sage: RBF(1).Ci()  # abs tol 5e-16
[0.337403922900968 +/- 3.25e-16]
>>> from sage.all import *
>>> RBF(Integer(1)).Ci()  # abs tol 5e-16
[0.337403922900968 +/- 3.25e-16]
cosh()[source]#

Return the hyperbolic cosine of this ball.

EXAMPLES:

sage: RBF(1).cosh()
[1.543080634815244 +/- ...e-16]
>>> from sage.all import *
>>> RBF(Integer(1)).cosh()
[1.543080634815244 +/- ...e-16]
cosh_integral()[source]#

Hyperbolic cosine integral

EXAMPLES:

sage: RBF(1).Chi()  # abs tol 1e-17
[0.837866940980208 +/- 4.72e-16]
>>> from sage.all import *
>>> RBF(Integer(1)).Chi()  # abs tol 1e-17
[0.837866940980208 +/- 4.72e-16]
cot()[source]#

Return the cotangent of this ball.

EXAMPLES:

sage: RBF(1).cot()
[0.642092615934331 +/- ...e-16]
sage: RBF(pi).cot()                                                         # needs sage.symbolic
nan
>>> from sage.all import *
>>> RBF(Integer(1)).cot()
[0.642092615934331 +/- ...e-16]
>>> RBF(pi).cot()                                                         # needs sage.symbolic
nan
coth()[source]#

Return the hyperbolic cotangent of this ball.

EXAMPLES:

sage: RBF(1).coth()
[1.313035285499331 +/- ...e-16]
sage: RBF(0).coth()
nan
>>> from sage.all import *
>>> RBF(Integer(1)).coth()
[1.313035285499331 +/- ...e-16]
>>> RBF(Integer(0)).coth()
nan
csc()[source]#

Return the cosecant of this ball.

EXAMPLES:

sage: RBF(1).csc()
[1.188395105778121 +/- ...e-16]
>>> from sage.all import *
>>> RBF(Integer(1)).csc()
[1.188395105778121 +/- ...e-16]
csch()[source]#

Return the hyperbolic cosecant of this ball.

EXAMPLES:

sage: RBF(1).csch()
[0.850918128239321 +/- ...e-16]
>>> from sage.all import *
>>> RBF(Integer(1)).csch()
[0.850918128239321 +/- ...e-16]
diameter()[source]#

Return the diameter of this ball.

EXAMPLES:

sage: RBF(1/3).diameter()
1.1102230e-16
sage: RBF(1/3).diameter().parent()
Real Field with 30 bits of precision
sage: RBF(RIF(1.02, 1.04)).diameter()
0.020000000
>>> from sage.all import *
>>> RBF(Integer(1)/Integer(3)).diameter()
1.1102230e-16
>>> RBF(Integer(1)/Integer(3)).diameter().parent()
Real Field with 30 bits of precision
>>> RBF(RIF(RealNumber('1.02'), RealNumber('1.04'))).diameter()
0.020000000
endpoints(rnd=None)[source]#

Return the endpoints of this ball, rounded outwards.

INPUT:

OUTPUT:

A pair of real numbers.

EXAMPLES:

sage: RBF(-1/3).endpoints()
(-0.333333333333334, -0.333333333333333)
>>> from sage.all import *
>>> RBF(-Integer(1)/Integer(3)).endpoints()
(-0.333333333333334, -0.333333333333333)

See also

lower(), upper()

erf()[source]#

Error function.

EXAMPLES:

sage: RBF(1/2).erf() # abs tol 1e-16
[0.520499877813047 +/- 6.10e-16]
>>> from sage.all import *
>>> RBF(Integer(1)/Integer(2)).erf() # abs tol 1e-16
[0.520499877813047 +/- 6.10e-16]
erfi()[source]#

Imaginary error function

EXAMPLES:

sage: RBF(1/2).erfi()
[0.614952094696511 +/- 2.22e-16]
>>> from sage.all import *
>>> RBF(Integer(1)/Integer(2)).erfi()
[0.614952094696511 +/- 2.22e-16]
exp()[source]#

Return the exponential of this ball.

EXAMPLES:

sage: RBF(1).exp()
[2.718281828459045 +/- ...e-16]
>>> from sage.all import *
>>> RBF(Integer(1)).exp()
[2.718281828459045 +/- ...e-16]
expm1()[source]#

Return exp(self) - 1, computed accurately when self is close to zero.

EXAMPLES:

sage: eps = RBF(1e-30)
sage: exp(eps) - 1
[+/- ...e-30]
sage: eps.expm1()
[1.000000000000000e-30 +/- ...e-47]
>>> from sage.all import *
>>> eps = RBF(RealNumber('1e-30'))
>>> exp(eps) - Integer(1)
[+/- ...e-30]
>>> eps.expm1()
[1.000000000000000e-30 +/- ...e-47]
floor()[source]#

Return the floor of this ball.

EXAMPLES:

sage: RBF(1000+1/3, rad=1.r).floor()
[1.00e+3 +/- 1.01]
>>> from sage.all import *
>>> RBF(Integer(1000)+Integer(1)/Integer(3), rad=1.).floor()
[1.00e+3 +/- 1.01]
gamma(a=None)[source]#

Image of this ball by the (upper incomplete) Euler Gamma function

For \(a\) real, return the upper incomplete Gamma function \(\Gamma(self,a)\).

For integer and rational arguments, gamma() may be faster.

EXAMPLES:

sage: RBF(1/2).gamma()
[1.772453850905516 +/- ...e-16]
sage: RBF(gamma(3/2, RBF(2).sqrt()))  # abs tol 2e-17
[0.37118875695353 +/- 3.00e-15]
sage: RBF(3/2).gamma_inc(RBF(2).sqrt())  # abs tol 2e-17
[0.37118875695353 +/- 3.00e-15]
>>> from sage.all import *
>>> RBF(Integer(1)/Integer(2)).gamma()
[1.772453850905516 +/- ...e-16]
>>> RBF(gamma(Integer(3)/Integer(2), RBF(Integer(2)).sqrt()))  # abs tol 2e-17
[0.37118875695353 +/- 3.00e-15]
>>> RBF(Integer(3)/Integer(2)).gamma_inc(RBF(Integer(2)).sqrt())  # abs tol 2e-17
[0.37118875695353 +/- 3.00e-15]

See also

gamma()

gamma_inc(a=None)[source]#

Image of this ball by the (upper incomplete) Euler Gamma function

For \(a\) real, return the upper incomplete Gamma function \(\Gamma(self,a)\).

For integer and rational arguments, gamma() may be faster.

EXAMPLES:

sage: RBF(1/2).gamma()
[1.772453850905516 +/- ...e-16]
sage: RBF(gamma(3/2, RBF(2).sqrt()))  # abs tol 2e-17
[0.37118875695353 +/- 3.00e-15]
sage: RBF(3/2).gamma_inc(RBF(2).sqrt())  # abs tol 2e-17
[0.37118875695353 +/- 3.00e-15]
>>> from sage.all import *
>>> RBF(Integer(1)/Integer(2)).gamma()
[1.772453850905516 +/- ...e-16]
>>> RBF(gamma(Integer(3)/Integer(2), RBF(Integer(2)).sqrt()))  # abs tol 2e-17
[0.37118875695353 +/- 3.00e-15]
>>> RBF(Integer(3)/Integer(2)).gamma_inc(RBF(Integer(2)).sqrt())  # abs tol 2e-17
[0.37118875695353 +/- 3.00e-15]

See also

gamma()

gamma_inc_lower(a)[source]#

Image of this ball by the lower incomplete Euler Gamma function

For \(a\) real, return the lower incomplete Gamma function of \(\Gamma(self,a)\).

EXAMPLES:

sage: RBF(gamma_inc_lower(1/2, RBF(2).sqrt()))
[1.608308637729248 +/- 8.14e-16]
sage: RealBallField(100)(7/2).gamma_inc_lower(5)
[2.6966551541863035516887949614 +/- 8.91e-29]
>>> from sage.all import *
>>> RBF(gamma_inc_lower(Integer(1)/Integer(2), RBF(Integer(2)).sqrt()))
[1.608308637729248 +/- 8.14e-16]
>>> RealBallField(Integer(100))(Integer(7)/Integer(2)).gamma_inc_lower(Integer(5))
[2.6966551541863035516887949614 +/- 8.91e-29]
identical(other)[source]#

Return True iff self and other are equal as balls, i.e. have both the same midpoint and radius.

Note that this is not the same thing as testing whether both self and other certainly represent the same 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 contains(), depending on the circumstance.

EXAMPLES:

sage: RBF(1).identical(RBF(3)-RBF(2))
True
sage: RBF(1, rad=0.25r).identical(RBF(1, rad=0.25r))
True
sage: RBF(1).identical(RBF(1, rad=0.25r))
False
>>> from sage.all import *
>>> RBF(Integer(1)).identical(RBF(Integer(3))-RBF(Integer(2)))
True
>>> RBF(Integer(1), rad=0.25).identical(RBF(Integer(1), rad=0.25))
True
>>> RBF(Integer(1)).identical(RBF(Integer(1), rad=0.25))
False
imag()[source]#

Return the imaginary part of this ball.

EXAMPLES:

sage: RBF(1/3).imag()
0
>>> from sage.all import *
>>> RBF(Integer(1)/Integer(3)).imag()
0
is_NaN()[source]#

Return True if this ball is not-a-number.

EXAMPLES:

sage: RBF(NaN).is_NaN()                                                     # needs sage.symbolic
True
sage: RBF(-5).gamma().is_NaN()
True
sage: RBF(infinity).is_NaN()
False
sage: RBF(42, rad=1.r).is_NaN()
False
>>> from sage.all import *
>>> RBF(NaN).is_NaN()                                                     # needs sage.symbolic
True
>>> RBF(-Integer(5)).gamma().is_NaN()
True
>>> RBF(infinity).is_NaN()
False
>>> RBF(Integer(42), rad=1.).is_NaN()
False
is_exact()[source]#

Return True iff the radius of this ball is zero.

EXAMPLES:

sage: RBF = RealBallField()
sage: RBF(1).is_exact()
True
sage: RBF(RIF(0.1, 0.2)).is_exact()
False
>>> from sage.all import *
>>> RBF = RealBallField()
>>> RBF(Integer(1)).is_exact()
True
>>> RBF(RIF(RealNumber('0.1'), RealNumber('0.2'))).is_exact()
False
is_finite()[source]#

Return True iff the midpoint and radius of this ball are both finite floating-point numbers, i.e. not infinities or NaN.

EXAMPLES:

sage: (RBF(2)^(2^1000)).is_finite()
True
sage: RBF(oo).is_finite()
False
>>> from sage.all import *
>>> (RBF(Integer(2))**(Integer(2)**Integer(1000))).is_finite()
True
>>> RBF(oo).is_finite()
False
is_infinity()[source]#

Return True if this ball contains or may represent a point at infinity.

This is the exact negation of is_finite(), used in comparisons with Sage symbolic infinities.

Warning

Contrary to the usual convention, a return value of True does not imply that all points of the ball satisfy the predicate. This is due to the way comparisons with symbolic infinities work in sage.

EXAMPLES:

sage: RBF(infinity).is_infinity()
True
sage: RBF(-infinity).is_infinity()
True
sage: RBF(NaN).is_infinity()                                                # needs sage.symbolic
True
sage: (~RBF(0)).is_infinity()
True
sage: RBF(42, rad=1.r).is_infinity()
False
>>> from sage.all import *
>>> RBF(infinity).is_infinity()
True
>>> RBF(-infinity).is_infinity()
True
>>> RBF(NaN).is_infinity()                                                # needs sage.symbolic
True
>>> (~RBF(Integer(0))).is_infinity()
True
>>> RBF(Integer(42), rad=1.).is_infinity()
False
is_negative_infinity()[source]#

Return True if this ball is the point -∞.

EXAMPLES:

sage: RBF(-infinity).is_negative_infinity()
True
>>> from sage.all import *
>>> RBF(-infinity).is_negative_infinity()
True
is_nonzero()[source]#

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: RBF = RealBallField()
sage: RBF(pi).is_nonzero()                                                  # needs sage.symbolic
True
sage: RBF(RIF(-0.5, 0.5)).is_nonzero()
False
>>> from sage.all import *
>>> RBF = RealBallField()
>>> RBF(pi).is_nonzero()                                                  # needs sage.symbolic
True
>>> RBF(RIF(-RealNumber('0.5'), RealNumber('0.5'))).is_nonzero()
False

See also

is_zero()

is_positive_infinity()[source]#

Return True if this ball is the point +∞.

EXAMPLES:

sage: RBF(infinity).is_positive_infinity()
True
>>> from sage.all import *
>>> RBF(infinity).is_positive_infinity()
True
is_zero()[source]#

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

EXAMPLES:

sage: RBF = RealBallField()
sage: RBF(0).is_zero()
True
sage: RBF(RIF(-0.5, 0.5)).is_zero()
False
>>> from sage.all import *
>>> RBF = RealBallField()
>>> RBF(Integer(0)).is_zero()
True
>>> RBF(RIF(-RealNumber('0.5'), RealNumber('0.5'))).is_zero()
False

See also

is_nonzero()

lambert_w()[source]#

Return the image of this ball by the Lambert W function.

EXAMPLES:

sage: RBF(1).lambert_w()
[0.5671432904097...]
>>> from sage.all import *
>>> RBF(Integer(1)).lambert_w()
[0.5671432904097...]
li()[source]#

Logarithmic integral

EXAMPLES:

sage: RBF(3).li()  # abs tol 1e-15
[2.16358859466719 +/- 4.72e-15]
>>> from sage.all import *
>>> RBF(Integer(3)).li()  # abs tol 1e-15
[2.16358859466719 +/- 4.72e-15]
log(base=None)[source]#

Return the logarithm of this ball.

INPUT:

  • base (optional, positive real ball or number) – if None, return the natural logarithm ln(self), otherwise, return the general logarithm ln(self)/ln(base)

EXAMPLES:

sage: RBF(3).log()
[1.098612288668110 +/- ...e-16]
sage: RBF(3).log(2)
[1.58496250072116 +/- ...e-15]
sage: log(RBF(5), 2)
[2.32192809488736 +/- ...e-15]

sage: RBF(-1/3).log()
nan
sage: RBF(3).log(-1)
nan
sage: RBF(2).log(0)
nan
>>> from sage.all import *
>>> RBF(Integer(3)).log()
[1.098612288668110 +/- ...e-16]
>>> RBF(Integer(3)).log(Integer(2))
[1.58496250072116 +/- ...e-15]
>>> log(RBF(Integer(5)), Integer(2))
[2.32192809488736 +/- ...e-15]

>>> RBF(-Integer(1)/Integer(3)).log()
nan
>>> RBF(Integer(3)).log(-Integer(1))
nan
>>> RBF(Integer(2)).log(Integer(0))
nan
log1p()[source]#

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

EXAMPLES:

sage: eps = RBF(1e-30)
sage: (1 + eps).log()
[+/- ...e-16]
sage: eps.log1p()
[1.00000000000000e-30 +/- ...e-46]
>>> from sage.all import *
>>> eps = RBF(RealNumber('1e-30'))
>>> (Integer(1) + eps).log()
[+/- ...e-16]
>>> eps.log1p()
[1.00000000000000e-30 +/- ...e-46]
log_gamma()[source]#

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

The complex branch structure is assumed, so if self <= 0, the result is an indeterminate interval.

EXAMPLES:

sage: RBF(1/2).log_gamma()
[0.572364942924700 +/- ...e-16]
>>> from sage.all import *
>>> RBF(Integer(1)/Integer(2)).log_gamma()
[0.572364942924700 +/- ...e-16]
log_integral()[source]#

Logarithmic integral

EXAMPLES:

sage: RBF(3).li()  # abs tol 1e-15
[2.16358859466719 +/- 4.72e-15]
>>> from sage.all import *
>>> RBF(Integer(3)).li()  # abs tol 1e-15
[2.16358859466719 +/- 4.72e-15]
log_integral_offset()[source]#

Offset logarithmic integral

EXAMPLES:

sage: RBF(3).Li()  # abs tol 1e-15
[1.11842481454970 +/- 7.61e-15]
>>> from sage.all import *
>>> RBF(Integer(3)).Li()  # abs tol 1e-15
[1.11842481454970 +/- 7.61e-15]
lower(rnd=None)[source]#

Return the right endpoint of this ball, rounded downwards.

INPUT:

OUTPUT:

A real number.

EXAMPLES:

sage: RBF(-1/3).lower()
-0.333333333333334
sage: RBF(-1/3).lower().parent()
Real Field with 53 bits of precision and rounding RNDD
>>> from sage.all import *
>>> RBF(-Integer(1)/Integer(3)).lower()
-0.333333333333334
>>> RBF(-Integer(1)/Integer(3)).lower().parent()
Real Field with 53 bits of precision and rounding RNDD

See also

upper(), endpoints()

max(*others)[source]#

Return a ball containing the maximum of this ball and the remaining arguments.

EXAMPLES:

sage: RBF(-1, rad=.5).max(0)
0

sage: RBF(0, rad=2.).max(RBF(0, rad=1.)).endpoints()
(-1.00000000465662, 2.00000000651926)

sage: RBF(-infinity).max(-3, 1/3)
[0.3333333333333333 +/- ...e-17]

sage: RBF('nan').max(0)
nan
>>> from sage.all import *
>>> RBF(-Integer(1), rad=RealNumber('.5')).max(Integer(0))
0

>>> RBF(Integer(0), rad=RealNumber('2.')).max(RBF(Integer(0), rad=RealNumber('1.'))).endpoints()
(-1.00000000465662, 2.00000000651926)

>>> RBF(-infinity).max(-Integer(3), Integer(1)/Integer(3))
[0.3333333333333333 +/- ...e-17]

>>> RBF('nan').max(Integer(0))
nan

See also

min()

mid()[source]#

Return the center of this ball.

EXAMPLES:

sage: RealBallField(16)(1/3).mid()
0.3333
sage: RealBallField(16)(1/3).mid().parent()
Real Field with 16 bits of precision
sage: RealBallField(16)(RBF(1/3)).mid().parent()
Real Field with 53 bits of precision
sage: RBF('inf').mid()
+infinity
>>> from sage.all import *
>>> RealBallField(Integer(16))(Integer(1)/Integer(3)).mid()
0.3333
>>> RealBallField(Integer(16))(Integer(1)/Integer(3)).mid().parent()
Real Field with 16 bits of precision
>>> RealBallField(Integer(16))(RBF(Integer(1)/Integer(3))).mid().parent()
Real Field with 53 bits of precision
>>> RBF('inf').mid()
+infinity
sage: b = RBF(2)^(2^1000)
sage: b.mid()
+infinity
>>> from sage.all import *
>>> b = RBF(Integer(2))**(Integer(2)**Integer(1000))
>>> b.mid()
+infinity

See also

rad(), squash()

min(*others)[source]#

Return a ball containing the minimum of this ball and the remaining arguments.

EXAMPLES:

sage: RBF(1, rad=.5).min(0)
0

sage: RBF(0, rad=2.).min(RBF(0, rad=1.)).endpoints()
(-2.00000000651926, 1.00000000465662)

sage: RBF(infinity).min(3, 1/3)
[0.3333333333333333 +/- ...e-17]

sage: RBF('nan').min(0)
nan
>>> from sage.all import *
>>> RBF(Integer(1), rad=RealNumber('.5')).min(Integer(0))
0

>>> RBF(Integer(0), rad=RealNumber('2.')).min(RBF(Integer(0), rad=RealNumber('1.'))).endpoints()
(-2.00000000651926, 1.00000000465662)

>>> RBF(infinity).min(Integer(3), Integer(1)/Integer(3))
[0.3333333333333333 +/- ...e-17]

>>> RBF('nan').min(Integer(0))
nan

See also

max()

nbits()[source]#

Return the minimum precision sufficient to represent this ball exactly.

In other words, return the number of bits needed to represent the absolute value of the mantissa of the midpoint of this ball. The result is 0 if the midpoint is a special value.

EXAMPLES:

sage: RBF(1/3).nbits()
53
sage: RBF(1023, .1).nbits()
10
sage: RBF(1024, .1).nbits()
1
sage: RBF(0).nbits()
0
sage: RBF(infinity).nbits()
0
>>> from sage.all import *
>>> RBF(Integer(1)/Integer(3)).nbits()
53
>>> RBF(Integer(1023), RealNumber('.1')).nbits()
10
>>> RBF(Integer(1024), RealNumber('.1')).nbits()
1
>>> RBF(Integer(0)).nbits()
0
>>> RBF(infinity).nbits()
0
overlaps(other)[source]#

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

If either self or other contains NaN, this method always returns nonzero (as a NaN could be anything, it could in particular contain any number that is included in the other operand).

EXAMPLES:

sage: RBF(pi).overlaps(RBF(pi) + 2**(-100))                                 # needs sage.symbolic
True
sage: RBF(pi).overlaps(RBF(3))                                              # needs sage.symbolic
False
>>> from sage.all import *
>>> RBF(pi).overlaps(RBF(pi) + Integer(2)**(-Integer(100)))                                 # needs sage.symbolic
True
>>> RBF(pi).overlaps(RBF(Integer(3)))                                              # needs sage.symbolic
False
polylog(s)[source]#

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

EXAMPLES:

sage: polylog(0, -1)                                                        # needs sage.symbolic
-1/2
sage: RBF(-1).polylog(0)
[-0.50000000000000 +/- ...e-16]
sage: polylog(1, 1/2)                                                       # needs sage.symbolic
-log(1/2)
sage: RBF(1/2).polylog(1)
[0.69314718055995 +/- ...e-15]
sage: RBF(1/3).polylog(1/2)
[0.44210883528067 +/- 6.7...e-15]
sage: RBF(1/3).polylog(RLF(pi))                                             # needs sage.symbolic
[0.34728895057225 +/- ...e-15]
>>> from sage.all import *
>>> polylog(Integer(0), -Integer(1))                                                        # needs sage.symbolic
-1/2
>>> RBF(-Integer(1)).polylog(Integer(0))
[-0.50000000000000 +/- ...e-16]
>>> polylog(Integer(1), Integer(1)/Integer(2))                                                       # needs sage.symbolic
-log(1/2)
>>> RBF(Integer(1)/Integer(2)).polylog(Integer(1))
[0.69314718055995 +/- ...e-15]
>>> RBF(Integer(1)/Integer(3)).polylog(Integer(1)/Integer(2))
[0.44210883528067 +/- 6.7...e-15]
>>> RBF(Integer(1)/Integer(3)).polylog(RLF(pi))                                             # needs sage.symbolic
[0.34728895057225 +/- ...e-15]
psi()[source]#

Compute the digamma function with argument self.

EXAMPLES:

sage: RBF(1).psi() # abs tol 1e-15
[-0.5772156649015329 +/- 4.84e-17]
>>> from sage.all import *
>>> RBF(Integer(1)).psi() # abs tol 1e-15
[-0.5772156649015329 +/- 4.84e-17]
rad()[source]#

Return the radius of this ball.

EXAMPLES:

sage: RBF(1/3).rad()
5.5511151e-17
sage: RBF(1/3).rad().parent()
Real Field with 30 bits of precision
>>> from sage.all import *
>>> RBF(Integer(1)/Integer(3)).rad()
5.5511151e-17
>>> RBF(Integer(1)/Integer(3)).rad().parent()
Real Field with 30 bits of precision
rad_as_ball()[source]#

Return an exact ball with center equal to the radius of this ball.

EXAMPLES:

sage: rad = RBF(1/3).rad_as_ball()
sage: rad
[5.55111512e-17 +/- ...e-26]
sage: rad.is_exact()
True
sage: rad.parent()
Real ball field with 30 bits of precision
>>> from sage.all import *
>>> rad = RBF(Integer(1)/Integer(3)).rad_as_ball()
>>> rad
[5.55111512e-17 +/- ...e-26]
>>> rad.is_exact()
True
>>> rad.parent()
Real ball field with 30 bits of precision

See also

squash(), rad()

real()[source]#

Return the real part of this ball.

EXAMPLES:

sage: RBF(1/3).real()
[0.3333333333333333 +/- 7.04e-17]
>>> from sage.all import *
>>> RBF(Integer(1)/Integer(3)).real()
[0.3333333333333333 +/- 7.04e-17]
rgamma()[source]#

Return the image of this ball by the function 1/Γ, avoiding division by zero at the poles of the gamma function.

EXAMPLES:

sage: RBF(-1).rgamma()
0
sage: RBF(3).rgamma()
0.5000000000000000
>>> from sage.all import *
>>> RBF(-Integer(1)).rgamma()
0
>>> RBF(Integer(3)).rgamma()
0.5000000000000000
rising_factorial(n)[source]#

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 real \(n\), it is a quotient of gamma functions.

EXAMPLES:

sage: RBF(1).rising_factorial(5)
120.0000000000000
sage: RBF(1/2).rising_factorial(1/3) # abs tol 1e-14
[0.636849884317974 +/- 8.98e-16]
>>> from sage.all import *
>>> RBF(Integer(1)).rising_factorial(Integer(5))
120.0000000000000
>>> RBF(Integer(1)/Integer(2)).rising_factorial(Integer(1)/Integer(3)) # abs tol 1e-14
[0.636849884317974 +/- 8.98e-16]
round()[source]#

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

EXAMPLES:

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

sage: b = RBF(pi.n(100))                                                    # needs sage.symbolic
sage: b.mid()                                                               # needs sage.symbolic
3.141592653589793238462643383
>>> from sage.all import *
>>> b = RBF(pi.n(Integer(100)))                                                    # needs sage.symbolic
>>> b.mid()                                                               # needs sage.symbolic
3.141592653589793238462643383

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

sage: b.round().mid()                                                       # needs sage.symbolic
3.14159265358979
>>> from sage.all import *
>>> b.round().mid()                                                       # needs sage.symbolic
3.14159265358979

See also

trim()

rsqrt()[source]#

Return the reciprocal square root of self.

At high precision, this is faster than computing a square root.

EXAMPLES:

sage: RBF(2).rsqrt()
[0.707106781186547 +/- ...e-16]
sage: RBF(0).rsqrt()
nan
>>> from sage.all import *
>>> RBF(Integer(2)).rsqrt()
[0.707106781186547 +/- ...e-16]
>>> RBF(Integer(0)).rsqrt()
nan
sec()[source]#

Return the secant of this ball.

EXAMPLES:

sage: RBF(1).sec()
[1.850815717680925 +/- ...e-16]
>>> from sage.all import *
>>> RBF(Integer(1)).sec()
[1.850815717680925 +/- ...e-16]
sech()[source]#

Return the hyperbolic secant of this ball.

EXAMPLES:

sage: RBF(1).sech()
[0.648054273663885 +/- ...e-16]
>>> from sage.all import *
>>> RBF(Integer(1)).sech()
[0.648054273663885 +/- ...e-16]
sin()[source]#

Return the sine of this ball.

EXAMPLES:

sage: RBF(pi).sin()                                                         # needs sage.symbolic
[+/- ...e-16]
>>> from sage.all import *
>>> RBF(pi).sin()                                                         # needs sage.symbolic
[+/- ...e-16]

See also

sinpi()

sin_integral()[source]#

Sine integral

EXAMPLES:

sage: RBF(1).Si() # abs tol 1e-15
[0.946083070367183 +/- 9.22e-16]
>>> from sage.all import *
>>> RBF(Integer(1)).Si() # abs tol 1e-15
[0.946083070367183 +/- 9.22e-16]
sinh()[source]#

Return the hyperbolic sine of this ball.

EXAMPLES:

sage: RBF(1).sinh()
[1.175201193643801 +/- ...e-16]
>>> from sage.all import *
>>> RBF(Integer(1)).sinh()
[1.175201193643801 +/- ...e-16]
sinh_integral()[source]#

Hyperbolic sine integral

EXAMPLES:

sage: RBF(1).Shi()
[1.05725087537573 +/- 2.77e-15]
>>> from sage.all import *
>>> RBF(Integer(1)).Shi()
[1.05725087537573 +/- 2.77e-15]
sqrt()[source]#

Return the square root of this ball.

EXAMPLES:

sage: RBF(2).sqrt()
[1.414213562373095 +/- ...e-16]
sage: RBF(-1/3).sqrt()
nan
>>> from sage.all import *
>>> RBF(Integer(2)).sqrt()
[1.414213562373095 +/- ...e-16]
>>> RBF(-Integer(1)/Integer(3)).sqrt()
nan
sqrt1pm1()[source]#

Return \(\sqrt{1+\mathrm{self}}-1\), computed accurately when self is close to zero.

EXAMPLES:

sage: eps = RBF(10^(-20))
sage: (1 + eps).sqrt() - 1
[+/- ...e-16]
sage: eps.sqrt1pm1()
[5.00000000000000e-21 +/- ...e-36]
>>> from sage.all import *
>>> eps = RBF(Integer(10)**(-Integer(20)))
>>> (Integer(1) + eps).sqrt() - Integer(1)
[+/- ...e-16]
>>> eps.sqrt1pm1()
[5.00000000000000e-21 +/- ...e-36]
sqrtpos()[source]#

Return the square root of this ball, assuming that it represents a nonnegative number.

Any negative numbers in the input interval are discarded.

EXAMPLES:

sage: RBF(2).sqrtpos()
[1.414213562373095 +/- ...e-16]
sage: RBF(-1/3).sqrtpos()
0
sage: RBF(0, rad=2.r).sqrtpos()
[+/- 1.42]
>>> from sage.all import *
>>> RBF(Integer(2)).sqrtpos()
[1.414213562373095 +/- ...e-16]
>>> RBF(-Integer(1)/Integer(3)).sqrtpos()
0
>>> RBF(Integer(0), rad=2.).sqrtpos()
[+/- 1.42]
squash()[source]#

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

EXAMPLES:

sage: mid = RealBallField(16)(1/3).squash()
sage: mid
[0.3333 +/- ...e-5]
sage: mid.is_exact()
True
sage: mid.parent()
Real ball field with 16 bits of precision
>>> from sage.all import *
>>> mid = RealBallField(Integer(16))(Integer(1)/Integer(3)).squash()
>>> mid
[0.3333 +/- ...e-5]
>>> mid.is_exact()
True
>>> mid.parent()
Real ball field with 16 bits of precision

See also

mid(), rad_as_ball()

tan()[source]#

Return the tangent of this ball.

EXAMPLES:

sage: RBF(1).tan()
[1.557407724654902 +/- ...e-16]
sage: RBF(pi/2).tan()                                                       # needs sage.symbolic
nan
>>> from sage.all import *
>>> RBF(Integer(1)).tan()
[1.557407724654902 +/- ...e-16]
>>> RBF(pi/Integer(2)).tan()                                                       # needs sage.symbolic
nan
tanh()[source]#

Return the hyperbolic tangent of this ball.

EXAMPLES:

sage: RBF(1).tanh()
[0.761594155955765 +/- ...e-16]
>>> from sage.all import *
>>> RBF(Integer(1)).tanh()
[0.761594155955765 +/- ...e-16]
trim()[source]#

Return a trimmed copy of this ball.

Round self to a number of bits equal to the accuracy() of self (as indicated by its radius), plus a few guard bits. The resulting ball is guaranteed to contain self, but is more economical if self has less than full accuracy.

EXAMPLES:

sage: b = RBF(0.11111111111111, rad=.001)
sage: b.mid()
0.111111111111110
sage: b.trim().mid()
0.111111104488373
>>> from sage.all import *
>>> b = RBF(RealNumber('0.11111111111111'), rad=RealNumber('.001'))
>>> b.mid()
0.111111111111110
>>> b.trim().mid()
0.111111104488373

See also

round()

union(other)[source]#

Return a ball containing the convex hull of self and other.

EXAMPLES:

sage: RBF(0).union(1).endpoints()
(-9.31322574615479e-10, 1.00000000093133)
>>> from sage.all import *
>>> RBF(Integer(0)).union(Integer(1)).endpoints()
(-9.31322574615479e-10, 1.00000000093133)
upper(rnd=None)[source]#

Return the right endpoint of this ball, rounded upwards.

INPUT:

OUTPUT:

A real number.

EXAMPLES:

sage: RBF(-1/3).upper()
-0.333333333333333
sage: RBF(-1/3).upper().parent()
Real Field with 53 bits of precision and rounding RNDU
>>> from sage.all import *
>>> RBF(-Integer(1)/Integer(3)).upper()
-0.333333333333333
>>> RBF(-Integer(1)/Integer(3)).upper().parent()
Real Field with 53 bits of precision and rounding RNDU

See also

lower(), endpoints()

zeta(a=None)[source]#

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

For a = 1 (or a = None), this computes the Riemann zeta function.

Otherwise, it computes the Hurwitz zeta function.

Use RealBallField.zeta() to compute the Riemann zeta function of a small integer without first converting it to a real ball.

EXAMPLES:

sage: RBF(-1).zeta()
[-0.0833333333333333 +/- ...e-17]
sage: RBF(-1).zeta(1)
[-0.0833333333333333 +/- ...e-17]
sage: RBF(-1).zeta(2)
[-1.083333333333333 +/- ...e-16]
>>> from sage.all import *
>>> RBF(-Integer(1)).zeta()
[-0.0833333333333333 +/- ...e-17]
>>> RBF(-Integer(1)).zeta(Integer(1))
[-0.0833333333333333 +/- ...e-17]
>>> RBF(-Integer(1)).zeta(Integer(2))
[-1.083333333333333 +/- ...e-16]
zetaderiv(k)[source]#

Return the image of this ball by the k-th derivative of the Riemann zeta function.

For a more flexible interface, see the low-level method _zeta_series of polynomials with complex ball coefficients.

EXAMPLES:

sage: RBF(1/2).zetaderiv(1)
[-3.92264613920915...]
sage: RBF(2).zetaderiv(3)
[-6.0001458028430...]
>>> from sage.all import *
>>> RBF(Integer(1)/Integer(2)).zetaderiv(Integer(1))
[-3.92264613920915...]
>>> RBF(Integer(2)).zetaderiv(Integer(3))
[-6.0001458028430...]
class sage.rings.real_arb.RealBallField(precision=53)[source]#

Bases: UniqueRepresentation, RealBallField

An approximation of the field of real numbers using mid-rad intervals, also known as balls.

INPUT:

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

EXAMPLES:

sage: RBF = RealBallField() # indirect doctest
sage: RBF(1)
1.000000000000000
>>> from sage.all import *
>>> RBF = RealBallField() # indirect doctest
>>> RBF(Integer(1))
1.000000000000000
sage: (1/2*RBF(1)) + AA(sqrt(2)) - 1 + polygen(QQ, 'x')                         # needs sage.symbolic
x + [0.914213562373095 +/- ...e-16]
>>> from sage.all import *
>>> (Integer(1)/Integer(2)*RBF(Integer(1))) + AA(sqrt(Integer(2))) - Integer(1) + polygen(QQ, 'x')                         # needs sage.symbolic
x + [0.914213562373095 +/- ...e-16]

See also

Element[source]#

alias of RealBall

algebraic_closure()[source]#

Return the complex ball field with the same precision.

EXAMPLES:

sage: from sage.rings.complex_arb import ComplexBallField
sage: RBF.complex_field()
Complex ball field with 53 bits of precision
sage: RealBallField(3).algebraic_closure()
Complex ball field with 3 bits of precision
>>> from sage.all import *
>>> from sage.rings.complex_arb import ComplexBallField
>>> RBF.complex_field()
Complex ball field with 53 bits of precision
>>> RealBallField(Integer(3)).algebraic_closure()
Complex ball field with 3 bits of precision
bell_number(n)[source]#

Return a ball enclosing the n-th Bell number.

EXAMPLES:

sage: [RBF.bell_number(n) for n in range(7)]
[1.000000000000000,
 1.000000000000000,
 2.000000000000000,
 5.000000000000000,
 15.00000000000000,
 52.00000000000000,
 203.0000000000000]
sage: RBF.bell_number(-1)
Traceback (most recent call last):
...
ValueError: expected a nonnegative index
sage: RBF.bell_number(10**20)
[5.38270113176282e+1794956117137290721328 +/- ...e+1794956117137290721313]
>>> from sage.all import *
>>> [RBF.bell_number(n) for n in range(Integer(7))]
[1.000000000000000,
 1.000000000000000,
 2.000000000000000,
 5.000000000000000,
 15.00000000000000,
 52.00000000000000,
 203.0000000000000]
>>> RBF.bell_number(-Integer(1))
Traceback (most recent call last):
...
ValueError: expected a nonnegative index
>>> RBF.bell_number(Integer(10)**Integer(20))
[5.38270113176282e+1794956117137290721328 +/- ...e+1794956117137290721313]
bernoulli(n)[source]#

Return a ball enclosing the n-th Bernoulli number.

EXAMPLES:

sage: [RBF.bernoulli(n) for n in range(4)]
[1.000000000000000, -0.5000000000000000, [0.1666666666666667 +/- ...e-17], 0]
sage: RBF.bernoulli(2**20)
[-1.823002872104961e+5020717 +/- ...e+5020701]
sage: RBF.bernoulli(2**1000)
Traceback (most recent call last):
...
ValueError: argument too large
>>> from sage.all import *
>>> [RBF.bernoulli(n) for n in range(Integer(4))]
[1.000000000000000, -0.5000000000000000, [0.1666666666666667 +/- ...e-17], 0]
>>> RBF.bernoulli(Integer(2)**Integer(20))
[-1.823002872104961e+5020717 +/- ...e+5020701]
>>> RBF.bernoulli(Integer(2)**Integer(1000))
Traceback (most recent call last):
...
ValueError: argument too large
catalan_constant()[source]#

Return a ball enclosing the Catalan constant.

EXAMPLES:

sage: RBF.catalan_constant()
[0.915965594177219 +/- ...e-16]
sage: RealBallField(128).catalan_constant()
[0.91596559417721901505460351493238411077 +/- ...e-39]
>>> from sage.all import *
>>> RBF.catalan_constant()
[0.915965594177219 +/- ...e-16]
>>> RealBallField(Integer(128)).catalan_constant()
[0.91596559417721901505460351493238411077 +/- ...e-39]
characteristic()[source]#

Real ball fields have characteristic zero.

EXAMPLES:

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

Return the complex ball field with the same precision.

EXAMPLES:

sage: from sage.rings.complex_arb import ComplexBallField
sage: RBF.complex_field()
Complex ball field with 53 bits of precision
sage: RealBallField(3).algebraic_closure()
Complex ball field with 3 bits of precision
>>> from sage.all import *
>>> from sage.rings.complex_arb import ComplexBallField
>>> RBF.complex_field()
Complex ball field with 53 bits of precision
>>> RealBallField(Integer(3)).algebraic_closure()
Complex ball field with 3 bits of precision
construction()[source]#

Return the construction of a real ball field as a completion of the rationals.

EXAMPLES:

sage: RBF = RealBallField(42)
sage: functor, base = RBF.construction()
sage: functor, base
(Completion[+Infinity, prec=42], Rational Field)
sage: functor(base) is RBF
True
>>> from sage.all import *
>>> RBF = RealBallField(Integer(42))
>>> functor, base = RBF.construction()
>>> functor, base
(Completion[+Infinity, prec=42], Rational Field)
>>> functor(base) is RBF
True
cospi(x)[source]#

Return a ball enclosing \(\cos(\pi x)\).

This works even if x itself is not a ball, and may be faster or more accurate where x is a rational number.

EXAMPLES:

sage: RBF.cospi(1)
-1.000000000000000
sage: RBF.cospi(1/3)
0.5000000000000000
>>> from sage.all import *
>>> RBF.cospi(Integer(1))
-1.000000000000000
>>> RBF.cospi(Integer(1)/Integer(3))
0.5000000000000000

See also

cos()

double_factorial(n)[source]#

Return a ball enclosing the n-th double factorial.

EXAMPLES:

sage: [RBF.double_factorial(n) for n in range(7)]
[1.000000000000000,
 1.000000000000000,
 2.000000000000000,
 3.000000000000000,
 8.000000000000000,
 15.00000000000000,
 48.00000000000000]
sage: RBF.double_factorial(2**20)
[1.448372990...e+2928836 +/- ...]
sage: RBF.double_factorial(2**1000)
Traceback (most recent call last):
...
ValueError: argument too large
sage: RBF.double_factorial(-1)
Traceback (most recent call last):
...
ValueError: expected a nonnegative index
>>> from sage.all import *
>>> [RBF.double_factorial(n) for n in range(Integer(7))]
[1.000000000000000,
 1.000000000000000,
 2.000000000000000,
 3.000000000000000,
 8.000000000000000,
 15.00000000000000,
 48.00000000000000]
>>> RBF.double_factorial(Integer(2)**Integer(20))
[1.448372990...e+2928836 +/- ...]
>>> RBF.double_factorial(Integer(2)**Integer(1000))
Traceback (most recent call last):
...
ValueError: argument too large
>>> RBF.double_factorial(-Integer(1))
Traceback (most recent call last):
...
ValueError: expected a nonnegative index
euler_constant()[source]#

Return a ball enclosing the Euler constant.

EXAMPLES:

sage: RBF.euler_constant() # abs tol 1e-15
[0.5772156649015329 +/- 9.00e-17]
sage: RealBallField(128).euler_constant()
[0.57721566490153286060651209008240243104 +/- ...e-39]
>>> from sage.all import *
>>> RBF.euler_constant() # abs tol 1e-15
[0.5772156649015329 +/- 9.00e-17]
>>> RealBallField(Integer(128)).euler_constant()
[0.57721566490153286060651209008240243104 +/- ...e-39]
fibonacci(n)[source]#

Return a ball enclosing the n-th Fibonacci number.

EXAMPLES:

sage: [RBF.fibonacci(n) for n in range(7)]
[0,
1.000000000000000,
1.000000000000000,
2.000000000000000,
3.000000000000000,
5.000000000000000,
8.000000000000000]
sage: RBF.fibonacci(-2)
-1.000000000000000
sage: RBF.fibonacci(10**20)
[3.78202087472056e+20898764024997873376 +/- ...e+20898764024997873361]
>>> from sage.all import *
>>> [RBF.fibonacci(n) for n in range(Integer(7))]
[0,
1.000000000000000,
1.000000000000000,
2.000000000000000,
3.000000000000000,
5.000000000000000,
8.000000000000000]
>>> RBF.fibonacci(-Integer(2))
-1.000000000000000
>>> RBF.fibonacci(Integer(10)**Integer(20))
[3.78202087472056e+20898764024997873376 +/- ...e+20898764024997873361]
gamma(x)[source]#

Return a ball enclosing the gamma function of x.

This works even if x itself is not a ball, and may be more efficient in the case where x is an integer or a rational number.

EXAMPLES:

sage: RBF.gamma(5)
24.00000000000000
sage: RBF.gamma(10**20)
[1.932849514310098...+1956570551809674817225 +/- ...]
sage: RBF.gamma(1/3)
[2.678938534707747 +/- ...e-16]
sage: RBF.gamma(-5)
nan
>>> from sage.all import *
>>> RBF.gamma(Integer(5))
24.00000000000000
>>> RBF.gamma(Integer(10)**Integer(20))
[1.932849514310098...+1956570551809674817225 +/- ...]
>>> RBF.gamma(Integer(1)/Integer(3))
[2.678938534707747 +/- ...e-16]
>>> RBF.gamma(-Integer(5))
nan

See also

gamma()

gens()[source]#

EXAMPLES:

sage: RBF.gens()
(1.000000000000000,)
sage: RBF.gens_dict()
{'1.000000000000000': 1.000000000000000}
>>> from sage.all import *
>>> RBF.gens()
(1.000000000000000,)
>>> RBF.gens_dict()
{'1.000000000000000': 1.000000000000000}
is_exact()[source]#

Real ball fields are not exact.

EXAMPLES:

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

Return a ball enclosing \(\log(2)\).

EXAMPLES:

sage: RBF.log2()
[0.6931471805599453 +/- ...e-17]
sage: RealBallField(128).log2()
[0.69314718055994530941723212145817656807 +/- ...e-39]
>>> from sage.all import *
>>> RBF.log2()
[0.6931471805599453 +/- ...e-17]
>>> RealBallField(Integer(128)).log2()
[0.69314718055994530941723212145817656807 +/- ...e-39]
maximal_accuracy()[source]#

Return the relative accuracy of exact elements measured in bits.

OUTPUT:

An integer.

EXAMPLES:

sage: RBF.maximal_accuracy()
9223372036854775807 # 64-bit
2147483647          # 32-bit
>>> from sage.all import *
>>> RBF.maximal_accuracy()
9223372036854775807 # 64-bit
2147483647          # 32-bit
pi()[source]#

Return a ball enclosing \(\pi\).

EXAMPLES:

sage: RBF.pi()
[3.141592653589793 +/- ...e-16]
sage: RealBallField(128).pi()
[3.1415926535897932384626433832795028842 +/- ...e-38]
>>> from sage.all import *
>>> RBF.pi()
[3.141592653589793 +/- ...e-16]
>>> RealBallField(Integer(128)).pi()
[3.1415926535897932384626433832795028842 +/- ...e-38]
prec()[source]#

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

EXAMPLES:

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

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

EXAMPLES:

sage: RealBallField().precision()
53
>>> from sage.all import *
>>> RealBallField().precision()
53
sinpi(x)[source]#

Return a ball enclosing \(\sin(\pi x)\).

This works even if x itself is not a ball, and may be faster or more accurate where x is a rational number.

EXAMPLES:

sage: RBF.sinpi(1)
0
sage: RBF.sinpi(1/3)
[0.866025403784439 +/- ...e-16]
sage: RBF.sinpi(1 + 2^(-100))
[-2.478279624546525e-30 +/- ...e-46]
>>> from sage.all import *
>>> RBF.sinpi(Integer(1))
0
>>> RBF.sinpi(Integer(1)/Integer(3))
[0.866025403784439 +/- ...e-16]
>>> RBF.sinpi(Integer(1) + Integer(2)**(-Integer(100)))
[-2.478279624546525e-30 +/- ...e-46]

See also

sin()

some_elements()[source]#

Real ball fields contain exact balls, inexact balls, infinities, and more.

EXAMPLES:

sage: RBF.some_elements()
[0, 1.000000000000000, [0.3333333333333333 +/- ...e-17],
[-4.733045976388941e+363922934236666733021124 +/- ...e+363922934236666733021108],
[+/- inf], [+/- inf], [+/- inf], nan]
>>> from sage.all import *
>>> RBF.some_elements()
[0, 1.000000000000000, [0.3333333333333333 +/- ...e-17],
[-4.733045976388941e+363922934236666733021124 +/- ...e+363922934236666733021108],
[+/- inf], [+/- inf], [+/- inf], nan]
zeta(s)[source]#

Return a ball enclosing the Riemann zeta function of s.

This works even if s itself is not a ball, and may be more efficient in the case where s is an integer.

EXAMPLES:

sage: RBF.zeta(3)
[1.202056903159594 +/- ...e-16]
sage: RBF.zeta(1)
nan
sage: RBF.zeta(1/2)
[-1.460354508809587 +/- ...e-16]
>>> from sage.all import *
>>> RBF.zeta(Integer(3))
[1.202056903159594 +/- ...e-16]
>>> RBF.zeta(Integer(1))
nan
>>> RBF.zeta(Integer(1)/Integer(2))
[-1.460354508809587 +/- ...e-16]

See also

zeta()

sage.rings.real_arb.create_RealBall(parent, serialized)[source]#

Create a RealBall from a serialized representation.