arb – real numbers

class flint.arb(mid=None, rad=None)

Represents a real number \(x\) by a midpoint \(m\) and a radius \(r\) such that \(x \in [m \pm r] = [m-r, m+r]\). The midpoint and radius are both floating-point numbers. The radius uses a fixed, implementation-defined precision (30 bits). The precision used for midpoints is controlled by ctx.prec (bits) or equivalently ctx.dps (digits).

The constructor accepts a midpoint mid and a radius rad, either of which defaults to zero if omitted. The arguments can be tuples \((a, b)\) representing exact floating-point data \(a 2^b\), integers, floating-point numbers, rational strings, or decimal strings. If the radius is nonzero, it might be rounded up to a slightly larger value than the exact value passed by the user.

>>> arb(10.25)
10.2500000000000
>>> print(1 / arb(4))  # exact
0.250000000000000
>>> print(1 / arb(3))  # approximate
[0.333333333333333 +/- 3.71e-16]
>>> print(arb("3.0"))
3.00000000000000
>>> print(arb("0.1"))
[0.100000000000000 +/- 2.23e-17]
>>> print(arb("1/10"))
[0.100000000000000 +/- 2.23e-17]
>>> print(arb("3.14159 +/- 0.00001"))
[3.1416 +/- 2.01e-5]
>>> ctx.dps = 50
>>> print(arb("1/3"))
[0.33333333333333333333333333333333333333333333333333 +/- 3.78e-51]
>>> ctx.default()

Converting to or from decimal results in some loss of accuracy. See arb.str() for details.

abs_lower(self)

Lower bound for the absolute value of self. The output is an arb holding an exact floating-point number that has been rounded down to the current precision.

>>> print(arb("-5 +/- 2").abs_lower().str(5, radius=False))
3.0000
abs_upper(self)

Upper bound for the absolute value of self. The output is an arb holding an exact floating-point number that has been rounded up to the current precision.

>>> print(arb("-5 +/- 2").abs_upper().str(5, radius=False))
7.0000
acos(s)

Inverse cosine \(\operatorname{acos}(s)\).

>>> from flint import showgood
>>> showgood(lambda: arb(0).acos(), dps=25)
1.570796326794896619231322
acosh(s)

Inverse hyperbolic cosine \(\operatorname{acosh}(s)\).

>>> from flint import showgood
>>> showgood(lambda: arb(2).acosh(), dps=25)
1.316957896924816708625046
agm(s, t=1)

The arithmetic-geometric mean \(M(s,t)\), or \(M(s) = M(s,1)\) if no extra parameter is passed.

>>> from flint import showgood
>>> showgood(lambda: arb(2).sqrt().agm(), dps=25)
1.198140234735592207439922
>>> showgood(lambda: arb(2).agm(100), dps=25)
29.64467336236643624632443
>>> showgood(lambda: arb(0).agm(0), dps=25)
0
>>> showgood(arb(0).agm, dps=25)
0
>>> showgood(arb(-2).agm, dps=25)   # not defined for negatives
Traceback (most recent call last):
  ...
ValueError: no convergence (maxprec=960, try higher maxprec)

This function is undefined for negative input. Use acb.agm() for the complex extension.

airy(s)

Computes the Airy function values \(\operatorname{Ai}(s)\), \(\operatorname{Ai}'(s)\), \(\operatorname{Bi}(s)\), \(\operatorname{Bi}'(s)\) simultaneously, returning a tuple.

>>> from flint import showgood
>>> showgood(lambda: arb(-1).airy(), dps=10)
(0.5355608833, -0.01016056712, 0.1039973895, 0.5923756264)
airy_ai(s, int derivative=0)

