acb – complex numbers¶
-
class
flint.
acb
(real=None, imag=None)¶ An acb represents a complex number by a rectangular enclosure consisting of arb balls for the real and imaginary parts.
>>> acb(2) 2.00000000000000 >>> acb(2+3j) 2.00000000000000 + 3.00000000000000j >>> acb("2 +/- 0.001", fmpq(2,3)) [2.00 +/- 1.01e-3] + [0.666666666666667 +/- 4.82e-16]j >>> acb(-1) ** 0.25 [0.707106781186547 +/- 6.14e-16] + [0.707106781186547 +/- 6.15e-16]j
-
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(acb(3, "-5 +/- 2").abs_lower().str(5, radius=False)) 4.2426
-
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(acb(3, "-5 +/- 2").abs_upper().str(5, radius=False)) 7.6158
-
acos
(s)¶ Inverse cosine \(\operatorname{acos}(s)\).
>>> showgood(lambda: acb(2).acos(), dps=25) 1.316957896924816708625046j
-
acosh
(s)¶ Inverse hyperbolic cosine \(\operatorname{acosh}(s)\).
>>> showgood(lambda: acb(2,3).acosh(), dps=25) 1.983387029916535432347077 + 1.000143542473797218521038j
-
agm
(s, t=None)¶ Arithmetic-geometric mean \(M(s,t)\), or \(M(s) = M(s,1)\) if no extra parameter is passed.
>>> showgood(lambda: acb(2).agm(), dps=25) 1.456791031046906869186432 >>> showgood(lambda: acb(1,1).agm(), dps=25) 1.049160528732780220531827 + 0.4781557460881612293261882j >>> showgood(lambda: (acb(-95,-65)/100).agm(acb(684,747)/1000), dps=25) -0.3711072435676023931065922 + 0.3199561471173686568561674j
-
airy
(s)¶ Computes the Airy function values \(\operatorname{Ai}(s)\), \(\operatorname{Ai}'(s)\), \(\operatorname{Bi}(s)\), \(\operatorname{Bi}'(s)\) simultaneously, returning a tuple.
>>> showgood(lambda: acb(-1+1j).airy(), dps=5) (0.82212 - 0.11997j, -0.37906 - 0.60450j, 0.21429 + 0.67392j, 0.83447 - 0.34653j)
-
airy_ai
(s, int derivative=0)¶ Airy function \(\operatorname{Ai}(s)\), or \(\operatorname{Ai}'(s)\) if derivative is 1.
>>> showgood(lambda: acb(-1+1j).airy_ai(), dps=25) 0.8221174265552725939610246 - 0.1199663426644243438939006j >>> showgood(lambda: acb(-1+1j).airy_ai(derivative=1), dps=25) -0.3790604792268334962367164 - 0.6045001308622460716372591j
-
airy_bi
(s, int derivative=0)¶ Airy function \(\operatorname{Bi}(s)\), or \(\operatorname{Bi}'(s)\) if derivative is 1.
>>> showgood(lambda: acb(-1+1j).airy_bi(), dps=25) 0.2142904015348735739780868 + 0.6739169237227052095951775j >>> showgood(lambda: acb(-1+1j).airy_bi(derivative=1), dps=25) 0.8344734885227826369040951 - 0.3465260632668285285537957j
-
arg
(self)¶ Complex argument (phase).
>>> showgood(lambda: acb("3.3").arg(), dps=25) 0 >>> showgood(lambda: acb(-1).arg(), dps=25) 3.141592653589793238462643 >>> acb(-1, "+/- 1").arg() [+/- 3.15]
-
asin
(s)¶ Inverse sine \(\operatorname{asin}(s)\).
>>> showgood(lambda: acb(2).asin(), dps=25) 1.570796326794896619231322 - 1.316957896924816708625046j
-
asinh
(s)¶ Inverse hyperbolic sine \(\operatorname{asinh}(s)\).
>>> showgood(lambda: acb(2,3).asinh(), dps=25) 1.968637925793096291788665 + 0.9646585044076027920454111j
-
atan
(s)¶ Computes the inverse tangent \(\operatorname{atan}(s)\).
>>> showgood(lambda: acb(1,2).atan(), dps=25) 1.338972522294493561124194 + 0.4023594781085250936501898j
-
atanh
(s)¶ Inverse hyperbolic tangent \(\operatorname{atanh}(s)\).
>>> showgood(lambda: acb(2,3).atanh(), dps=25) 0.1469466662255297520474328 + 1.338972522294493561124194j
-
barnes_g
(s)¶ Barnes G-function \(G(s)\).
>>> acb(8).barnes_g() 24883200.0000000 >>> showgood(lambda: acb(2+3j).barnes_g(), dps=25) -0.1781021386408216960641890 + 0.04504542715447837909120582j
-
bernoulli_poly
(s, ulong n)¶ Returns the value of the Bernoulli polynomial \(B_n(s)\).
>>> showgood(lambda: acb(0.25+0.25j).bernoulli_poly(5), dps=25) -0.05859375000000000000000000 + 0.006510416666666666666666667j
-
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.
>>> showgood(lambda: acb(5).bessel_i(1), dps=25) 24.33564214245052719914305 >>> showgood(lambda: acb(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.
>>> showgood(lambda: acb(5).bessel_j(1), dps=25) -0.3275791375914652220377343 >>> showgood(lambda: acb(2+3j).bessel_j(1+2j), dps=25) 0.5041904509946947234759103 - 0.1765180072689450645147231j
-
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.
>>> showgood(lambda: acb(5).bessel_k(1), dps=25) 0.004044613445452164208365022 >>> showgood(lambda: acb(5).bessel_k(1, scaled=True), dps=25) 0.6002738587883125829360457 >>> showgood(lambda: acb(2+3j).bessel_k(1+2j), dps=25) -0.09884736370006798963642719 - 0.02870616366668971734065520j
-
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.
>>> showgood(lambda: acb(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.
>>> showgood(lambda: acb(2,3).beta_lower(1, 2.5), dps=25) 0.2650137734913866999568502 - 7.111836702381954625752124j >>> showgood(lambda: acb(2,3).beta_lower(1, 2.5, regularized=True), dps=25) 0.6625344337284667498921254 - 17.77959175595488656438031j
-
chebyshev_t
(s, n)¶ Chebyshev function of the first kind \(T_n(s)\).
>>> showgood(lambda: (acb(1)/3).chebyshev_t(3), dps=25) -0.8518518518518518518518519
-
chebyshev_u
(s, n)¶ Chebyshev function of the second kind \(U_n(s)\).
>>> showgood(lambda: (acb(1)/3).chebyshev_u(3), dps=25) -1.037037037037037037037037
-
chi
(s)¶ Hyperbolic cosine integral \(\operatorname{Chi}(s)\).
>>> showgood(lambda: acb(10).chi(), dps=25) 1246.114486042454414726558
-
ci
(s)¶ Cosine integral \(\operatorname{Ci}(s)\).
>>> showgood(lambda: acb(10).ci(), dps=25) -0.04545643300445537263453283
-
complex_rad
(self)¶ Returns an acb representing the radii of the real and imaginary parts of self together a single complex number.
>>> print(acb("1 +/- 0.3", "2 +/- 0.4").complex_rad().str(5, radius=False)) 0.30000 + 0.40000j
-
contains
(self, other)¶
-
cos
(s)¶ Cosine function \(\cos(s)\).
>>> showgood(lambda: acb(1,2).cos(), dps=25) 2.032723007019665529436343 - 3.051897799151800057512116j
-
cos_pi
(s)¶ Cosine function \(\cos(\pi s)\).
>>> showgood(lambda: acb(1,2).cos_pi(), dps=25) -267.7467614837482222459319
-
cosh
(s)¶ Hyperbolic cosine function \(\cosh(s)\).
>>> showgood(lambda: acb(2,3).cosh(), dps=25) -3.724545504915322565473971 + 0.5118225699873846088344638j
-
cot
(s)¶ Cotangent function \(\cot(s)\).
>>> showgood(lambda: acb(1,2).cot(), dps=25) 0.03279775553375259406276455 - 0.9843292264581910294718882j
-
cot_pi
(s)¶ Cotangent function \(\cot(\pi s)\).
>>> showgood(lambda: acb(1,2).cot_pi(), dps=25) -1.000006974709035616233122j
-
coth
(s)¶ Hyperbolic cotangent function \(\coth(s)\).
>>> showgood(lambda: acb(2,3).coth(), dps=25) 1.035746637764995396112759 + 0.01060478347033710175031690j
-
csc
(s)¶ Cosecant function \(\sec(s)\).
>>> showgood(lambda: acb(2,3).csc(), dps=25) 0.09047320975320743980579048 + 0.04120098628857412646300981j
-
csch
(s)¶ Hyperbolic cosecant function \(\operatorname{csch}(s)\).
>>> showgood(lambda: acb(2,3).csch(), dps=25) -0.2725486614629401995124985 - 0.04030057885689152187513248j
-
csgn
(self)¶ Complex sign function defined as a piecewise extension of the real sign function.
>>> showgood(lambda: acb(2,3).csgn(), dps=25) 1.000000000000000000000000 >>> showgood(lambda: acb(-1).csgn(), dps=25) -1.000000000000000000000000
-
static
dft
(vec, bool inverse=False)¶ Returns the discrete Fourier transform (DFT) of a given iterable vec. The output is a list of acb entries. If inverse is True, computes the inverse transform instead.
>>> for c in acb.dft(acb.dft(range(1,12)), inverse=True): ... print(c) ... [1.000000000000 +/- 1.06e-13] + [+/- 8.66e-14]j [2.000000000000 +/- 1.25e-13] + [+/- 1.17e-13]j [3.000000000000 +/- 1.29e-13] + [+/- 1.21e-13]j [4.000000000000 +/- 1.34e-13] + [+/- 1.22e-13]j [5.000000000000 +/- 1.40e-13] + [+/- 1.24e-13]j [6.000000000000 +/- 1.50e-13] + [+/- 1.26e-13]j [7.000000000000 +/- 1.43e-13] + [+/- 1.26e-13]j [8.000000000000 +/- 1.47e-13] + [+/- 1.23e-13]j [9.000000000000 +/- 1.47e-13] + [+/- 1.22e-13]j [10.00000000000 +/- 1.50e-13] + [+/- 1.20e-13]j [11.000000000000 +/- 1.43e-13] + [+/- 1.22e-13]j >>> sum(acb.dft(acb.dft(range(1,10001)), inverse=True)).contains(50005000) True
-
digamma
(s)¶ Digamma function \(\psi(s)\).
>>> showgood(lambda: acb(1,2).digamma(), dps=25) 0.7145915153739775266568699 + 1.320807282642230228386088j
-
dirichlet_eta
(s)¶ Dirichlet eta function \(\eta(s)\).
>>> showgood(lambda: acb(1).dirichlet_eta(), dps=25) 0.6931471805599453094172321 >>> showgood(lambda: acb(0).dirichlet_eta(), dps=25) 0.5000000000000000000000000 >>> showgood(lambda: acb(0.5+10000j).dirichlet_eta(), dps=25) -0.08240476492345768707759740 - 0.4458305781469329610368705j
-
ei
(s)¶ Exponential integral \(\operatorname{Ei}(s)\).
>>> showgood(lambda: acb(10).ei(), dps=25) 2492.228976241877759138440
-
elliptic_e
(m)¶ - Complete elliptic integral of the second kind \(E(m)\).
>>> showgood(lambda: 2 * acb(0).elliptic_e(), dps=25) 3.141592653589793238462643 >>> showgood(lambda: acb(100+50j).elliptic_e(), dps=25) 2.537235535915583230880553 - 10.10997759792420947130321j >>> showgood(lambda: (1 - acb("1e-100")).elliptic_e(), dps=25) 1.000000000000000000000000
-
static
elliptic_e_inc
(phi, m, bool pi=False)¶ Incomplete elliptic integral \(E(\phi,m)\).
>>> showgood(lambda: acb.elliptic_e_inc(2, 0.75), dps=25) 1.443433069099461566497882 >>> showgood(lambda: acb.elliptic_e_inc(2, 0.75, pi=True), dps=25) 4.844224110273838099214252
-
static
elliptic_f
(phi, m, bool pi=False)¶ Incomplete elliptic integral \(F(\phi,m)\).
>>> showgood(lambda: acb.elliptic_f(2, 0.75), dps=25) 2.952569673655779502336268 >>> showgood(lambda: acb.elliptic_f(2, 0.75, pi=True), dps=25) 8.626062589998572941754700
-
elliptic_inv_p
(self, tau)¶ Inverse Weierstrass elliptic function \(\sigma^{-1}(z, \tau)\) where z is given by self.
>>> showgood(lambda: acb("1/3","1/5").elliptic_p(1j).elliptic_inv_p(1j), dps=25) 0.3333333333333333333333333 + 0.2000000000000000000000000j
-
elliptic_invariants
(tau)¶ Returns the lattice invariants \(g_2, g_3\) given the lattice parameter \(\tau\).
>>> showgood(lambda: acb(0.5+1j).elliptic_invariants(), dps=25) (72.64157667926127619414883, 536.6642788346023116199232)
-
elliptic_k
(m)¶ Complete elliptic integral of the first kind \(K(m)\).
>>> showgood(lambda: 2 * acb(0).elliptic_k(), dps=25) 3.141592653589793238462643 >>> showgood(lambda: acb(100+50j).elliptic_k(), dps=25) 0.2052037361984861505113972 + 0.3158446040520529200980591j >>> showgood(lambda: (1 - acb("1e-100")).elliptic_k(), dps=25) 116.5155490108221748197340
-
elliptic_p
(self, tau)¶ Weierstrass elliptic function \(\wp(z, \tau)\) where z is given by self.
>>> showgood(lambda: acb("1/3","1/5").elliptic_p(1j), dps=25) 3.686380646078879780287811 - 4.591498371497259422829963j >>> showgood(lambda: acb("1/3","6/5").elliptic_p(1j), dps=25) 3.686380646078879780287811 - 4.591498371497259422829963j
-
static
elliptic_pi
(n, m)¶ Complete elliptic integral \(\Pi(n,m)\).
>>> showgood(lambda: acb.elliptic_pi(0.25, 0.125), dps=25) 1.879349451879603826415650
-
static
elliptic_pi_inc
(n, phi, m, bool pi=False)¶ Incomplete elliptic integral \(\Pi(n,\phi,m)\).
>>> showgood(lambda: acb.elliptic_pi_inc(0.25, 0.5, 0.125), dps=25) 0.5128718023282086251399279 >>> showgood(lambda: acb.elliptic_pi_inc(0.25, 0.5, 0.125, pi=True), dps=25) 1.879349451879603826415650
-
static
elliptic_rc
(x, y)¶ Carlson incomplete elliptic integral \(R_C(x,y)\).
>>> showgood(lambda: acb.elliptic_rc(1, 2+3j), dps=25) 0.5952169239306156198937085 - 0.2387981909090509407085918j
-
static
elliptic_rd
(x, y, z)¶ Carlson incomplete elliptic integral \(R_D(x,y,z)\).
>>> showgood(lambda: acb.elliptic_rd(1, 2, 1+2j), dps=25) 0.2043722510302629773408686 - 0.3559745898273715996616328j
-
static
elliptic_rf
(x, y, z)¶ Carlson incomplete elliptic integral \(R_F(x,y,z)\).
>>> showgood(lambda: acb.elliptic_rf(1, 2+3j, 3+4j), dps=25) 0.5577655465453921610327459 - 0.2202042457195556054465308j
-
static
elliptic_rg
(x, y, z)¶ Carlson incomplete elliptic integral \(R_G(x,y,z)\).
>>> showgood(lambda: acb.elliptic_rg(1, 2, 1+2j), dps=25) 1.206557168056722980475052 + 0.2752176688707739710008275j
-
static
elliptic_rj
(x, y, z, p)¶ Carlson incomplete elliptic integral \(R_J(x,y,z,p)\).
>>> showgood(lambda: acb.elliptic_rj(1, 2, 1+2j, 2+3j), dps=25) 0.1604659632144333057202530 - 0.2502751672723324201692394j
-
elliptic_roots
(tau)¶ Returns the elliptic roots \(e_1, e_2, e_3\) given the lattice parameter \(\tau\).
>>> showgood(lambda: acb(0.5+1j).elliptic_roots(), dps=10) (6.285388119, -3.142694059 + 3.386618325j, -3.142694059 - 3.386618325j)
-
elliptic_sigma
(self, tau)¶ Weierstrass elliptic function \(\sigma(z, \tau)\) where z is given by self.
>>> showgood(lambda: acb("1/3","1/5").elliptic_sigma(1j), dps=25) 0.3396549497136886526515370 + 0.1970690762350931272896772j
-
elliptic_zeta
(self, tau)¶ Weierstrass elliptic function \(\zeta(z, \tau)\) where z is given by self.
>>> showgood(lambda: acb("1/3","1/5").elliptic_zeta(1j), dps=25) 2.219680339508418716159086 - 1.504947925755241700681002j
-
erf
(s)¶ Error function \(\operatorname{erf}(s)\).
>>> showgood(lambda: acb(2+3j).erf() - 1, dps=25) -21.82946142761456838910309 + 8.687318271470163144428079j >>> showgood(lambda: acb("77.7").erf() - 1, dps=25, maxdps=10000) -7.929310690520378873143053e-2625
-
erfc
(s)¶ Complementary error function \(\operatorname{erfc}(s)\).
>>> showgood(lambda: acb("77.7").erfc(), dps=25) 7.929310690520378873143053e-2625
-
erfi
(s)¶ Imaginary error function \(\operatorname{erfi}(s)\).
>>> showgood(lambda: acb(10).erfi(), dps=25) 1.524307422708669699360547e+42
-
exp
(s)¶ Exponential function \(\exp(s)\).
>>> showgood(lambda: acb(1,2).exp(), dps=25) -1.131204383756813638431255 + 2.471726672004818927616931j
-
exp_pi_i
(s)¶ Exponential function of modified argument \(\exp(\pi i s)\).
>>> showgood(lambda: acb(1,2).exp_pi_i(), dps=25) -0.001867442731707988814430213 >>> showgood(lambda: acb(1.5,2.5).exp_pi_i(), dps=25) -0.0003882032039267662472325299j >>> showgood(lambda: acb(1.25,2.25).exp_pi_i(), dps=25) -0.0006020578259597635239581705 - 0.0006020578259597635239581705j
-
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.
>>> showgood(lambda: acb(2+3j).expint(1+2j), dps=25) -0.01442661495527080336037156 + 0.01942348372986687164972794j
-
expm1
(s)¶ Exponential function \(\exp(s)-1\), computed accurately for small s.
>>> showgood(lambda: acb("1e-10000").expm1(), dps=25) 1.000000000000000000000000e-10000
-
fresnel_c
(s, bool normalized=True)¶ Fresnel cosine integral \(C(s)\), optionally not normalized.
>>> showgood(lambda: acb(3).fresnel_c(), dps=25) 0.6057207892976856295561611 >>> showgood(lambda: acb(3).fresnel_c(normalized=False), dps=25) 0.7028635577302687301744099
-
fresnel_s
(s, bool normalized=True)¶ Fresnel sine integral \(S(s)\), optionally not normalized.
>>> showgood(lambda: acb(3).fresnel_s(), dps=25) 0.4963129989673750360976123 >>> showgood(lambda: acb(3).fresnel_s(normalized=False), dps=25) 0.7735625268937690171497722
-
gamma
(s)¶ Gamma function \(\Gamma(s)\).
>>> showgood(lambda: acb(1,2).gamma(), dps=25) 0.1519040026700361374481610 + 0.01980488016185498197191013j
-
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.
>>> showgood(lambda: acb(2+3j).gamma_lower(2.5), dps=25) 1.646077010134876664349297 + 1.140585862703100904414519j >>> showgood(lambda: acb(2+3j).gamma_lower(2.5, regularized=1), dps=25) 1.238266003780709160156983 + 0.8580088838385611055576971j >>> showgood(lambda: acb(2+3j).gamma_lower(2.5, regularized=2), dps=25) -0.01687942194354244633506487 - 0.05864800341005793099108467j
-
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.
>>> showgood(lambda: acb(2+3j).gamma_upper(1+2j), dps=25) 0.02614303924198793235765248 - 0.0007536537278463329391666679j >>> showgood(lambda: acb(2+3j).gamma_upper(1+2j, regularized=True), dps=25) 0.1685897763996036749499208 - 0.02694171301624093921683609j
-
gegenbauer_c
(s, n, m)¶ Gegenbauer function \(C_n^{m}(s)\).
>>> showgood(lambda: (acb(1)/3).gegenbauer_c(5, 0.25), dps=25) 0.1321855709876543209876543
-
hermite_h
(s, n)¶ Hermite function \(H_n(s)\).
>>> showgood(lambda: (acb(1)/3).hermite_h(5), dps=25) 34.20576131687242798353909
-
hypgeom
(self, a, b, bool regularized=False, long n=-1)¶ 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.
>>> showgood(lambda: acb.pi().hypgeom([1+1j, 2-2j], [3, fmpq(1,3)]), dps=25) # 2F2 144.9760711583421645394627 - 51.06535684838559608699106j >>> showgood(lambda: acb.pi().hypgeom([1+1j, 2-2j], [3, fmpq(1,3)], regularized=True), dps=25) 27.05849150326959272369764 - 9.530893707861133993129862j
The optional parameter n, if nonnegative, controls the number of terms to add in the hypergeometric series. This is just a tuning parameter: a rigorous error bound is computed regardless of n.
-
hypgeom_0f1
(self, a, bool regularized=False)¶ Hypergeometric function \({}_0F_1(a,z)\) where the argument z is given by self Optionally the regularized version of this function can be computed.
>>> showgood(lambda: acb(-5).hypgeom_0f1(2.5), dps=25) 0.003114611044402738470826907 >>> showgood(lambda: acb(-5).hypgeom_0f1(2.5, regularized=True), dps=25) 0.002342974810739764377177885
-
hypgeom_1f1
(self, a, b, bool regularized=False)¶ Kummer’s confluent hypergeometric function \({}_1F_1(a,b,z)\) where z is given by self.. Optionally, computes the regularized version.
>>> showgood(lambda: acb(40000+50000j).hypgeom_1f1(2+3j, 3+4j), dps=25) 3.730925582634533963357515e+17366 + 3.199717318207534638202987e+17367j >>> showgood(lambda: acb(40000+50000j).hypgeom_1f1(2+3j, 3+4j) / acb(3+4j).gamma(), dps=25) -1.846160890579724375436801e+17368 + 2.721369772032882467996588e+17367j >>> showgood(lambda: acb(10).hypgeom_1f1(5, -3, regularized=True), dps=25) 832600407043.6938843410086 >>> showgood(lambda: acb(10).hypgeom_1f1(-5,-6), dps=25) 403.7777777777777777777778 >>> showgood(lambda: acb(10).hypgeom_1f1(-5,-5,regularized=True), dps=25) 0 >>> showgood(lambda: acb(10).hypgeom_1f1(-5,-6,regularized=True), dps=25) 0 >>> showgood(lambda: acb(10).hypgeom_1f1(-5,-4,regularized=True), dps=25) -100000.0000000000000000000
-
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.
>>> showgood(lambda: acb(5).hypgeom_2f1(1,2,3), dps=25) -0.5109035488895912495067571 - 0.2513274122871834590770115j >>> showgood(lambda: acb(5).hypgeom_2f1(1,2,3,regularized=True), dps=25) -0.2554517744447956247533786 - 0.1256637061435917295385057j
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.
>>> showgood(lambda: acb("11/10").hypgeom_2f1(arb(2).sqrt(), 0.5, arb(2).sqrt()+1.5, abc=True), dps=25) 1.801782659480054173351004 - 0.3114019850045849100659641j >>> showgood(lambda: acb("11/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, long n=-1, bool asymp=False)¶ Tricomi’s confluent hypergeometric function \(U(a,b,z)\) where z is given by self.
If asymp is set to True the asymptotic series is forced. If \(|z|\) is small, the attainable accuracy is then limited. The optional parameter n, if nonnegative, controls the number of terms to add in the asymptotic series. This is just a tuning parameter: a rigorous error bound is computed regardless of n.
>>> showgood(lambda: acb(400+500j).hypgeom_u(1+1j, 2+3j), dps=25) 0.001836433961463105354717547 - 0.003358699641979853540147122j >>> showgood(lambda: acb(-30).hypgeom_u(1+1j, 2+3j), dps=25) 0.7808944974399200669072087 - 0.2674783064947089569672470j >>> print(acb(-30).hypgeom_u(1+1j, 2+3j, n=0, asymp=True)) [+/- 3.41] + [+/- 3.41]j >>> print(acb(-30).hypgeom_u(1+1j, 2+3j, n=30, asymp=True)) [0.7808944974 +/- 7.46e-11] + [-0.2674783065 +/- 3.31e-11]j >>> print(acb(-30).hypgeom_u(1+1j, 2+3j, n=60, asymp=True)) [0.78089 +/- 8.14e-6] + [-0.26748 +/- 5.34e-6]j
-
imag
¶
-
static
integral
(func, a, b, params=None, rel_tol=None, abs_tol=None, deg_limit=None, eval_limit=None, depth_limit=None, use_heap=None, verbose=None)¶ Computes the integral \(\int_a^b f(x) dx\) where the integrand f is defined by func.
>>> showgood(lambda: acb.integral(lambda x, _: x.sin(), 0, arb.pi()), dps=25) 2.000000000000000000000000 >>> showgood(lambda: acb.integral(lambda x, _: (x + x.sin()).gamma(), 1, 1+1j), dps=25) -0.2732681890680958866139676 + 0.7064496061603478580993410j
The function func takes two parameters as input: the argument x and a boolean flag analytic. If analytic is False, func should return \(f(x)\), and in this case there are no restrictions on f; for instance, the integrand can be discontinuous on x. If analytic is True, func must verify that f is analytic on x, and if not, return a non-finite ball instead of \(f(x)\).
The analytic flag can be ignored for meromorphic functions, because evaluation at poles automatically leads to non-finite balls. However, it must be checked for functions that are non-analytic in some regions (for instance, functions with branch cuts such as \(\sqrt{x}\) and \(\log(x)\)), even if the integration path does not touch any non-analytic points. Some methods have an analytic option built-in, so the user simply has to forward this flag:
>>> showgood(lambda: acb.integral(lambda x, _: x.sqrt(), 1, 4), dps=25) # WRONG!!! 4.669414894781006338774348 >>> showgood(lambda: acb.integral(lambda x, a: x.sqrt(analytic=a), 1, 4), dps=25) # correct 4.666666666666666666666667
The following works without handling the analytic flag, because the integrand is meromorphic:
>>> showgood(lambda: acb.integral(lambda x, _: x.sech(), -1000, 1000), dps=25) 3.141592653589793238462643
The options rel_tol and abs_tol specify the relative and absolute tolerance goal for the integration. Both default to \(2^{-p}\) where p is the current precision.
The options deg_limit, eval_limit, depth_limit and use_heap allow control over the amount of work done before aborting; see the documentation for acb_calc_integrate for details.
-
jacobi_p
(s, n, a, b)¶ Jacobi polynomial (or Jacobi function) \(P_n^{a,b}(s)\).
>>> showgood(lambda: (acb(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)\).
>>> showgood(lambda: (acb(1)/3).laguerre_l(5, 0.25), dps=25) 0.03871323490012002743484225
-
lambertw
(s, branch=0, bool left=False, bool middle=False)¶ Lambert W function, \(W_k(s)\) where k is given by branch.
>>> showgood(lambda: acb(1).lambertw(), dps=25) 0.5671432904097838729999687 >>> showgood(lambda: acb(1).lambertw(-1), dps=25) -1.533913319793574507919741 - 4.375185153061898385470907j >>> showgood(lambda: acb(1).lambertw(-2), dps=25) -2.401585104868002884174140 - 10.77629951611507089849710j
The branch cuts follow Corless et al. by default. This function allows selecting alternative branch cuts in order to support continuous analytic continuation on complex intervals.
If left is set, computes \(W_{\mathrm{left}|k}(z)\) which corresponds to \(W_k(z)\) in the upper half plane and \(W_{k+1}(z)\) in the lower half plane, connected continuously to the left of the branch points. In other words, the branch cut on \((-\infty,0)\) is rotated counterclockwise to \((0,+\infty)\). (For \(k = -1\) and \(k = 0\), there is also a branch cut on \((-1/e,0)\), continuous from below instead of from above to maintain counterclockwise continuity.)
If middle is set, computes \(W_{\mathrm{middle}}(z)\) which corresponds to \(W_{-1}(z)\) in the upper half plane and \(W_{1}(z)\) in the lower half plane, connected continuously through \((-1/e,0)\) with branch cuts on \((-\infty,-1/e)\) and \((0,+\infty)\). \(W_{\mathrm{middle}}(z)\) extends the real analytic function \(W_{-1}(x)\) defined on \((-1/e,0)\) to a complex analytic function, whereas the standard branch \(W_{-1}(z)\) has a branch cut along the real segment.
>>> acb(-5,"+/- 1e-20").lambertw() [0.844844605432170 +/- 6.18e-16] + [+/- 1.98]j >>> acb(-5,"+/- 1e-20").lambertw(left=True) [0.844844605432170 +/- 7.85e-16] + [1.97500875488903 +/- 4.05e-15]j >>> acb(-5,"+/- 1e-20").lambertw(-1, left=True) [0.844844605432170 +/- 7.85e-16] + [-1.97500875488903 +/- 4.05e-15]j >>> acb(-0.25,"+/- 1e-20").lambertw(middle=True) [-2.15329236411035 +/- 1.48e-15] + [+/- 4.87e-16]j
-
legendre_p
(s, n, m=0, int type=2)¶ Legendre function of the first kind \(P_n^m(z)\).
>>> showgood(lambda: (acb(1)/3).legendre_p(5), dps=25) 0.3333333333333333333333333 >>> showgood(lambda: (acb(1)/3).legendre_p(5, 1.5), dps=25) -2.372124991643971726805456 >>> showgood(lambda: (acb(3)).legendre_p(5, 1.5, type=2), dps=25) -12091.31397811720689120900 - 12091.31397811720689120900j >>> showgood(lambda: (acb(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).
-
legendre_q
(s, n, m=0, int type=2)¶ Legendre function of the second kind \(Q_n^m(z)\).
>>> showgood(lambda: (acb(1)/3).legendre_q(5), dps=25) 0.1655245300933242182362054 >>> showgood(lambda: (acb(1)/3).legendre_q(5, 1.5), dps=25) -6.059967350218583975575616 >>> showgood(lambda: (acb(3)).legendre_q(5, 1.5, type=2), dps=25) -18992.99179585607569083751 + 18992.99179585607569083751j >>> showgood(lambda: (acb(3)).legendre_q(5, 1.5, type=3), dps=25) -0.0003010942389043591251820434j
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)\). The function is defined to be continuous away from the negative half-axis and thus differs from \(\log(\Gamma(s))\) in general.
>>> showgood(lambda: acb(1,2).lgamma(), dps=25) -1.876078786430929341229996 + 0.1296463163097883113837075j >>> showgood(lambda: (acb(0,10).lgamma() - acb(0,10).gamma().log()).imag / arb.pi(), dps=25) 4.000000000000000000000000
-
li
(s, bool offset=False)¶ Logarithmic integral \(\operatorname{li}(s)\), optionally the offset logarithmic integral \(\operatorname{Li}(s) = \operatorname{li}(s) - \operatorname{li}(2)\).
>>> showgood(lambda: acb(10).li(), dps=25) 6.165599504787297937522982 >>> showgood(lambda: acb(10).li(offset=True), dps=25) 5.120435724669805152678393
-
log
(s, bool analytic=False)¶ Natural logarithm \(\log(s)\). The analytic flag allows verifying that the branch cut is not touched; this is useful for numerical integration.
>>> showgood(lambda: acb(1,2).log(), dps=25) 0.8047189562170501873003797 + 1.107148717794090503017065j >>> showgood(lambda: acb(-5).log(), dps=25) 1.609437912434100374600759 + 3.141592653589793238462643j
-
log1p
(s)¶ Computes \(\log(1+s)\), accurately for small s.
>>> showgood(lambda: acb(1,2).log1p(), dps=25) 1.039720770839917964125848 + 0.7853981633974483096156608j >>> showgood(lambda: acb(0,"1e-100000000000000000").log1p(), dps=25) 5.000000000000000000000000e-200000000000000001 + 1.000000000000000000000000e-100000000000000000j
-
log_barnes_g
(s)¶ Logarithmic Barnes G-function \(\log G(s)\). Like the logarithmic gamma function, continuous analytic continuation is implied.
>>> showgood(lambda: acb(2+3j).log_barnes_g(), dps=25) -1.694395396880976849503750 - 3.389316783507118550918827j >>> showgood(lambda: acb(2+3j).barnes_g().log(), dps=25) -1.694395396880976849503750 + 2.893868523672467926006460j
-
log_sin_pi
(s)¶ Logarithmic sine function with analytic continuation defined to be consistent with the functional equation of the log gamma function.
>>> showgood(lambda: acb(5+2j).log_sin_pi(), dps=25) 5.590034639271204166020651 - 14.13716694115406957308190j >>> showgood(lambda: acb(5+2j).sin_pi().log(), dps=25) 5.590034639271204166020651 - 1.570796326794896619231322j
-
mid
(self)¶ Returns an exact acb representing the midpoint of self:
>>> acb("1 +/- 0.3", "2 +/- 0.4").mid() 1.00000000000000 + 2.00000000000000j
-
modular_delta
(tau)¶ Modular discriminant \(\Delta(\tau)\).
>>> showgood(lambda: acb(0.25+5j).modular_delta(), dps=25) 1.237896015010281684100435e-26 + 2.271101068324093838679275e-14j
-
modular_eta
(tau)¶ Dedekind eta function \(\eta(\tau)\).
>>> showgood(lambda: acb(1+1j).modular_eta(), dps=25) 0.7420487758365647263392722 + 0.1988313702299107190516142j
-
modular_j
(tau)¶ Modular j-invariant \(j(\tau)\).
>>> showgood(lambda: (1 + acb(-163).sqrt()/2).modular_j(), dps=25) 262537412640769488.0000000 >>> showgood(lambda: acb(3.25+100j).modular_j() - 744, dps=25, maxdps=2000) -3.817428033843319485125773e-539 - 7.503618895582604309279721e+272j
-
modular_lambda
(tau)¶ Modular lambda function \(\lambda(\tau)\).
>>> showgood(lambda: acb(0.25+5j).modular_lambda(), dps=25) 1.704995415668039343330405e-6 + 1.704992508662079437786598e-6j
-
modular_theta
(z, tau)¶ Computes the Jacobi theta functions \(\theta_1(z,\tau)\), \(\theta_2(z,\tau)\), \(\theta_3(z,\tau)\), \(\theta_4(z,\tau)\), returning all four values as a tuple. We define the theta functions with a factor \(\pi\) for z included; for example \(\theta_3(z,\tau) = 1 + 2 \sum_{n=1}^{\infty} q^{n^2} \cos(2n\pi z)\).
>>> for i in range(4): ... showgood(lambda: acb(1+1j).modular_theta(1.25+3j)[i], dps=25) ... 1.820235910124989594900076 - 1.216251950154477951760042j -1.220790267576967690128359 - 1.827055516791154669091679j 0.9694430387796704100046143 - 0.03055696120816803328582847j 1.030556961196006476576271 + 0.03055696120816803328582847j
-
overlaps
(self, other)¶
-
static
pi
()¶ Returns tthe constant \(\pi\) as an acb.
>>> showgood(lambda: acb.pi(), dps=25) 3.141592653589793238462643
-
polygamma
(self, s)¶ Polygamma function \(\psi_s(z)\) where z is given by self.
>>> showgood(lambda: acb(2+3j).polygamma(2), dps=25) 0.05267618908093586035755719 + 0.07303622933440580692454450j
-
polylog
(self, s)¶ Computes the polylogarithm \(\operatorname{Li}_s(z)\) where the argument z is given by self and the order s is given as an extra parameter.
>>> showgood(lambda: acb(3).polylog(2), dps=25) 2.320180423313098396406194 - 3.451392295223202661433821j >>> showgood(lambda: acb(2,3).polylog(1+2j), dps=25) -6.854984251829306907116775 + 7.375884252099214498443165j
-
pow
(s, t, bool analytic=False)¶ Power \(s^t\). The analytic flag allows verifying that the branch cut is not touched; this is useful for numerical integration.
>>> acb.integral(lambda z, a: z.pow(acb("1/3")), -5-1j, -5+1j) # WRONG!!! [+/- 5.03e-15] + [1.81137435753228 +/- 7.32e-15]j >>> acb.integral(lambda z, a: z.pow(acb("1/3"), analytic=a), -5-1j, -5+1j) [+/- 2.66e-14] + [1.8108516218463 +/- 3.58e-14]j
-
rad
(self)¶ Returns an upper bound for the radius (magnitude of error) of self as an arb.
>>> print(acb("1 +/- 0.3", "2 +/- 0.4").rad().str(5, radius=False)) 0.50000
-
real
¶
-
real_abs
(s, bool analytic=False)¶ Absolute value of a real variable, extended to a piecewise complex analytic function. This function is useful for integration.
>>> f = lambda z, a: (z**3).real_abs(analytic=a) >>> showgood(lambda: acb.integral(f, -1, 1), dps=25) 0.5000000000000000000000000
-
real_ceil
(s, bool analytic=False)¶ Ceiling function of a real variable, extended to a piecewise complex analytic function. This function is useful for integration.
>>> showgood(lambda: acb.integral(lambda x, a: x.real_ceil(analytic=a), 0, 100), dps=15) 5050.00000000000
-
real_floor
(s, bool analytic=False)¶ Floor function of a real variable, extended to a piecewise complex analytic function. This function is useful for integration.
>>> showgood(lambda: acb.integral(lambda x, a: x.real_floor(analytic=a), 0, 100), dps=15) 4950.00000000000
-
real_heaviside
(s, bool analytic=False)¶ Heaviside step function of a real variable, extended to a piecewise complex analytic function. This function is useful for integration.
>>> acb(-5+2j).real_heaviside(analytic=True) 0 >>> acb(5+2j).real_heaviside(analytic=True) 1.00000000000000 >>> acb(0).real_heaviside(analytic=True) nan + nanj >>> acb(0).real_heaviside() 0.500000000000000
-
real_max
(s, t, bool analytic=False)¶ Maximum value of two real variables, extended to a piecewise complex analytic function. This function is useful for integration.
>>> f = lambda x, a: x.sin().real_max(x.cos(), analytic=a) >>> showgood(lambda: acb.integral(f, 0, acb.pi())) 2.41421356237310
-
real_min
(s, t, bool analytic=False)¶ Minimum value of two real variables, extended to a piecewise complex analytic function. This function is useful for integration.
>>> f = lambda x, a: x.sin().real_min(x.cos(), analytic=a) >>> showgood(lambda: acb.integral(f, 0, acb.pi())) -0.414213562373095
-
real_sgn
(s, bool analytic=False)¶ Sign function of a real variable, extended to a piecewise complex analytic function. This function is useful for integration.
>>> acb(-5+2j).real_sgn() -1.00000000000000 >>> acb(0).real_sgn() 0 >>> acb(0).real_sgn(analytic=True) nan + nanj
-
real_sqrt
(s, bool analytic=False)¶ Square root of a real variable assumed to be nonnegative, extended to a piecewise complex analytic function. This function is useful for integration.
>>> acb.integral(lambda x, a: x.sqrt(analytic=a), 0, 1) [0.6666666666667 +/- 4.19e-14] + [+/- 1.12e-16]j >>> acb.integral(lambda x, a: x.real_sqrt(analytic=a), 0, 1) [0.6666666666667 +/- 4.19e-14]
-
rel_accuracy_bits
(self)¶
-
repr
(self)¶
-
rgamma
(s)¶ Reciprocal gamma function \(1/\Gamma(s)\), avoiding division by zero at the poles of the gamma function.
>>> showgood(lambda: acb(1,2).rgamma(), dps=25) 6.473073626019134501563613 - 0.8439438407732021454882999j >>> print(acb(0).rgamma()) 0 >>> print(acb(-1).rgamma()) 0
-
rising
(s, n)¶ Rising factorial \((s)_n\).
>>> showgood(lambda: acb(1,2).rising(5), dps=25) -540.0000000000000000000000 - 100.0000000000000000000000j >>> showgood(lambda: acb(1,2).rising(2+3j), dps=25) 0.3898076751098812033498554 - 0.01296289149274537721607465j
-
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.
>>> showgood(lambda: acb(1,2).rising2(5), dps=25) (-540.0000000000000000000000 - 100.0000000000000000000000j, -666.0000000000000000000000 + 420.0000000000000000000000j)
-
root
(s, ulong n)¶ Principal n-th root of s.
>>> showgood(lambda: acb(-1).root(3), dps=25) 0.5000000000000000000000000 + 0.8660254037844386467637232j >>> showgood(lambda: acb(10,11).root(3), dps=25) 2.364674532267765964410078 + 0.6739866114818132500777466j
-
rsqrt
(s, bool analytic=False)¶ Reciprocal square root \(1/\sqrt{s}\). The analytic flag allows verifying that the branch cut is not touched; this is useful for numerical integration.
>>> showgood(lambda: acb(1,2).rsqrt(), dps=25) 0.5688644810057831072783079 - 0.3515775842541429284870573j
-
sec
(s)¶ Secant function \(\sec(s)\).
>>> showgood(lambda: acb(2,3).sec(), dps=25) -0.04167496441114427004834991 + 0.09061113719623759652966120j
-
sech
(s)¶ Hyperbolic secant function \(\operatorname{sech}(s)\).
>>> showgood(lambda: acb(2,3).sech(), dps=25) -0.2635129751583893096436042 - 0.03621163655876852087145690j
-
sgn
(self)¶ Complex sign function.
>>> showgood(lambda: acb(-1).sgn(), dps=25) -1.000000000000000000000000 >>> showgood(lambda: acb(5,5).sgn(), dps=25) 0.7071067811865475244008444 + 0.7071067811865475244008444j >>> showgood(lambda: acb(0).sgn(), dps=25) 0
-
shi
(s)¶ Hyperbolic sine integral \(\operatorname{Shi}(s)\).
>>> showgood(lambda: acb(10).shi(), dps=25) 1246.114490199423344411882
-
si
(s)¶ Sine integral \(\operatorname{Si}(s)\).
>>> showgood(lambda: acb(10).si(), dps=25) 1.658347594218874049330972
-
sin
(s)¶ Sine function \(\sin(s)\).
>>> showgood(lambda: acb(1,2).sin(), dps=25) 3.165778513216168146740735 + 1.959601041421605897070352j
-
sin_cos
(s)¶ Computes \(\sin(s)\) and \(\cos(s)\) simultaneously.
>>> showgood(lambda: acb(1,2).sin_cos(), dps=15) (3.16577851321617 + 1.95960104142161j, 2.03272300701967 - 3.05189779915180j)
-
sin_cos_pi
(s)¶ Computes \(\sin(\pi s)\) and \(\cos(\pi s)\) simultaneously.
>>> showgood(lambda: acb(1,2).sin_cos_pi(), dps=25) (-267.7448940410165142571174j, -267.7467614837482222459319)
-
sin_pi
(s)¶ Sine function \(\sin(\pi s)\).
>>> showgood(lambda: acb(1,2).sin_pi(), dps=25) -267.7448940410165142571174j
-
sinc
(s)¶ Sinc function, \(\operatorname{sinc}(x) = \sin(x)/x\).
>>> showgood(lambda: acb(2,3).sinc(), dps=25) 0.4463290318402435457438585 - 2.753947027743647493993194j
-
sinc_pi
(s)¶ Normalized sinc function, \(\operatorname{sinc}(\pi x) = \sin(\pi x)/(\pi x)\).
>>> showgood(lambda: acb(2,3).sinc_pi(), dps=25) 455.1212281938610260024088 + 303.4141521292406840016059j
-
sinh
(s)¶ Hyperbolic sine function \(\sinh(s)\).
>>> showgood(lambda: acb(2,3).sinh(), dps=25) -3.590564589985779952012565 + 0.5309210862485198052670401j
-
sinh_cosh
(s)¶ Computes \(\sinh(s)\) and \(\cosh(s)\) simultaneously.
>>> showgood(lambda: acb(1,2).sinh_cosh(), dps=15) (-0.489056259041294 + 1.40311925062204j, -0.642148124715520 + 1.06860742138278j)
-
static
spherical_y
(long n, long m, theta, phi)¶ Spherical harmonic \(Y_n^m(\theta, \phi)\). The present implementation only supports integer n and m.
>>> showgood(lambda: acb.spherical_y(5, 3, 0.25, 0.75), dps=25) 0.02451377199072374024317003 - 0.03036343496553117039110087j
-
sqrt
(s, bool analytic=False)¶ Square root \(\sqrt{s}\).
>>> showgood(lambda: acb(1,2).sqrt(), dps=25) 1.272019649514068964252422 + 0.7861513777574232860695586j
The analytic flag allows verifying that the branch cut is not touched; this is useful for numerical integration.
>>> showgood(lambda: acb.integral(lambda z, a: z.sqrt(), 0, 1).real, dps=25) # WRONG!!! 0.6738873386790491615691993 >>> showgood(lambda: acb.integral(lambda z, a: z.sqrt(analytic=a), 0, 1).real, dps=25) 0.6666666666666666666666667
-
static
stieltjes
(n, a=1)¶ Generalized Stieltjes constant \(\gamma_n(a)\).
>>> showgood(lambda: acb.stieltjes(1), dps=25) -0.07281584548367672486058638 >>> showgood(lambda: acb.stieltjes(10**10), dps=25) 7.588362123713105194822403e+12397849705 >>> showgood(lambda: acb.stieltjes(1000, 2+3j), dps=25) -1.206122870741999199264747e+494 - 1.389205283963836265123849e+494j
-
str
(self, *args, **kwargs)¶
-
tan
(s)¶ Tangent function \(\tan(s)\).
>>> showgood(lambda: acb(1,2).tan(), dps=25) 0.03381282607989669028437056 + 1.014793616146633568117054j
-
tan_pi
(s)¶ Tangent function \(\tan(\pi s)\).
>>> showgood(lambda: acb(1,2).tan_pi(), dps=25) 0.9999930253396106106051072j
-
tanh
(s)¶ Hyperbolic tangent function \(\tanh(s)\).
>>> showgood(lambda: acb(2,3).tanh(), dps=25) 0.9653858790221331242784803 - 0.009884375038322493720314034j
-
unique_fmpz
(self)¶ If self represents exactly one integer, returns this value as an fmpz; otherwise returns None.
>>> acb("5 +/- 0.1").unique_fmpz() 5 >>> acb("5 +/- 0.9", "10 +/- 11").unique_fmpz() 5 >>> acb("5 +/- 0.9", 11).unique_fmpz() >>> >>> acb("5.1 +/- 0.9").unique_fmpz() >>>
-
zeta
(s, a=None)¶ Riemann zeta function \(\zeta(s)\), or the Hurwitz zeta function \(\zeta(s,a)\) if a second parameter is passed.
>>> showgood(lambda: acb(0.5,1000).zeta(), dps=25) 0.3563343671943960550744025 + 0.9319978312329936651150604j >>> showgood(lambda: acb(1,2).zeta(acb(2,3)), dps=25) -2.953059572088556722876240 + 3.410962524512050603254574j
-