Airy function \(\operatorname{Ai}(s)\), or \(\operatorname{Ai}'(s)\) if derivative is 1.

>>> from flint import showgood
>>> showgood(lambda: arb(-1).airy_ai(), dps=25)
0.5355608832923521187995166
>>> showgood(lambda: arb(-1).airy_ai(derivative=1), dps=25)
-0.01016056711664520939504547
static airy_ai_zero(n, int derivative=0)

For positive integer n, returns the zero \(a_n\) of the Airy function \(\operatorname{Ai}(s)\), or the corresponding zero \(a'_n\) of \(\operatorname{Ai}'(s)\) if derivative is 1.

>>> from flint import showgood
>>> showgood(lambda: arb.airy_ai_zero(1), dps=25)
-2.338107410459767038489197
>>> showgood(lambda: arb.airy_ai_zero(1000), dps=25)
-281.0315196125215528353364
>>> showgood(lambda: arb.airy_ai_zero(1, derivative=1), dps=25)
-1.018792971647471089017325
airy_bi(s, int derivative=0)

Airy function \(\operatorname{Bi}(s)\), or \(\operatorname{Bi}'(s)\) if derivative is 1.

>>> from flint import showgood
>>> showgood(lambda: arb(-1).airy_bi(), dps=25)
0.1039973894969446118886900
>>> showgood(lambda: arb(-1).airy_bi(derivative=1), dps=25)
0.5923756264227923508167792
static airy_bi_zero(n, int derivative=0)

For positive integer n, returns the zero \(b_n\) of the Airy function \(\operatorname{Bi}(s)\), or the corresponding zero \(b'_n\) of \(\operatorname{Bi}'(s)\) if derivative is 1.

>>> from flint import showgood
>>> showgood(lambda: arb.airy_bi_zero(1), dps=25)
-1.173713222709127924919980
>>> showgood(lambda: arb.airy_bi_zero(1000), dps=25)
-280.9378112034152401578834
>>> showgood(lambda: arb.airy_bi_zero(1, derivative=1), dps=25)
-2.294439682614123246622459
asin(s)

Inverse sine \(\operatorname{asin}(s)\).

>>> from flint import showgood
>>> showgood(lambda: arb(1).asin(), dps=25)
1.570796326794896619231322
asinh(s)

Inverse hyperbolic sine \(\operatorname{asinh}(s)\).

>>> from flint import showgood
>>> showgood(lambda: arb(1000).asinh(), dps=25)
7.600902709541988611523290
atan(s)

Inverse tangent \(\operatorname{atan}(s)\).

>>> from flint import showgood
>>> showgood(lambda: arb(1).atan(), dps=25)
0.7853981633974483096156608
static atan2(s, t)

Two-argument inverse tangent \(\operatorname{atan2}(s,t)\).

>>> from flint import showgood
>>> showgood(lambda: arb.atan2(-10,-5), dps=25)
-2.034443935795702735445578
atanh(s)

Inverse hyperbolic tangent \(\operatorname{atanh}(s)\).

>>> from flint import showgood
>>> showgood(lambda: arb("0.99").atanh(), dps=25)
2.646652412362246197705061
backlund_s(x)

Backlund S function related to the Riemann zeta function.

>>> from flint import showgood
>>> showgood(lambda: arb(123).backlund_s(), dps=25)
0.4757920863536796196115749
static bell_number(n)

Computes the Bell number \(B_n\) as an arb, where n is a given integer.

>>> from flint import showgood
>>> showgood(lambda: arb.bell_number(1), dps=25)
1.000000000000000000000000
>>> showgood(lambda: arb.bell_number(10), dps=25)
115975.0000000000000000000
>>> showgood(lambda: arb.bell_number(10**20), dps=25)
5.382701131762816107395343e+1794956117137290721328
static bernoulli(n)

Computes the Bernoulli number \(B_n\) as an arb, where n is a given integer.

>>> from flint import showgood
>>> showgood(lambda: arb.bernoulli(1), dps=25)
-0.5000000000000000000000000
>>> showgood(lambda: arb.bernoulli(2), dps=25)
0.1666666666666666666666667
>>> showgood(lambda: arb.bernoulli(10**7), dps=25)
-4.983176441416329828329241e+57675260
>>> showgood(lambda: arb.bernoulli(10**30), dps=25)
-5.048207707665410387507573e+28767525649738633122783863898083
bernoulli_poly(s, ulong n)

Returns the value of the Bernoulli polynomial \(B_n(s)\).

>>> from flint import showgood
>>> showgood(lambda: arb("3.141").bernoulli_poly(5), dps=25)
113.5165117028492010000000
bessel_i(self, n, bool scaled=False)

Bessel function \(I_n(z)\), where the argument z is given by self and the order n is passed as an extra parameter. Optionally a scaled Bessel function can be computed.

>>> from flint import showgood
>>> showgood(lambda: arb(5).bessel_i(1), dps=25)
24.33564214245052719914305
>>> showgood(lambda: arb(5).bessel_i(1, scaled=True), dps=25)
0.1639722669445423569261229
bessel_j(self, n)

Bessel function \(J_n(z)\), where the argument z is given by self and the order n is passed as an extra parameter.

>>> from flint import showgood
>>> showgood(lambda: arb(5).bessel_j(1), dps=25)
-0.3275791375914652220377343
bessel_k(self, n, bool scaled=False)

Bessel function \(K_n(z)\), where the argument z is given by self and the order n is passed as an extra parameter. Optionally a scaled Bessel function can be computed.

>>> from flint import showgood
>>> showgood(lambda: arb(5).bessel_k(1), dps=25)
0.004044613445452164208365022
>>> showgood(lambda: arb(5).bessel_k(1, scaled=True), dps=25)
0.6002738587883125829360457
bessel_y(self, n)

Bessel function \(Y_n(z)\), where the argument z is given by self and the order n is passed as an extra parameter.

>>> from flint import showgood
>>> showgood(lambda: arb(5).bessel_y(1), dps=25)
0.1478631433912268448010507
beta_lower(self, a, b, int regularized=0)

Lower incomplete beta function \(B(a,b;z)\). The argument z is given by self and the parameters a and b are passed as extra function arguments. Optionally the regularized version of this function can be computed.

>>> from flint import showgood
>>> showgood(lambda: arb("0.9").beta_lower(1, 2.5), dps=25)
0.3987350889359326482672004
>>> showgood(lambda: arb("0.9").beta_lower(1, 2.5, regularized=True), dps=25)
0.9968377223398316206680011
bin(s, ulong k)

Binomial coefficient \({s \choose k}\). Currently k is limited to an integer; this restriction will be removed in the future by using the gamma function.

>>> from flint import arb, ctx, showgood
>>> ctx.prec = 53
>>> print(arb(10).bin(5))
252.000000000000
>>> showgood(lambda: arb.pi().bin(100), dps=25)
5.478392395095119521549286e-9
static bin_uiui(ulong n, ulong k)

Binomial coefficient \({n \choose k}\).

>>> from flint import arb, ctx
>>> ctx.prec = 53
>>> print(arb.bin_uiui(10, 5))
252.000000000000
bits(self)

Returns number of bits needed to represent absolute value of mantissa of the midpoint; returns 0 if midpoint is special value.

>>> arb("2047/2048").bits()
11
ceil(s)

Ceiling function \(\lceil s \rceil\).

>>> print(arb.pi().ceil())
4.00000000000000
>>> print((arb.pi() - arb.pi()).ceil().str(more=True))
[0.500000000000000 +/- 0.501]
chebyshev_t(s, n)

Chebyshev function of the first kind \(T_n(s)\).

>>> from flint import showgood
>>> showgood(lambda: (arb(1)/3).chebyshev_t(3), dps=25)
-0.8518518518518518518518519
chebyshev_u(s, n)

Chebyshev function of the second kind \(U_n(s)\).

>>> from flint import showgood
>>> showgood(lambda: (arb(1)/3).chebyshev_u(3), dps=25)
-1.037037037037037037037037
chi(s)

Hyperbolic cosine integral \(\operatorname{Chi}(s)\).

>>> from flint import showgood
>>> showgood(lambda: arb(3).chi(), dps=25)
4.960392094765609760297918
ci(s)

Cosine integral \(\operatorname{Ci}(s)\).

>>> from flint import showgood
>>> showgood(lambda: arb(3).ci(), dps=25)
0.1196297860080003276264723
static const_catalan()

Catalan’s constant.

>>> from flint import showgood
>>> showgood(arb.const_catalan, dps=25)
0.9159655941772190150546035
static const_e()

The constant \(e\).

>>> from flint import showgood
>>> showgood(arb.const_e, dps=25)
2.718281828459045235360287
static const_euler()

Euler’s constant \(\gamma\).

>>> from flint import showgood
>>> showgood(arb.const_euler, dps=25)
0.5772156649015328606065121
static const_glaisher()

Glaisher’s constant.

>>> from flint import showgood
>>> showgood(arb.const_glaisher, dps=25)
1.282427129100622636875343
static const_khinchin()

Khinchin’s constant \(K\).

>>> from flint import showgood
>>> showgood(arb.const_khinchin, dps=25)
2.685452001065306445309715
static const_log10()

The constant \(\log(10)\).

>>> from flint import showgood
>>> showgood(arb.const_log10, dps=25)
2.302585092994045684017991
static const_log2()

The constant \(\log(2)\).

>>> from flint import showgood
>>> showgood(arb.const_log2, dps=25)
0.6931471805599453094172321
static const_sqrt_pi()

The constant \(\sqrt{\pi}\).

>>> from flint import showgood
>>> showgood(arb.const_sqrt_pi, dps=25)
1.772453850905516027298167
contains(self, other)
contains_integer(self)
contains_interior(self, other)
cos(s)

Cosine function \(\cos(s)\).

>>> from flint import showgood
>>> showgood(lambda: arb(1).cos(), dps=25)
0.5403023058681397174009366
cos_pi(s)

Cosine function \(\cos(\pi s)\).

>>> from flint import showgood
>>> showgood(lambda: arb(0.75).cos_pi(), dps=25)
-0.7071067811865475244008444
>>> showgood(lambda: arb(0.5).cos_pi(), dps=25)
0
static cos_pi_fmpq(fmpq s)

Returns the algebraic cosine value \(\cos(\pi s)\).

>>> from flint import fmpq
>>> from flint import showgood
>>> showgood(lambda: arb.cos_pi_fmpq(fmpq(3,4)), dps=25)
-0.7071067811865475244008444
cosh(s)

Hyperbolic cosine \(\cosh(s)\).

>>> from flint import showgood
>>> showgood(lambda: arb(1).cosh(), dps=25)
1.543080634815243778477906
cot(s)

Cotangent function \(\cot(s)\).

>>> from flint import showgood
>>> showgood(lambda: arb(1).cot(), dps=25)
0.6420926159343307030064200
cot_pi(s)

Cotangent function \(\cot(\pi s)\).

>>> from flint import showgood
>>> showgood(lambda: arb(0.125).cot_pi(), dps=25)
2.414213562373095048801689
>>> showgood(lambda: arb(0.5).cot_pi(), dps=25)
0
coth(s)

Hyperbolic cotangent \(\coth(s)\).

>>> from flint import showgood
>>> showgood(lambda: arb(1).coth(), dps=25)
1.313035285499331303636161
coulomb(self, l, eta)

Computes the Coulomb wave functions \(F_{\ell}(\eta,z)\), \(G_{\ell}(\eta,z)\), where z is given by self. Both function values are computed simultaneously and a tuple is returned.

>>> from flint import showgood
>>> showgood(lambda: arb(1).coulomb(0.5, 0.25), dps=10)
(0.4283180781, 1.218454487)
coulomb_f(self, l, eta)

Regular Coulomb wave function \(F_{\ell}(\eta,z)\) where z is given by self.

>>> from flint import showgood
>>> showgood(lambda: arb(1).coulomb_f(0.5, 0.25), dps=25)
0.4283180781043541845555944
coulomb_g(self, l, eta)

Irregular Coulomb wave function \(G_{\ell}(\eta,z)\) where z is given by self.

>>> from flint import showgood
>>> showgood(lambda: arb(1).coulomb_g(0.5, 0.25), dps=25)
1.218454487206367973745641
csc(s)

Cosecant function \(\operatorname{csc}(s)\).

>>> from flint import showgood
>>> showgood(lambda: arb(1).csc(), dps=25)
1.188395105778121216261599
csch(s)

Hyperbolic cosecant \(\operatorname{csch}(s)\).

>>> from flint import showgood
>>> showgood(lambda: arb(1).csch(), dps=25)
0.8509181282393215451338428
digamma(s)

Digamma function \(\psi(s)\).

>>> from flint import showgood
>>> showgood(lambda: arb(1).digamma(), dps=25)
-0.5772156649015328606065121
ei(s)

Exponential integral \(\operatorname{Ei}(s)\).

>>> from flint import showgood
>>> showgood(lambda: arb(3).ei(), dps=25)
9.933832570625416558008336
erf(s)

Error function \(\operatorname{erf}(s)\).

>>> from flint import showgood
>>> showgood(lambda: arb(3).erf(), dps=25)
0.9999779095030014145586272
erfc(s)

Complementary error function \(\operatorname{erfc}(s)\).

>>> from flint import showgood
>>> showgood(lambda: arb(3).erfc(), dps=25)
2.209049699858544137277613e-5
erfcinv(s)

Inverse complementary error function \(\operatorname{erfcinv}(s)\).

>>> from flint import showgood
>>> showgood(lambda: arb(.25).erfcinv(), dps=25)
0.8134198475976185416902894
erfi(s)

Imaginary error function \(\operatorname{erfi}(s)\).

>>> from flint import showgood
>>> showgood(lambda: arb(3).erfi(), dps=25)
1629.994622601565651061648
erfinv(s)

Inverse error function \(\operatorname{erfinv}(s)\).

>>> from flint import showgood
>>> showgood(lambda: arb(.25).erfinv(), dps=25)
0.2253120550121781047250140
exp(s)

Exponential function \(\exp(s)\).

>>> from flint import showgood
>>> showgood(lambda: arb(1).exp(), dps=25)
2.718281828459045235360287
expint(self, s)

Generalized exponential integral \(E_s(z)\). The argument z is given by self and the order s is passed as an extra parameter.

>>> from flint import showgood
>>> showgood(lambda: arb(5).expint(2), dps=25)
0.0009964690427088381099832386
expm1(s)

Exponential function \(\exp(s) - 1\), computed accurately for small s.

>>> from flint import showgood
>>> showgood(lambda: (arb(10) ** -8).expm1(), dps=25)
1.000000005000000016666667e-8
fac(s)

Factorial, equivalent to \(s! = \Gamma(s+1)\).

>>> from flint import showgood
>>> showgood(lambda: arb(5).fac(), dps=25)
120.0000000000000000000000
>>> showgood(lambda: arb("0.1").fac(), dps=25)
0.9513507698668731836292487
static fac_ui(ulong n)

Factorial \(n!\), given an integer.

>>> from flint import arb, ctx
>>> ctx.prec = 53
>>> print(arb.fac_ui(10))
3628800.00000000
>>> from flint import showgood
>>> showgood(lambda: arb.fac_ui(10**9).log(), dps=25)
19723265848.22698260792313
static fib(n)

Computes the Fibonacci number \(F_n\) as an arb, where n is a given integer.

>>> from flint import arb, ctx
>>> ctx.prec = 53
>>> print(arb.fib(10))
55.0000000000000
>>> from flint import showgood
>>> showgood(lambda: arb.fib(10**100).log(), dps=25)
4.812118250596034474977589e+99
floor(s)

Floor function \(\lfloor s \rfloor\).

>>> from flint import arb, ctx
>>> ctx.prec = 53
>>> print(arb.pi().floor())
3.00000000000000
>>> print((arb.pi() - arb.pi()).floor().str(more=True))
[-0.500000000000000 +/- 0.501]
fmpq(self)
fmpz(self)
fresnel_c(s, bool normalized=True)

Fresnel cosine integral \(C(s)\), optionally not normalized.

>>> from flint import showgood
>>> showgood(lambda: arb(3).fresnel_c(), dps=25)
0.6057207892976856295561611
>>> showgood(lambda: arb(3).fresnel_c(normalized=False), dps=25)
0.7028635577302687301744099
fresnel_s(s, bool normalized=True)

Fresnel sine integral \(S(s)\), optionally not normalized.

>>> from flint import showgood
>>> showgood(lambda: arb(3).fresnel_s(), dps=25)
0.4963129989673750360976123
>>> showgood(lambda: arb(3).fresnel_s(normalized=False), dps=25)
0.7735625268937690171497722
gamma(s)

Gamma function \(\Gamma(s)\).

>>> from flint import showgood
>>> showgood(lambda: arb(10).gamma(), dps=25)
362880.0000000000000000000
>>> showgood(lambda: arb(-2.5).gamma(), dps=25)
-0.9453087204829418812256893
>>> showgood(lambda: (arb.pi() ** 10).gamma(), dps=25)
1.705646271897306403570389e+424898
>>> showgood(lambda: arb(0).gamma(), dps=25)  # pole
Traceback (most recent call last):
  ...
ValueError: no convergence (maxprec=960, try higher maxprec)
static gamma_fmpq(fmpq s)

Computes the gamma function \(\Gamma(s)\) of a given fmpq s, exploiting the fact that s is an exact rational number to improve performance.

>>> from flint import fmpq
>>> from flint import showgood
>>> showgood(lambda: arb.gamma_fmpq(fmpq(1,4)), dps=25)
3.625609908221908311930685
gamma_lower(self, s, int regularized=0)

Lower incomplete gamma function \(\gamma(s,z)\). The argument z is given by self and the order s is passed as an extra parameter. Optionally the regularized versions \(P(s,z)\) and \(\gamma^{*}(s,z) = z^{-s} P(s,z)\) can be computed.

>>> from flint import showgood
>>> showgood(lambda: arb(5).gamma_lower(2.5), dps=25)
1.229327136861979598135209
>>> showgood(lambda: arb(5).gamma_lower(2.5, regularized=1), dps=25)
0.9247647538534878212779231
>>> showgood(lambda: arb(5).gamma_lower(2.5, regularized=2), dps=25)
0.01654269482249807489997922
gamma_upper(self, s, int regularized=0)

Upper incomplete gamma function \(\Gamma(s,z)\). The argument z is given by self and the order s is passed as an extra parameter. Optionally the regularized version \(Q(s,z)\) can be computed.

>>> from flint import showgood
>>> showgood(lambda: arb(5).gamma_upper(2.5), dps=25)
0.1000132513171574223384170
>>> showgood(lambda: arb(5).gamma_upper(2.5, regularized=1), dps=25)
0.07523524614651217872207687
gegenbauer_c(s, n, m)

Gegenbauer function \(C_n^{m}(s)\).

>>> from flint import showgood
>>> showgood(lambda: (arb(1)/3).gegenbauer_c(5, 0.25), dps=25)
0.1321855709876543209876543
static gram_point(n)

Returns the n-th Gram point.

>>> from flint import showgood
>>> showgood(lambda: arb.gram_point(-1), dps=25)
9.666908056130192141261536
>>> showgood(lambda: arb.gram_point(0), dps=25)
17.84559954041086081682634
>>> showgood(lambda: arb.gram_point(10**30), dps=25)
9.829776286927442475869051e+28
hermite_h(s, n)

Hermite function \(H_n(s)\).

>>> from flint import showgood
>>> showgood(lambda: (arb(1)/3).hermite_h(5), dps=25)
34.20576131687242798353909
hypgeom(self, a, b, bool regularized=False)

Generalized hypergeometric function \({}_pF_q(a;b;z)\). The argument z is given by self and a and b are additional lists of complex numbers defining the parameters. Optionally the regularized hypergeometric function can be computed.

>>> from flint import showgood
>>> showgood(lambda: arb(5).hypgeom([1,2,3],[5,4.5,6]), dps=25)
1.301849229178968153998454
>>> showgood(lambda: arb(5).hypgeom([1,2,3],[5,4.5,6],regularized=True), dps=25)
3.886189282817193519132054e-5
hypgeom_0f1(self, a, bool regularized=False)

The hypergeometric function \({}_0F_1(a,z)\) where the argument z is given by self Optionally the regularized version of this function can be computed.

>>> from flint import showgood
>>> showgood(lambda: arb(-5).hypgeom_0f1(2.5), dps=25)
0.003114611044402738470826907
>>> showgood(lambda: arb(-5).hypgeom_0f1(2.5, regularized=True), dps=25)
0.002342974810739764377177885
hypgeom_1f1(self, a, b, bool regularized=False)

The hypergeometric function \({}_1F_1(a,b,z)\) where the argument z is given by self Optionally the regularized version of this function can be computed.

>>> from flint import showgood
>>> showgood(lambda: arb(-5).hypgeom_1f1(1, 2.5), dps=25)
0.2652884733179767675077230
>>> showgood(lambda: arb(-5).hypgeom_1f1(1, 2.5, regularized=True), dps=25)
0.1995639910417191573048664
hypgeom_2f1(self, a, b, c, bool regularized=False, bool ab=False, bool ac=False, bc=False, abc=False)

The hypergeometric function \({}_2F_1(a,b,c,z)\) where the argument z is given by self Optionally the regularized version of this function can be computed.

>>> from flint import showgood
>>> showgood(lambda: arb(-5).hypgeom_2f1(1,2,3), dps=25)
0.2566592424617555999350018
>>> showgood(lambda: arb(-5).hypgeom_2f1(1,2,3,regularized=True), dps=25)
0.1283296212308777999675009

The flags ab, ac, bc, abc can be used to specify whether the parameter differences \(a-b\), \(a-c\), \(b-c\) and \(a+b-c\) represent exact integers, even if the input intervals are inexact. If the parameters are exact, these flags are not needed.

>>> from flint import showgood
>>> showgood(lambda: arb("9/10").hypgeom_2f1(arb(2).sqrt(), 0.5, arb(2).sqrt()+1.5, abc=True), dps=25)
1.447530478120770807945697
>>> showgood(lambda: arb("9/10").hypgeom_2f1(arb(2).sqrt(), 0.5, arb(2).sqrt()+1.5), dps=25)
Traceback (most recent call last):
  ...
ValueError: no convergence (maxprec=960, try higher maxprec)
hypgeom_u(self, a, b)

The hypergeometric function \(U(a,b,z)\) where the argument z is given by self

>>> from flint import showgood
>>> showgood(lambda: arb(5).hypgeom_u(1, 2.5), dps=25)
0.2184157028890778783289036
imag
intersection(s, t)

Returns a ball containing the intersection of s and t. If s and t are non-overlapping, raises ValueError.

>>> arb("10 +/- 8.001").intersection(arb("0 +/- 2.001"))
[2.00 +/- 1.01e-3]
>>> arb(2).intersection(3)
Traceback (most recent call last):
  ...
ValueError: empty intersection
is_exact(self) bool
is_finite(self) bool
is_integer(self) bool
is_nan(self) bool
is_zero(self) bool
jacobi_p(s, n, a, b)

Jacobi polynomial (or Jacobi function) \(P_n^{a,b}(s)\).

>>> from flint import showgood
>>> showgood(lambda: (arb(1)/3).jacobi_p(5, 0.25, 0.5), dps=25)
0.4131944444444444444444444
laguerre_l(s, n, m=0)

Laguerre function \(L_n^{m}(s)\).

>>> from flint import showgood
>>> showgood(lambda: (arb(1)/3).laguerre_l(5, 0.25), dps=25)
0.03871323490012002743484225
lambertw(s, int branch=0)

Lambert W function, \(W_k(s)\). Either the principal branch (\(k = 0\)) or the \(k = -1\) branch can be computed by specifying branch.

>>> from flint import showgood
>>> showgood(lambda: arb(10).lambertw(), dps=25)
1.745528002740699383074301
>>> showgood(lambda: arb("-0.1").lambertw(-1), dps=25)
-3.577152063957297218409392
legendre_p(s, n, m=0, int type=2)

Legendre function of the first kind \(P_n^m(z)\).

>>> from flint import showgood
>>> showgood(lambda: (arb(1)/3).legendre_p(5), dps=25)
0.3333333333333333333333333
>>> showgood(lambda: (arb(1)/3).legendre_p(5, 1.5), dps=25)
-2.372124991643971726805456
>>> showgood(lambda: (arb(3)).legendre_p(5, 1.5, type=3), dps=25)
17099.70021476473458984981

The optional parameter type can be 2 or 3, and selects between two different branch cut conventions (see Mathematica and mpmath).

static legendre_p_root(ulong n, ulong k, bool weight=False)

Returns the index-k zero of the Legendre polynomial \(P_n(x)\). The zeros are indexed in decreasing order.

If weight is True, returns a tuple (x, w) where x is the zero and w is the corresponding weight for Gauss-Legendre quadrature on \((-1,1)\).

>>> from flint import showgood
>>> for k in range(5):
...     showgood(lambda: arb.legendre_p_root(5,k), dps=25)
...
0.9061798459386639927976269
0.5384693101056830910363144
0
-0.5384693101056830910363144
-0.9061798459386639927976269
>>> for k in range(3):
...     showgood(lambda: arb.legendre_p_root(3,k,weight=True), dps=15)
...
(0.774596669241483, 0.555555555555556)
(0, 0.888888888888889)
(-0.774596669241483, 0.555555555555556)
legendre_q(s, n, m=0, int type=2)

Legendre function of the second kind \(Q_n^m(z)\).

>>> from flint import showgood
>>> showgood(lambda: (arb(1)/3).legendre_q(5), dps=25)
0.1655245300933242182362054
>>> showgood(lambda: (arb(1)/3).legendre_q(5, 1.5), dps=25)
-6.059967350218583975575616

The optional parameter type can be 2 or 3, and selects between two different branch cut conventions (see Mathematica and mpmath).

lgamma(s)

Logarithmic gamma function \(\log \Gamma(s)\).

>>> from flint import showgood
>>> showgood(lambda: arb(100).lgamma(), dps=25)
359.1342053695753987760440

This function is undefined for negative \(s\). Use acb.lgamma() for the complex extension.

li(s, bool offset=False)

Logarithmic integral \(\operatorname{li}(s)\), optionally the offset logarithmic integral \(\operatorname{Li}(s) = \operatorname{li}(s) - \operatorname{li}(2)\).

>>> from flint import showgood
>>> showgood(lambda: arb(10).li(), dps=25)
6.165599504787297937522982
>>> showgood(lambda: arb(10).li(offset=True), dps=25)
5.120435724669805152678393
log(s)

Natural logarithm \(\log(s)\).

>>> from flint import showgood
>>> showgood(lambda: arb(2).log(), dps=25)
0.6931471805599453094172321
>>> showgood(lambda: arb(100).exp().log(), dps=25)
100.0000000000000000000000
>>> showgood(lambda: arb(-1).sqrt(), dps=25)
Traceback (most recent call last):
  ...
ValueError: no convergence (maxprec=960, try higher maxprec)

This function is undefined for negative input. Use acb.log() for the complex extension.

log1p(s)

Natural logarithm \(\log(1+s)\), computed accurately for small s.

>>> from flint import showgood
>>> from flint import acb
>>> showgood(lambda: acb(1).log1p(), dps=25)
0.6931471805599453094172321
>>> showgood(lambda: arb("1e-100000000000000000").log1p(), dps=25)
1.000000000000000000000000e-100000000000000000

This function is undefined for \(s \le -1\). Use acb.log1p() for the complex extension.

log_base(s, ulong b)

Returns \(\log_b(s)\), computed exactly when possible.

>>> from flint import arb, ctx
>>> ctx.prec = 53
>>> arb(2048).log_base(2)
11.0000000000000
lower(self)

Lower bound for self (towards \(-\infty\)). The output is an arb holding an exact floating-point number that has been rounded down to the current precision.

>>> print(arb("-5 +/- 2").lower().str(5, radius=False))
-7.0000
man_exp(self)

Decomposes self into an integer mantissa and an exponent, returning an fmpz pair. Requires that self is exact and finite.

>>> from flint import arb, ctx
>>> ctx.prec = 53
>>> arb("1.1").mid().man_exp()
(4953959590107545, -52)
>>> arb("1.1").rad().man_exp()
(1, -52)
>>> arb(0).man_exp()
(0, 0)
>>> arb("1.1").man_exp()
Traceback (most recent call last):
  ...
ValueError: man_exp requires an exact, finite value
>>> arb("+inf").man_exp()
Traceback (most recent call last):
  ...
ValueError: man_exp requires an exact, finite value
max(s, t)

Maximum value of s and t.

>>> arb(2).max(arb("3 +/- 1.1"))
[+/- 4.11]
>>> arb(4).max(arb("3 +/- 1.1"))
[4e+0 +/- 0.101]
mid(self)

Returns the midpoint of self as an exact arb:

>>> from flint import arb, ctx
>>> ctx.prec = 53
>>> arb("1 +/- 0.3").mid()
1.00000000000000
mid_rad_10exp(self, long n=0)

Returns an fmpz triple (mid, rad, exp) where the larger of mid and rad has n digits plus a few digits (n defaults to the current precision), such that self is contained in \([\operatorname{mid} \pm \operatorname{rad}] 10^{\operatorname{exp}}\).

>>> (arb(1) / 3).mid_rad_10exp(10)
(333333333333333, 2, -15)
>>> (arb(1) / 3).mid_rad_10exp(20)
(3333333333333333148296162, 555111516, -25)
>>> arb(0,1e-100).mid_rad_10exp(10)
(0, 100000000376832, -114)
>>> arb(-1,1e100).mid_rad_10exp()
(0, 10000000083585662976, 81)
min(s, t)

Minimum value of s and t.

>>> from flint import arb, ctx
>>> ctx.prec = 53
>>> arb(2).min(3)
2.00000000000000
>>> arb(2).min(arb("3 +/- 1.1"))
[2e+0 +/- 0.101]
static nan()
neg(self, bool exact=False)
static neg_inf()
nonnegative_part(self)
overlaps(self, other)
static partitions_p(n)

Number of partitions of the integer n, evaluated as an arb.

>>> from flint import showgood
>>> showgood(lambda: arb.partitions_p(10), dps=25)
42.00000000000000000000000
>>> showgood(lambda: arb.partitions_p(100), dps=25)
190569292.0000000000000000
>>> showgood(lambda: arb.partitions_p(10**50), dps=25)
3.285979358867807890529967e+11140086280105007830283557
static pi()

Returns the constant \(\pi\) as an arb.

>>> from flint import showgood
>>> showgood(arb.pi, dps=25)
3.141592653589793238462643
polylog(self, s)

Polylogarithm \(\operatorname{Li}_s(z)\) where the argument z is given by self and the order s is given as an extra parameter.

>>> from flint import showgood
>>> showgood(lambda: arb(-1).polylog(2), dps=25)
-0.8224670334241132182362076
>>> showgood(lambda: arb(-3).polylog(1.75), dps=25)
-1.813689945878060161612620
>>> showgood(lambda: arb(-2.5).polylog(4.75), dps=25)
-2.322090601785704585092044

This function is undefined for some combinations of \(s, z\). Use acb.polylog() for the complex extension.

static pos_inf()
rad(self)

Returns the radius of self as an exact arb:

>>> print(arb("1 +/- 0.3").rad().str(5, radius=False))
0.30000
real
rel_accuracy_bits(self)
rel_one_accuracy_bits(self)
repr(self)
rgamma(s)

Reciprocal gamma function \(1/\Gamma(s)\), avoiding division by zero at the poles of the gamma function.

>>> from flint import showgood
>>> showgood(lambda: arb(1.5).rgamma(), dps=25)
1.128379167095512573896159
>>> print(arb(0).rgamma())
0
>>> print(arb(-1).rgamma())
0
>>> print(arb(-3,1e-10).rgamma())
[+/- 6.01e-10]
rising(s, n)

Rising factorial \((s)_n\).

>>> from flint import showgood
>>> showgood(lambda: arb.pi().rising(0), dps=25)
1.000000000000000000000000
>>> showgood(lambda: arb.pi().rising(10), dps=25)
299606572.3661012684972888
>>> showgood(lambda: arb.pi().rising(0.5), dps=25)
1.703592785410167015590330
rising2(s, ulong n)

Computes the rising factorial \((s)_n\) where n is an unsigned integer, along with the first derivative with respect to \((s)_n\).

The current implementation does not use the gamma function, so n should be moderate.

>>> from flint import arb, ctx
>>> ctx.prec = 53
>>> u, v = arb(3).rising2(5)
>>> print(u); print(v)
2520.00000000000
2754.00000000000
static rising_fmpq_ui(fmpq s, ulong n)

Computes the rising factorial \((s)_n\) where s is a rational number and n is an unsigned integer. The current implementation does not use the gamma function, so n should be moderate.

>>> from flint import fmpq
>>> from flint import showgood
>>> showgood(lambda: arb.rising_fmpq_ui(fmpq(-1,3), 100), dps=25)
-4.960517984074284420131903e+154
root(s, ulong n)

Principal n-th root of s.

>>> from flint import showgood
>>> showgood(lambda: arb(10).root(3), dps=25)
2.154434690031883721759294
rsqrt(s)

Reciprocal square root \(1/\sqrt{s}\).

>>> from flint import showgood
>>> showgood(lambda: arb(3).rsqrt(), dps=25)
0.5773502691896257645091488
>>> showgood(lambda: arb(0).rsqrt(), dps=25)
Traceback (most recent call last):
  ...
ValueError: no convergence (maxprec=960, try higher maxprec)
>>> showgood(lambda: arb(-1).rsqrt(), dps=25)
Traceback (most recent call last):
  ...
ValueError: no convergence (maxprec=960, try higher maxprec)

This function is undefined for negative input. Use acb.rsqrt() for the complex extension.

sec(s)

Secant function \(\operatorname{sec}(s)\).

>>> from flint import showgood
>>> showgood(lambda: arb(1).sec(), dps=25)
1.850815717680925617911753
sech(s)

Hyperbolic secant \(\operatorname{sech}(s)\).

>>> from flint import showgood
>>> showgood(lambda: arb(1).sech(), dps=25)
0.6480542736638853995749774
sgn(self)

Sign function, returning an arb.

>>> from flint import arb, ctx
>>> ctx.prec = 53
>>> arb(-3).sgn()
-1.00000000000000
>>> arb(0).sgn()
0
>>> arb("0 +/- 1").sgn()
[+/- 1.01]
shi(s)

Hyperbolic sine integral \(\operatorname{Shi}(s)\).

>>> from flint import showgood
>>> showgood(lambda: arb(3).shi(), dps=25)
4.973440475859806797710418
si(s)

Sine integral \(\operatorname{Si}(s)\).

>>> from flint import showgood
>>> showgood(lambda: arb(3).si(), dps=25)
1.848652527999468256397730
sin(s)

Sine function \(\sin(s)\).

>>> from flint import showgood
>>> showgood(lambda: arb(1).sin(), dps=25)
0.8414709848078965066525023
sin_cos(s)

Computes \(\sin(s)\) and \(\cos(s)\) simultaneously.

>>> from flint import showgood
>>> showgood(lambda: arb(1).sin_cos(), dps=25)
(0.8414709848078965066525023, 0.5403023058681397174009366)
sin_cos_pi(s)

Computes \(\sin(\pi s)\) and \(\cos(\pi s)\) simultaneously.

>>> from flint import showgood
>>> showgood(lambda: arb(0.75).sin_cos_pi(), dps=25)
(0.7071067811865475244008444, -0.7071067811865475244008444)
static sin_cos_pi_fmpq(fmpq s)

Computes \(\sin(\pi s)\) and \(\cos(\pi s)\) simultaneously.

>>> from flint import fmpq
>>> from flint import showgood
>>> showgood(lambda: arb.sin_cos_pi_fmpq(fmpq(3,4)), dps=25)
(0.7071067811865475244008444, -0.7071067811865475244008444)
sin_pi(s)

Sine function \(\sin(\pi s)\).

>>> from flint import showgood
>>> showgood(lambda: arb(0.75).sin_pi(), dps=25)
0.7071067811865475244008444
>>> showgood(lambda: arb(1).sin_pi(), dps=25)
0
static sin_pi_fmpq(fmpq s)

Returns the algebraic sine value \(\sin(\pi s)\).

>>> from flint import fmpq
>>> from flint import showgood
>>> showgood(lambda: arb.sin_pi_fmpq(fmpq(3,4)), dps=25)
0.7071067811865475244008444
sinc(s)

Sinc function, \(\operatorname{sinc}(x) = \sin(x)/x\).

>>> from flint import showgood
>>> showgood(lambda: arb(3).sinc(), dps=25)
0.04704000268662240736691493
sinc_pi(s)

Normalized sinc function, \(\operatorname{sinc}(\pi x) = \sin(\pi x)/(\pi x)\).

>>> from flint import showgood
>>> showgood(lambda: arb(1.5).sinc_pi(), dps=25)
-0.2122065907891937810251784
>>> showgood(lambda: arb(2).sinc_pi(), dps=25)
0
sinh(s)

Hyperbolic sine \(\sinh(s)\).

>>> from flint import showgood
>>> showgood(lambda: arb(1).sinh(), dps=25)
1.175201193643801456882382
sinh_cosh(s)

Computes \(\sinh(s)\) and \(\cosh(s)\) simultaneously.

>>> from flint import showgood
>>> showgood(lambda: arb(1).sinh_cosh(), dps=25)
(1.175201193643801456882382, 1.543080634815243778477906)
sqrt(s)

Square root \(\sqrt{s}\).

>>> from flint import showgood
>>> showgood(lambda: arb(3).sqrt(), dps=25)
1.732050807568877293527446
>>> showgood(lambda: arb(0).sqrt(), dps=25)
0
>>> showgood(lambda: arb(-1).sqrt(), dps=25)
Traceback (most recent call last):
  ...
ValueError: no convergence (maxprec=960, try higher maxprec)

This function is undefined for negative input. Use acb.sqrt() for the complex extension.

str(self, long n=0, bool radius=True, bool more=False, long condense=0)

Produces a human-readable decimal representation of self, with up to n printed digits (which defaults to the current precision) for the midpoint. The output can be parsed by the arb constructor.

Since the internal representation is binary, conversion to decimal (and back from decimal) is generally inexact. Binary-decimal-binary roundtrips may result in significantly larger intervals, and should therefore be done sparingly.

>>> from flint import arb, ctx
>>> ctx.prec = 53
>>> print(arb.pi().str())
[3.14159265358979 +/- 3.34e-15]
>>> print(arb.pi().str(5))
[3.1416 +/- 7.35e-6]
>>> print(arb.pi().str(5, radius=False))
3.1416

By default, the output is truncated so that all displayed digits are guaranteed to be correct, up to adding or subtracting 1 in the last displayed digit (as a special case, if the output ends with a string of 0s, the correct decimal expansion to infinite precision could have a string of 9s).

>>> print((arb(1) - arb("1e-10")).str(5))
[1.0000 +/- 4e-10]
>>> print((arb(1) - arb("1e-10")).str(10))
[0.9999999999 +/- 3e-15]

To force more digits, set more to True.

>>> print(arb("0.1").str(30))
[0.100000000000000 +/- 2.23e-17]
>>> print(arb("0.1").str(30, more=True))
[0.0999999999999999916733273153113 +/- 1.39e-17]

Note that setting more to True results in a smaller printed radius, since there is less error from the conversion back to decimal.

>>> x = arb.pi().sin()
>>> print(x.str())
[+/- 3.46e-16]
>>> print(x.str(more=True))
[1.22460635382238e-16 +/- 2.23e-16]

The error indicated in the output may be much larger than the actual error in the internal representation of self. For example, if self is accurate to 1000 digits and printing is done at 10-digit precision, the output might only appear as being accurate to 10 digits. It is even possible for self to be exact and have an inexact decimal representation.

The condense option can be passed to show only leading and trailing digits of the fractional, integer and exponent parts of the output.

>>> ctx.dps = 1000
>>> print(arb.pi().str(condense=10))
[3.1415926535{...979 digits...}9216420199 +/- 9.28e-1001]
>>> print(arb.fac_ui(300).str(condense=10))
3060575122{...595 digits...}0000000000.0000000000{...365 digits...}0000000000
>>> print(arb(10**100).exp().str(condense=5))
[1.53837{...989 digits...}96534e+43429{...90 digits...}17483 +/- 4.84e+43429{...90 digits...}16483]
>>> ctx.default()
tan(s)

Tangent function \(\tan(s)\).

>>> from flint import showgood
>>> showgood(lambda: arb(1).tan(), dps=25)
1.557407724654902230506975
tan_pi(s)

Tangent function \(\tan(\pi s)\).

>>> from flint import showgood
>>> showgood(lambda: arb(0.125).tan_pi(), dps=25)
0.4142135623730950488016887
>>> showgood(lambda: arb(1).tan_pi(), dps=25)
0
tanh(s)

Hyperbolic tangent \(\tanh(s)\).

>>> from flint import showgood
>>> showgood(lambda: arb(1).tanh(), dps=25)
0.7615941559557648881194583
union(s, t)

Returns a ball containing the union of s and t.

>>> from flint import arb, ctx
>>> ctx.prec = 53
>>> x = arb(3).union(5); x.lower(); x.upper()
[2.99999999813735 +/- 4.86e-15]
[5.00000000186265 +/- 4.86e-15]
unique_fmpz(self)

If self represents exactly one integer, returns this value as an fmpz; otherwise returns None.

>>> arb("5 +/- 0.1").unique_fmpz()
5
>>> arb("5 +/- 0.9").unique_fmpz()
5
>>> arb("5.1 +/- 0.9").unique_fmpz()
>>>
upper(self)

Upper bound for self (towards \(+\infty\)). The output is an arb holding an exact floating-point number that has been rounded up to the current precision.

>>> print(arb("-5 +/- 2").upper().str(5, radius=False))
-3.0000
zeta(s, a=None)

Riemann zeta function \(\zeta(s)\) or the Hurwitz zeta function \(\zeta(s,a)\) if a second parameter is passed.

>>> from flint import showgood
>>> showgood(lambda: arb(4.25).zeta(), dps=25)
1.066954190711214532450370
>>> showgood(lambda: arb(4.25).zeta(2.75), dps=25)
0.01991885526414599096374229

This function is undefined for some combinations of \(s, a\). Use acb.zeta() for the complex extension.

zeta_nzeros(x)

Number of zeros of the Riemann zeta function with positive imaginary part between 0 and x.

>>> from flint import arb, ctx
>>> ctx.prec = 53
>>> arb("-5").zeta_nzeros()
0
>>> arb("14").zeta_nzeros()
0
>>> arb("15").zeta_nzeros()
1.00000000000000
>>> arb("14.1 +/- 0.1").zeta_nzeros()
[+/- 1.01]
>>> arb("100").zeta_nzeros()
29.0000000000000
>>> arb("1e6").zeta_nzeros()
1747146.00000000