Python interface (python-flint)¶
Introduction¶
Examples:
>>> from flint2 import *
>>> QQ.bernoulli(50)
495057205241079648212477525/66
>>> _.factor()
(1, [5, 417202699, 47464429777438199, 2, 3, 11], [2, 1, 1, -1, -1, -1])
>>> sign * (primes ** exponents).product()
495057205241079648212477525/66
Types, parents and coercions¶
>>> ZZ(5)
5
>>> _.parent()
Integer ring (fmpz)
>>> QQ(5)
5
>>> _.parent()
Rational field (fmpq)
>>> ZZ(10) / ZZ(6)
Traceback (most recent call last):
...
FlintDomainError: x / y is not an element of {Integer ring (fmpz)} for {x = 10}, {y = 6}
>>> x = QQ(1) / 2; x ** x
Traceback (most recent call last):
...
FlintDomainError: x ** y is not an element of {Rational field (fmpq)} for {x = 1/2}, {y = 1/2}
>>> ZZ(10) / QQ(6)
5/3
>>> x = QQbar(1) / 2; x ** x
Root a = 0.707107 of 2*a^2-1
Real and complex numbers¶
>>> RR.zeta(2)
[1.644934066848226 +/- 4.57e-16]
>>> RR.prec = 128
>>> RR.zeta(2)
[1.64493406684822643647241516664602518922 +/- 2.88e-39]
>>> RR.prec = 53 # restore default
API documentation¶
- class flint_ctypes.FlintDomainError¶
Raised when an operation does not have a well-defined result in the target domain.
- class flint_ctypes.FlintUnableError¶
Raised when an operation cannot be performed because the algorithm is not implemented or there is insufficient precision, memory, etc.
- class flint_ctypes.fmpz_poly_struct[source]¶
- alloc¶
Structure/Union member
- coeffs¶
Structure/Union member
- length¶
Structure/Union member
- class flint_ctypes.fmpq_poly_struct[source]¶
- alloc¶
Structure/Union member
- coeffs¶
Structure/Union member
- den¶
Structure/Union member
- length¶
Structure/Union member
- class flint_ctypes.qqbar_struct[source]¶
- enclosure¶
Structure/Union member
- poly¶
Structure/Union member
- class flint_ctypes.nmod_poly_struct[source]¶
- alloc¶
Structure/Union member
- coeffs¶
Structure/Union member
- length¶
Structure/Union member
- n¶
Structure/Union member
- ninv¶
Structure/Union member
- nnorm¶
Structure/Union member
- class flint_ctypes.fq_struct[source]¶
- alloc¶
Structure/Union member
- coeffs¶
Structure/Union member
- length¶
Structure/Union member
- class flint_ctypes.fq_nmod_struct[source]¶
- alloc¶
Structure/Union member
- coeffs¶
Structure/Union member
- length¶
Structure/Union member
- n¶
Structure/Union member
- ninv¶
Structure/Union member
- nnorm¶
Structure/Union member
- class flint_ctypes.gr_vec_struct[source]¶
- alloc¶
Structure/Union member
- entries¶
Structure/Union member
- length¶
Structure/Union member
- class flint_ctypes.gr_poly_struct[source]¶
- alloc¶
Structure/Union member
- coeffs¶
Structure/Union member
- length¶
Structure/Union member
- class flint_ctypes.gr_series_struct[source]¶
- alloc¶
Structure/Union member
- coeffs¶
Structure/Union member
- error¶
Structure/Union member
- length¶
Structure/Union member
- class flint_ctypes.psl2z_struct[source]¶
- a¶
Structure/Union member
- b¶
Structure/Union member
- c¶
Structure/Union member
- d¶
Structure/Union member
- class flint_ctypes.dirichlet_char_struct[source]¶
- log¶
Structure/Union member
- n¶
Structure/Union member
- class flint_ctypes.gr_mat_struct[source]¶
- c¶
Structure/Union member
- entries¶
Structure/Union member
- r¶
Structure/Union member
- rows¶
Structure/Union member
- class flint_ctypes.LogicContext(value)[source]¶
Handle the result of predicates (experimental):
>>> a = (RR_arb(1) / 3) * 3 >>> a [1.00000000000000 +/- 3.89e-16] >>> with strict_logic: ... a == 1 ... Traceback (most recent call last): ... Undecidable: unable to decide x == y for x = [1.00000000000000 +/- 3.89e-16], y = 1.000000000000000 over Real numbers (arb, prec = 53) >>> with pessimistic_logic: ... a == 1 ... False >>> with optimistic_logic: ... a == 1 ... True
- class flint_ctypes.gr_ctx[source]¶
-
- property prec¶
- gen()[source]¶
Generator of this domain.
>>> QQbar.i() Root a = 1.00000*I of a^2+1 >>> QQ.i() Traceback (most recent call last): ... FlintDomainError: i is not an element of {Rational field (fmpq)}
- i()[source]¶
Imaginary unit as an element of this domain.
>>> QQbar.i() Root a = 1.00000*I of a^2+1 >>> QQ.i() Traceback (most recent call last): ... FlintDomainError: i is not an element of {Rational field (fmpq)}
- pi()[source]¶
The number pi as an element of this domain.
>>> RR.pi() [3.141592653589793 +/- 3.39e-16] >>> QQbar.pi() Traceback (most recent call last): ... FlintDomainError: pi is not an element of {Complex algebraic numbers (qqbar)}
- euler()[source]¶
Euler’s constant as an element of this domain.
>>> RR.euler() [0.5772156649015329 +/- 9.00e-17]
We do not know whether Euler’s constant is rational:
>>> QQ.euler() Traceback (most recent call last): ... FlintUnableError: failed to compute euler in {Rational field (fmpq)}
- catalan()[source]¶
Catalan’s constant as an element of this domain.
>>> RR.catalan() [0.915965594177219 +/- 1.23e-16]
- khinchin()[source]¶
Khinchin’s constant as an element of this domain.
>>> RR.khinchin() [2.685452001065306 +/- 6.82e-16]
- glaisher()[source]¶
Khinchin’s constant as an element of this domain.
>>> RR.glaisher() [1.282427129100623 +/- 6.02e-16]
- fresnel_s(x, normalized=False)[source]¶
>>> RR.fresnel_s(1) [0.3102683017233811 +/- 2.67e-18] >>> RR.fresnel_s(1, normalized=True) [0.4382591473903548 +/- 9.24e-17]
- fresnel_c(x, normalized=False)[source]¶
>>> RR.fresnel_c(1) [0.904524237900272 +/- 1.46e-16] >>> RR.fresnel_c(1, normalized=True) [0.779893400376823 +/- 3.59e-16]
- fresnel(x, normalized=False)[source]¶
>>> RR.fresnel(1) ([0.3102683017233811 +/- 2.67e-18], [0.904524237900272 +/- 1.46e-16]) >>> RR.fresnel(1, normalized=True) ([0.4382591473903548 +/- 9.24e-17], [0.779893400376823 +/- 3.59e-16])
- gamma_upper(x, y, regularized=0)[source]¶
>>> RR.gamma_upper(3, 4) [0.476206611107089 +/- 5.30e-16] >>> RR.gamma_upper(3, 4, regularized=True) [0.2381033055535443 +/- 8.24e-17]
- gamma_lower(x, y, regularized=0)[source]¶
>>> RR.gamma_lower(3, 4) [1.52379338889291 +/- 2.89e-15] >>> RR.gamma_lower(3, 4, regularized=True) [0.76189669444646 +/- 6.52e-15]
- beta_lower(a, b, x, regularized=0)[source]¶
>>> RR.beta_lower(2, 3, 0.5) [0.0572916666666667 +/- 6.08e-17] >>> RR.beta_lower(2, 3, 0.5, regularized=True) [0.687500000000000 +/- 5.24e-16]
- log_integral(x, offset=False)[source]¶
>>> RR.log_integral(2) [1.04516378011749 +/- 4.01e-15] >>> RR.log_integral(2, offset=True) 0
- bessel_i(x, y, scaled=False)[source]¶
>>> RR.bessel_i(2, 3) [2.24521244092995 +/- 1.88e-15] >>> RR.bessel_i(2, 3, scaled=True) [0.111782545296958 +/- 2.09e-16]
- bessel_k(x, y, scaled=False)[source]¶
>>> RR.bessel_k(2, 3) [0.06151045847174 +/- 8.87e-15] >>> RR.bessel_k(2, 3, scaled=True) [1.235470584796 +/- 5.14e-13]
- coulomb_hpos(x, y, z)[source]¶
>>> CC.coulomb_hpos(2, 3, 4) ([5.371722466 +/- 6.15e-10] + [0.101631502833431 +/- 8.03e-16]*I)
- coulomb_hneg(x, y, z)[source]¶
>>> CC.coulomb_hneg(2, 3, 4) ([5.371722466 +/- 6.15e-10] + [-0.101631502833431 +/- 8.03e-16]*I)
- chebyshev_t(n, x)[source]¶
Chebyshev polynomial of the first kind.
>>> [ZZ.chebyshev_t(n, 2) for n in range(5)] [1, 2, 7, 26, 97] >>> RR.chebyshev_t(0.5, 0.75) [0.935414346693485 +/- 5.18e-16] >>> ZZx.chebyshev_t(4, [0, 1]) 1 - 8*x^2 + 8*x^4
- chebyshev_u(n, x)[source]¶
Chebyshev polynomial of the second kind.
>>> [ZZ.chebyshev_u(n, 2) for n in range(5)] [1, 4, 15, 56, 209] >>> RR.chebyshev_u(0.5, 0.75) [1.33630620956212 +/- 2.68e-15] >>> ZZx.chebyshev_u(4, [0, 1]) 1 - 12*x^2 + 16*x^4
- jacobi_p(n, a, b, x)[source]¶
Jacobi polynomial.
>>> RR.jacobi_p(3, 1, 2, 4) [602.500000000000 +/- 3.28e-13]
- gegenbauer_c(n, m, x)[source]¶
Gegenbauer polynomial.
>>> RR.gegenbauer_c(3, 2, 4) [2000.00000000000 +/- 3.60e-12]
- laguerre_l(n, m, x)[source]¶
Associated Laguerre polynomial (or Laguerre function).
>>> RR.laguerre_l(3, 2, 4) [-0.66666666666667 +/- 5.71e-15]
- hermite_h(n, x)[source]¶
Hermite polynomial (Hermite function).
>>> RR.hermite_h(3, 4) 464.0000000000000
- spherical_y(n, m, theta, phi)[source]¶
Spherical harmonic.
>>> CC.spherical_y(4, 3, 0.5, 0.75) ([0.076036396941350 +/- 2.18e-16] + [-0.094180781089734 +/- 4.96e-16]*I)
- legendre_p_root(n, k, weight=False)[source]¶
Root of Legendre polynomial. With weight=True, also returns the corresponding weight for Gauss-Legendre quadrature.
>>> RR.legendre_p_root(5, 1) [0.538469310105683 +/- 1.15e-16] >>> RR.legendre_p(5, 0, RR.legendre_p_root(5, 1)) [+/- 8.15e-16] >>> RR.legendre_p_root(5, 1, weight=True) ([0.538469310105683 +/- 1.15e-16], [0.4786286704993664 +/- 7.10e-17])
- hypgeom_0f1(a, z, regularized=False)[source]¶
Hypergeometric function 0F1, optionally regularized.
>>> RR.hypgeom_0f1(3, 4) [3.21109468764205 +/- 5.00e-15] >>> RR.hypgeom_0f1(3, 4, regularized=True) [1.60554734382103 +/- 5.20e-15] >>> CC.hypgeom_0f1(1, 2+2j) ([2.435598449671389 +/- 7.27e-16] + [4.43452765355337 +/- 4.91e-15]*I)
- hypgeom_1f1(a, b, z, regularized=False)[source]¶
Hypergeometric function 1F1, optionally regularized.
>>> RR.hypgeom_1f1(3, 4, 5) [60.504568913851 +/- 3.82e-13] >>> RR.hypgeom_1f1(3, 4, 5, regularized=True) [10.0840948189752 +/- 3.31e-14]
- hypgeom_u(a, b, z)[source]¶
Hypergeometric function U.
>>> RR.hypgeom_u(1, 2, 3) [0.3333333333333333 +/- 7.04e-17]
- hypgeom_2f1(a, b, c, z, regularized=False)[source]¶
Hypergeometric function 2F1, optionally regularized.
>>> RR.hypgeom_2f1(1, 2, 3, -4) [0.29882026094574 +/- 8.48e-15] >>> RR.hypgeom_2f1(1, 2, 3, -4, regularized=True) [0.14941013047287 +/- 4.24e-15]
- hypgeom_pfq(a, b, z, regularized=False)[source]¶
Generalized hypergeometric function, optionally regularized.
>>> RR.hypgeom_pfq([1,2], [3,4], 0.5) [1.09002619782383 +/- 4.32e-15] >>> RR.hypgeom_pfq([1, 2], [3, 4], 0.5, regularized=True) [0.090835516485319 +/- 2.36e-16] >>> CC.hypgeom_pfq([1,2], [3,4], 0.5+0.5j) ([1.08239550393928 +/- 2.16e-15] + [0.096660812453003 +/- 5.55e-16]*I)
- fac(x)[source]¶
Factorial.
>>> ZZ.fac(10) 3628800 >>> ZZ.fac(-1) Traceback (most recent call last): ... FlintDomainError: fac(x) is not an element of {Integer ring (fmpz)} for {x = -1}
Real and complex factorials extend using the gamma function:
>>> RR.fac(10**20) [1.93284951431010e+1956570551809674817245 +/- 3.03e+1956570551809674817230] >>> RR.fac(0.5) [0.886226925452758 +/- 1.78e-16] >>> CC.fac(1+1j) ([0.652965496420167 +/- 6.21e-16] + [0.343065839816545 +/- 5.38e-16]*I)
Factorials mod N:
>>> ZZmod(10**7 + 19).fac(10**7) 2343096
- fac_vec(length)[source]¶
Vector of factorials.
>>> ZZ.fac_vec(10) [1, 1, 2, 6, 24, 120, 720, 5040, 40320, 362880] >>> QQ.fac_vec(10) / 3 [1/3, 1/3, 2/3, 2, 8, 40, 240, 1680, 13440, 120960] >>> ZZmod(7).fac_vec(10) [1, 1, 2, 6, 3, 1, 6, 0, 0, 0] >>> sum(RR.fac_vec(100)) [9.427862397658e+155 +/- 3.19e+142]
- rfac(x)[source]¶
Reciprocal factorial.
>>> QQ.rfac(5) 1/120 >>> ZZ.rfac(-2) 0 >>> ZZ.rfac(2) Traceback (most recent call last): ... FlintDomainError: rfac(x) is not an element of {Integer ring (fmpz)} for {x = 2} >>> RR.rfac(0.5) [1.128379167095513 +/- 7.02e-16]
- rfac_vec(length)[source]¶
Vector of reciprocal factorials.
>>> QQ.rfac_vec(8) [1, 1, 1/2, 1/6, 1/24, 1/120, 1/720, 1/5040] >>> ZZmod(7).rfac_vec(7) [1, 1, 4, 6, 5, 1, 6] >>> ZZmod(7).rfac_vec(8) Traceback (most recent call last): ... FlintDomainError: rfac_vec(length) is not an element of {Integers mod 7 (_gr_nmod)} for {length = 8} >>> sum(RR.rfac_vec(20)) [2.71828182845904 +/- 8.66e-15]
- rising(x, n)[source]¶
Rising factorial.
>>> [ZZ.rising(3, k) for k in range(5)] [1, 3, 12, 60, 360] >>> ZZx.rising(ZZx([0,1]), 5) 24*x + 50*x^2 + 35*x^3 + 10*x^4 + x^5 >>> RR.rising(1, 10**7) [1.202423400515903e+65657059 +/- 5.57e+65657043]
- falling(x, n)[source]¶
Falling factorial.
>>> [ZZ.falling(3, k) for k in range(5)] [1, 3, 6, 6, 0] >>> ZZx.falling(ZZx([0,1]), 5) 24*x - 50*x^2 + 35*x^3 - 10*x^4 + x^5 >>> RR.log(RR.falling(RR.pi(), 10**7)) [151180898.7174084 +/- 9.72e-8]
- bin(x, y)[source]¶
Binomial coefficient.
>>> [ZZ.bin(5, k) for k in range(7)] [1, 5, 10, 10, 5, 1, 0] >>> RR.bin(100000, 50000) [2.52060836892200e+30100 +/- 5.36e+30085] >>> ZZmod(1000).bin(10000, 3000) 200 >>> ZZp64.bin(100000, 50000) 5763493550349629692 >>> ZZp64.bin(10**30, 2) 998763921924463582
- bin_vec(n, length=None)[source]¶
Vector of binomial coefficients, optionally truncated to specified length.
>>> ZZ.bin_vec(8) [1, 8, 28, 56, 70, 56, 28, 8, 1] >>> ZZmod(5).bin_vec(8) [1, 3, 3, 1, 0, 1, 3, 3, 1] >>> ZZ.bin_vec(0) [1] >>> ZZ.bin_vec(1000, 3) [1, 1000, 499500] >>> ZZ.bin_vec(4, 8) [1, 4, 6, 4, 1, 0, 0, 0] >>> QQ.bin_vec(QQ(1)/2, 5) [1, 1/2, -1/8, 1/16, -5/128]
- doublefac(x)[source]¶
Double factorial (semifactorial).
>>> [ZZ.doublefac(n) for n in range(10)] [1, 1, 2, 3, 8, 15, 48, 105, 384, 945] >>> RR.doublefac(2.5) [2.40706945611604 +/- 5.54e-15] >>> CC.doublefac(1+1j) ([0.250650779545753 +/- 7.56e-16] + [0.100474421235437 +/- 4.14e-16]*I)
- harmonic(x)[source]¶
Harmonic numbers.
>>> [QQ.harmonic(n) for n in range(6)] [0, 1, 3/2, 11/6, 25/12, 137/60] >>> RR.harmonic(10**9) [21.30048150234794 +/- 8.48e-15] >>> ZZp64.harmonic(1000) 6514760847963681162
- beta(x, y)[source]¶
Beta function.
>>> RR.beta(3, 4.5) [0.01243201243201243 +/- 6.93e-18] >>> CC.beta(1j, 1+1j) ([-1.18807306241087 +/- 5.32e-15] + [-1.31978426013907 +/- 4.09e-15]*I)
- barnes_g(x)[source]¶
Barnes G-function.
>>> RR.barnes_g(7) 34560.00000000000 >>> CC.barnes_g(1+2j) ([0.54596949228965 +/- 7.69e-15] + [-3.98421873125106 +/- 8.76e-15]*I)
- log_barnes_g(x)[source]¶
Logarithmic Barnes G-function.
>>> RR.log_barnes_g(100) [15258.0613921488 +/- 3.87e-11] >>> CC.log_barnes_g(10+20j) ([-452.057343313397 +/- 6.85e-13] + [121.014356688943 +/- 2.52e-13]*I)
- stieltjes(n, a=1)[source]¶
Stieltjes constant.
>>> CC.stieltjes(1) [-0.0728158454836767 +/- 2.78e-17] >>> CC.stieltjes(1, a=0.5) [-1.353459680804942 +/- 7.22e-16]
- polygamma(s, z)[source]¶
Polygamma function.
>>> CC.polygamma(2, 3) [-0.1541138063191886 +/- 7.16e-17]
- lerch_phi(z, s, a)[source]¶
>>> CC.lerch_phi(2, 3, 4) ([-0.00213902437921 +/- 1.70e-15] + [-0.04716836434127 +/- 5.28e-15]*I)
- dirichlet_eta(x)[source]¶
Dirichlet eta function.
>>> CC.dirichlet_eta(1) [0.6931471805599453 +/- 6.93e-17] >>> CC.dirichlet_eta(2) [0.822467033424113 +/- 2.36e-16]
- riemann_xi(x)[source]¶
Riemann xi function.
>>> s = 2+3j; CC.riemann_xi(s); CC.riemann_xi(1-s) ([0.41627125989962 +/- 4.65e-15] + [0.08882330496564 +/- 1.43e-15]*I) ([0.41627125989962 +/- 4.65e-15] + [0.08882330496564 +/- 1.43e-15]*I)
- bernoulli(n)[source]¶
Bernoulli number \(B_n\) as an element of this domain.
>>> QQ.bernoulli(10) 5/66 >>> RR.bernoulli(10) [0.0757575757575757 +/- 5.97e-17]
>>> ZZ.bernoulli(0) 1 >>> ZZ.bernoulli(1) Traceback (most recent call last): ... FlintDomainError: bernoulli(n) is not an element of {Integer ring (fmpz)} for {n = 1}
Huge Bernoulli numbers can be computed numerically:
>>> RR.bernoulli(10**20) [-1.220421181609039e+1876752564973863312289 +/- 4.69e+1876752564973863312273] >>> RF.bernoulli(10**20) -1.220421181609039e+1876752564973863312289 >>> QQ.bernoulli(10**20) Traceback (most recent call last): ... FlintUnableError: failed to compute bernoulli(n) in {Rational field (fmpq)} for {n = 100000000000000000000}
- bernoulli_vec(length)[source]¶
Vector of Bernoulli numbers.
>>> QQ.bernoulli_vec(12) [1, -1/2, 1/6, 0, -1/30, 0, 1/42, 0, -1/30, 0, 5/66, 0] >>> CC_ca.bernoulli_vec(5) [1, -0.500000 {-1/2}, 0.166667 {1/6}, 0, -0.0333333 {-1/30}] >>> sum(RR.bernoulli_vec(100)) [1.127124216595034e+76 +/- 6.74e+60] >>> sum(RF.bernoulli_vec(100)) 1.127124216595034e+76 >>> sum(CC.bernoulli_vec(100)) [1.127124216595034e+76 +/- 6.74e+60]
- eulernum(n)[source]¶
Euler number \(E_n\) as an element of this domain.
>>> ZZ.eulernum(10) -50521 >>> RR.eulernum(10) -50521.00000000000
Huge Euler numbers can be computed numerically:
>>> RR.eulernum(10**20) [4.346791453661149e+1936958564106659551331 +/- 8.35e+1936958564106659551315] >>> RF.eulernum(10**20) 4.346791453661149e+1936958564106659551331 >>> ZZ.eulernum(10**20) Traceback (most recent call last): ... FlintUnableError: failed to compute eulernum(n) in {Integer ring (fmpz)} for {n = 100000000000000000000}
- eulernum_vec(length)[source]¶
Vector of Euler numbers.
>>> ZZ.eulernum_vec(12) [1, 0, -1, 0, 5, 0, -61, 0, 1385, 0, -50521, 0] >>> QQ.eulernum_vec(12) / 3 [1/3, 0, -1/3, 0, 5/3, 0, -61/3, 0, 1385/3, 0, -50521/3, 0] >>> sum(RR.eulernum_vec(100)) [-7.23465655613392e+134 +/- 3.20e+119] >>> sum(RF.eulernum_vec(100)) -7.234656556133921e+134
- fib(n)[source]¶
Fibonacci number \(F_n\) as an element of this domain.
>>> ZZ.fib(10) 55 >>> RR.fib(10) 55.00000000000000 >>> ZZ.fib(-10) -55
Huge Fibonacci numbers can be computed numerically and in modular arithmetic:
>>> RR.fib(10**20) [3.78202087472056e+20898764024997873376 +/- 4.02e+20898764024997873361] >>> RF.fib(10**20) 3.782020874720557e+20898764024997873376 >>> F = FiniteField_fq(17, 1) >>> n = 10**20; F.fib(n); F.fib(n-1) + F.fib(n-2) 13 13
- fib_vec(length)[source]¶
Vector of Fibonacci numbers.
>>> ZZ.fib_vec(10) [0, 1, 1, 2, 3, 5, 8, 13, 21, 34] >>> QQ.fib_vec(10) / 3 [0, 1/3, 1/3, 2/3, 1, 5/3, 8/3, 13/3, 7, 34/3] >>> sum(RR.fib_vec(100)) [5.7314784401...e+20 +/- ...] >>> sum(RF.fib_vec(100)) 5.731478440138172e+20
- stirling_s1u(n, k)[source]¶
Unsigned Stirling number of the first kind.
>>> ZZ.stirling_s1u(5, 2) 50 >>> QQ.stirling_s1u(5, 2) 50 >>> ZZ.stirling_s1u(50, 21) 33187391298039120738041153829116024033357291261862000 >>> RR.stirling_s1u(50, 21) [3.318739129803912e+52 +/- 8.66e+36]
- stirling_s1(n, k)[source]¶
Signed Stirling number of the first kind.
>>> ZZ.stirling_s1(5, 2) -50 >>> QQ.stirling_s1(5, 2) -50 >>> RR.stirling_s1(5, 2) -50.00000000000000
- stirling_s2(n, k)[source]¶
Stirling number of the second kind.
>>> ZZ.stirling_s2(5, 2) 15 >>> QQ.stirling_s2(5, 2) 15 >>> RR.stirling_s2(5, 2) 15.00000000000000 >>> RR.stirling_s2(50, 20) [7.59792160686099e+45 +/- 5.27e+30]
- stirling_s1u_vec(n, length=None)[source]¶
Vector of unsigned Stirling numbers of the first kind, optionally truncated to specified length.
>>> ZZ.stirling_s1u_vec(5) [0, 24, 50, 35, 10, 1] >>> QQ.stirling_s1u_vec(5) / 3 [0, 8, 50/3, 35/3, 10/3, 1/3] >>> RR.stirling_s1u_vec(5, 3) [0, 24.00000000000000, 50.00000000000000]
- stirling_s1_vec(n, length=None)[source]¶
Vector of signed Stirling numbers of the first kind, optionally truncated to specified length.
>>> ZZ.stirling_s1_vec(5) [0, 24, -50, 35, -10, 1] >>> QQ.stirling_s1_vec(5) / 3 [0, 8, -50/3, 35/3, -10/3, 1/3] >>> RR.stirling_s1_vec(5, 3) [0, 24.00000000000000, -50.00000000000000]
- stirling_s2_vec(n, length=None)[source]¶
Vector of Stirling numbers of the second kind, optionally truncated to specified length.
>>> ZZ.stirling_s2_vec(5) [0, 1, 15, 25, 10, 1] >>> QQ.stirling_s2_vec(5) / 3 [0, 1/3, 5, 25/3, 10/3, 1/3] >>> RR.stirling_s2_vec(5, 3) [0, 1.000000000000000, 15.00000000000000]
- bellnum(n)[source]¶
Bell number \(E_n\) as an element of this domain.
>>> ZZ.bellnum(10) 115975 >>> RR.bellnum(10) 115975.0000000000 >>> ZZp64.bellnum(10000) 355901145009109239 >>> ZZmod(1000).bellnum(10000) 635
Huge Bell numbers can be computed numerically:
>>> RR.bellnum(10**20) [5.38270113176282e+1794956117137290721328 +/- 5.44e+1794956117137290721313] >>> ZZ.bellnum(10**20) Traceback (most recent call last): ... FlintUnableError: failed to compute bellnum(n) in {Integer ring (fmpz)} for {n = 100000000000000000000}
- bellnum_vec(length)[source]¶
Vector of Bell numbers.
>>> ZZ.bellnum_vec(10) [1, 1, 2, 5, 15, 52, 203, 877, 4140, 21147] >>> QQ.bellnum_vec(10) / 3 [1/3, 1/3, 2/3, 5/3, 5, 52/3, 203/3, 877/3, 1380, 7049] >>> RR.bellnum_vec(100).sum() [1.67618752079292e+114 +/- 4.30e+99] >>> RF.bellnum_vec(100).sum() 1.676187520792924e+114 >>> ZZmod(10000).bellnum_vec(10000).sum() 7337
- partitions(n)[source]¶
Partition function \(p(n)\) as an element of this domain.
>>> ZZ.partitions(10) 42 >>> QQ.partitions(10) / 5 42/5 >>> RR.partitions(10) 42.00000000000000 >>> RR.partitions(10**20) [1.838176508344883e+11140086259 +/- 8.18e+11140086243]
- partitions_vec(length)[source]¶
Vector of partition numbers.
>>> ZZ.partitions_vec(10) [1, 1, 2, 3, 5, 7, 11, 15, 22, 30] >>> QQ.partitions_vec(10) / 3 [1/3, 1/3, 2/3, 1, 5/3, 7/3, 11/3, 5, 22/3, 10] >>> ZZmod(10).partitions_vec(10) [1, 1, 2, 3, 5, 7, 1, 5, 2, 0] >>> sum(ZZmod(10).partitions_vec(100)) 6 >>> sum(RR.partitions_vec(100)) 1452423276.000000
- zeta_zero(n)[source]¶
Zero of the Riemann zeta function.
>>> CC.zeta_zero(1) (0.5000000000000000 + [14.13472514173469 +/- 4.71e-15]*I) >>> CC.zeta_zero(2) (0.5000000000000000 + [21.02203963877155 +/- 6.02e-15]*I)
- zeta_zeros(num, start=1)[source]¶
Zeros of the Riemann zeta function.
>>> [x.im() for x in CC.zeta_zeros(4)] [[14.13472514173469 +/- 4.71e-15], [21.02203963877155 +/- 6.02e-15], [25.01085758014569 +/- 7.84e-15], [30.42487612585951 +/- 5.96e-15]] >>> [x.im() for x in CC.zeta_zeros(2, start=100)] [[236.5242296658162 +/- 3.51e-14], [237.7698204809252 +/- 5.29e-14]]
- zeta_nzeros(t)[source]¶
Number of zeros of Riemann zeta function up to given height.
>>> CC.zeta_nzeros(100) 29.00000000000000
- dirichlet_l(s, chi)[source]¶
Dirichlet L-function with character chi.
>>> CC.dirichlet_l(2, DirichletGroup(1)(1)) [1.644934066848226 +/- 4.57e-16] >>> RR.dirichlet_l(2, DirichletGroup(4)(3)) [0.915965594177219 +/- 2.68e-16] >>> CC.dirichlet_l(2+3j, DirichletGroup(7)(3)) ([1.273313649440491 +/- 9.69e-16] + [-0.074323294425594 +/- 6.96e-16]*I)
- hardy_theta(s, chi=None)[source]¶
Hardy theta function.
>>> CC.hardy_theta(10) [-3.06707439628989 +/- 6.66e-15] >>> CC.hardy_theta(10, DirichletGroup(4)(3)) [4.64979557270698 +/- 4.41e-15]
- hardy_z(s, chi=None)[source]¶
Hardy Z-function.
>>> CC.hardy_z(2) [-0.539633125646145 +/- 8.59e-16] >>> CC.hardy_z(2, DirichletGroup(4)(3)) [1.15107760668266 +/- 5.01e-15]
- dirichlet_chi(n, chi)[source]¶
Value of the Dirichlet character chi(n).
>>> chi = DirichletGroup(5)(3) >>> [CC.dirichlet_chi(n, chi) for n in range(5)] [0, 1.000000000000000, -1.000000000000000*I, 1.000000000000000*I, -1.000000000000000]
- modular_lambda(tau)[source]¶
Modular lambda function lambda(tau).
>>> CC.modular_lambda(1j) [0.50000000000000 +/- 2.16e-15]
- modular_delta(tau)[source]¶
Modular discriminant delta(tau).
>>> CC.modular_delta(1j) [0.0017853698506421 +/- 6.01e-17]
- dedekind_eta(tau)[source]¶
Dedekind eta function eta(tau).
>>> CC.dedekind_eta(1j) [0.768225422326057 +/- 9.03e-16]
- hilbert_class_poly(D, x)[source]¶
Hilbert class polynomial H_D(x) evaluated at x.
>>> ZZx.hilbert_class_poly(-20, ZZx.gen()) -681472000 - 1264000*x + x^2 >>> CC.hilbert_class_poly(-20, 1+1j) (-682736000.0000000 - 1263998.000000000*I) >>> ZZx.hilbert_class_poly(-21, ZZx.gen()) Traceback (most recent call last): ... FlintDomainError: hilbert_class_poly(D, x) is not an element of {Ring of polynomials over Integer ring (fmpz)} for {D = -21}, {x = x}
- eisenstein_g(n, tau)[source]¶
Eisenstein series G_n(tau).
>>> CC.eisenstein_g(2, 1j) [3.14159265358979 +/- 8.71e-15] >>> CC.eisenstein_g(4, 1j); RR.gamma(0.25)**8 / (960 * RR.pi()**2) [3.1512120021539 +/- 3.41e-14] [3.15121200215390 +/- 7.72e-15]
- eisenstein_e(n, tau)[source]¶
Eisenstein series E_n(tau).
>>> CC.eisenstein_e(2, 1j) [0.95492965855137 +/- 3.85e-15] >>> CC.eisenstein_e(4, 1j); 3*RR.gamma(0.25)**8/(64*RR.pi()**6) [1.4557628922687 +/- 1.32e-14] [1.45576289226871 +/- 3.76e-15]
- eisenstein_g_vec(tau, n)[source]¶
Vector of Eisenstein series [G_4(tau), G_6(tau), …]. Note that G_2(tau) is omitted.
>>> CC.eisenstein_g_vec(1j, 3) [[3.1512120021539 +/- 3.41e-14], [+/- 4.40e-14], [4.255773035365 +/- 2.12e-13]]
- agm(x, y=None)[source]¶
Arithmetic-geometric mean.
>>> RR.agm(2) [1.45679103104691 +/- 3.98e-15] >>> RR.agm(2, 3) [2.47468043623630 +/- 4.68e-15]
- jacobi_theta(z, tau)[source]¶
Simultaneous computation of the four Jacobi theta functions.
>>> CC.jacobi_theta(0.125, 1j) ([0.347386687929454 +/- 3.21e-16], [0.843115469091413 +/- 8.18e-16], [1.061113709291166 +/- 5.74e-16], [0.938886290708834 +/- 3.52e-16])
- jacobi_theta_1(z, tau)[source]¶
Jacobi theta function.
>>> CC.jacobi_theta_1(0.125, 1j) [0.347386687929454 +/- 3.21e-16]
- jacobi_theta_2(z, tau)[source]¶
Jacobi theta function.
>>> CC.jacobi_theta_2(0.125, 1j) [0.843115469091413 +/- 8.18e-16]
- jacobi_theta_3(z, tau)[source]¶
Jacobi theta function.
>>> CC.jacobi_theta_3(0.125, 1j) [1.061113709291166 +/- 5.74e-16]
- jacobi_theta_4(z, tau)[source]¶
Jacobi theta function.
>>> CC.jacobi_theta_4(0.125, 1j) [0.938886290708834 +/- 3.52e-16]
- elliptic_invariants(tau)[source]¶
>>> g2, g3 = CC.elliptic_invariants(1j) >>> CC.weierstrass_p_prime(0.25, 1j)**2; 4*CC.weierstrass_p(0.25, 1j)**3 - g2*CC.weierstrass_p(0.25, 1j) - g3 [15152.862386715 +/- 7.03e-10] [15152.86238672 +/- 5.09e-9]
- elliptic_roots(tau)[source]¶
>>> e1, e2, e3 = CC.elliptic_roots(1j) >>> g2, g3 = CC.elliptic_invariants(1j) >>> 4*e1**3 - g2*e1 - g3 [+/- 3.12e-11] >>> 4*e2**3 - g2*e2 - g3 [+/- 8.29e-12] >>> 4*e3**3 - g2*e3 - g3 [+/- 3.14e-11]
- class flint_ctypes.gr_elem(val=None, context=None, random=False)[source]¶
Base class for elements.
- __init__(val=None, context=None, random=False)[source]¶
>>> ZZ(QQ(1)) 1 >>> ZZ(QQ(1) / 3) Traceback (most recent call last): ... FlintDomainError: 1/3 is not defined in Integer ring (fmpz)
- parent()[source]¶
Return the parent object of this element.
>>> ZZ(0).parent() Integer ring (fmpz) >>> ZZ(0).parent() is ZZ True
- is_invertible()[source]¶
Return whether self has a multiplicative inverse in its domain.
>>> >>> ZZ(3).is_invertible() False >>> ZZ(-1).is_invertible() True
- divides(other)[source]¶
Return whether self divides other.
>>> ZZ(5).divides(10) True >>> ZZ(5).divides(12) False
- factor()[source]¶
Returns a factorization of self as a tuple (prefactor, factors, exponents).
>>> ZZ(-120).factor() (-1, [2, 3, 5], [3, 1, 1])
- is_square()[source]¶
Return whether self is a perfect square in its domain.
>>> ZZ(3).is_square() False >>> ZZ(4).is_square() True >>> QQbar(3).is_square() True
- inv()[source]¶
Multiplicative inverse of this element.
>>> QQ(3).inv() 1/3 >>> QQ(0).inv() Traceback (most recent call last): ... FlintDomainError: inv(x) is not an element of {Rational field (fmpq)} for {x = 0}
- sqrt()[source]¶
Square root of this element.
>>> ZZ(4).sqrt() 2 >>> ZZ(2).sqrt() Traceback (most recent call last): ... FlintDomainError: sqrt(x) is not an element of {Integer ring (fmpz)} for {x = 2} >>> QQbar(2).sqrt() Root a = 1.41421 of a^2-2 >>> (QQ(25)/16).sqrt() 5/4 >>> QQbar(-1).sqrt() Root a = 1.00000*I of a^2+1 >>> RR(-1).sqrt() Traceback (most recent call last): ... FlintDomainError: sqrt(x) is not an element of {Real numbers (arb, prec = 53)} for {x = -1.000000000000000} >>> RF(-1).sqrt() nan
- floor()[source]¶
Floor function: closest integer in the direction of \(-\infty\).
>>> (QQ(3) / 2).floor() 1 >>> (QQ(3) / 2).ceil() 2 >>> (QQ(3) / 2).nint() 2 >>> (QQ(3) / 2).trunc() 1
- ceil()[source]¶
Ceiling function: closest integer in the direction of \(+\infty\).
>>> (QQ(3) / 2).ceil() 2
- trunc()[source]¶
Truncate to integer: closest integer in the direction of zero.
>>> (QQ(3) / 2).trunc() 1
- nint()[source]¶
Nearest integer function: nearest integer, rounding to even on a tie.
>>> (QQ(3) / 2).nint() 2
- conj()[source]¶
Complex conjugate.
>>> QQbar.i().conj() Root a = -1.00000*I of a^2+1 >>> CC(-2).log().conj() ([0.693147180559945 +/- 4.12e-16] + [-3.141592653589793 +/- 3.39e-16]*I) >>> QQ(3).conj() 3
- im()[source]¶
Imaginary part.
>>> QQ(1).im() 0 >>> (QQbar(-1) ** (QQ(1) / 3)).im() Root a = 0.866025 of 4*a^2-3
- csgn()[source]¶
Real-valued extension of the sign function: gives the sign of the real part when nonzero, and the sign of the imaginary part when on the imaginary axis.
>>> QQbar(-10).sqrt().csgn() 1 >>> (-QQbar(-10).sqrt()).csgn() -1
- mul_2exp(other)[source]¶
Exact multiplication by a dyadic number \(2^y\).
>>> QQ(3).mul_2exp(5) 96 >>> QQ(3).mul_2exp(-5) 3/32 >>> ZZ(100).mul_2exp(-2) 25 >>> ZZ(100).mul_2exp(-3) Traceback (most recent call last): ... FlintDomainError: mul_2exp(x, y) is not an element of {Integer ring (fmpz)} for {x = 100}, {y = -3}
- exp()[source]¶
Exponential function.
>>> RR(1).exp() [2.718281828459045 +/- 5.41e-16] >>> RR_ca(1).exp() 2.71828 {a where a = 2.71828 [Exp(1)]} >>> QQ(0).exp() 1 >>> QQ(1).exp() Traceback (most recent call last): ... FlintUnableError: failed to compute exp(x) in {Rational field (fmpq)} for {x = 1}
- expm1()[source]¶
Exponential function minus 1.
>>> RR("1e-10").expm1() [1.000000000050000e-10 +/- 3.86e-26] >>> CC(RR("1e-10")).expm1() [1.000000000050000e-10 +/- 3.86e-26] >>> RF("1e-10").expm1() 1.000000000050000e-10 >>> CF(RF("1e-10")).expm1() 1.000000000050000e-10 >>> QQ(0).expm1() 0 >>> QQ(1).expm1() Traceback (most recent call last): ... FlintUnableError: failed to compute expm1(x) in {Rational field (fmpq)} for {x = 1}
- exp2()[source]¶
Exponential function with base 2.
>>> QQ(5).exp2() 32 >>> RF(0.5).exp2() 1.414213562373095
- exp10()[source]¶
Exponential function with base 10.
>>> QQ(5).exp2() 32 >>> RF(0.5).exp10() 3.162277660168380
- log()[source]¶
Natural logarithm.
>>> QQ(1).log() 0 >>> QQ(2).log() Traceback (most recent call last): ... FlintUnableError: failed to compute log(x) in {Rational field (fmpq)} for {x = 2} >>> RF(2).log() 0.6931471805599453
- log1p()[source]¶
Natural logarithm with one added to the argument.
>>> QQ(0).log1p() 0 >>> RF(-0.5).log1p() -0.6931471805599453 >>> RR(1).log1p() [0.693147180559945 +/- 4.12e-16]
- exp_pi_i()[source]¶
\(\exp(\pi i x)\) evaluated at self.
>>> (QQbar(1) / 3).exp_pi_i() Root a = 0.500000 + 0.866025*I of a^2-a+1 >>> (QQbar(2).sqrt()).exp_pi_i() Traceback (most recent call last): ... FlintDomainError: exp_pi_i(x) is not an element of {Complex algebraic numbers (qqbar)} for {x = Root a = 1.41421 of a^2-2}
- log_pi_i()[source]¶
\(\log(x) / (\pi i)\) evaluated at self.
>>> (QQbar(-1) ** (QQbar(7) / 5)).log_pi_i() -3/5 >>> (QQbar(1) / 2).log_pi_i() Traceback (most recent call last): ... FlintDomainError: log_pi_i(x) is not an element of {Complex algebraic numbers (qqbar)} for {x = 1/2}
- sin_pi()[source]¶
\(\sin(\pi x)\) evaluated at self.
>>> (QQbar(1) / 3).sin_pi() Root a = 0.866025 of 4*a^2-3
- tan_pi()[source]¶
\(\tan(\pi x)\) evaluated at self.
>>> (QQbar(1) / 3).tan_pi() Root a = 1.73205 of a^2-3
- cot_pi()[source]¶
\(\cot(\pi x)\) evaluated at self.
>>> (QQbar(1) / 3).cot_pi() Root a = 0.577350 of 3*a^2-1
- class flint_ctypes.gr_poly(val=None, context=None, random=False)[source]¶
- __init__(val=None, context=None, random=False)[source]¶
>>> ZZ(QQ(1)) 1 >>> ZZ(QQ(1) / 3) Traceback (most recent call last): ... FlintDomainError: 1/3 is not defined in Integer ring (fmpz)
- is_monic()[source]¶
>>> RRx([2,3,4]).is_monic() False >>> RRx([2,3,1]).is_monic() True >>> RRx([]).is_monic() False
- monic()[source]¶
Return self rescaled to a monic polynomial.
>>> f = RRx([1,RR.pi()]) >>> f.monic() [0.318309886183791 +/- 4.43e-16] + 1.000000000000000*x >>> RRx([]).monic() # the zero polynomial cannot be made monic Traceback (most recent call last): ... ValueError >>> (f - f).monic() # unknown whether it is the zero polynomial Traceback (most recent call last): ... NotImplementedError
- roots(domain=None)[source]¶
Computes the roots in the coefficient ring of this polynomial, returning a tuple (
roots,multiplicities). If the ring is not algebraically closed, the sum of multiplicities can be smaller than the degree of the polynomial. Ifdomainis given, returns roots in that ring instead.>>> (ZZx([3,2]) * ZZx([15,1])**2 * ZZx([-10,1])).roots() ([10, -15], [1, 2]) >>> ZZx([1]).roots() ([], [])
We consider roots of the zero polynomial to be ill-defined:
>>> ZZx([]).roots() Traceback (most recent call last): ... ValueError
We construct an integer polynomial with rational, real algebraic and complex algebraic roots and extract its roots over different domains:
>>> f = ZZx([-2,0,1]) * ZZx([1, 0, 1]) * ZZx([3, 2])**2 >>> f.roots() # integer roots (there are none) ([], []) >>> f.roots(domain=QQ) # rational roots ([-3/2], [2]) >>> f.roots(domain=AA) # real algebraic roots ([Root a = 1.41421 of a^2-2, Root a = -1.41421 of a^2-2, -3/2], [1, 1, 2]) >>> f.roots(domain=QQbar) # complex algebraic roots ([Root a = 1.00000*I of a^2+1, Root a = -1.00000*I of a^2+1, Root a = 1.41421 of a^2-2, Root a = -1.41421 of a^2-2, -3/2], [1, 1, 1, 1, 2]) >>> f.roots(domain=RR) # real ball roots ([[-1.414213562373095 +/- 4.89e-17], [1.414213562373095 +/- 4.89e-17], -1.500000000000000], [1, 1, 2]) >>> f.roots(domain=CC) # complex ball roots ([[-1.414213562373095 +/- 4.89e-17], [1.414213562373095 +/- 4.89e-17], 1.000000000000000*I, -1.000000000000000*I, -1.500000000000000], [1, 1, 1, 1, 2]) >>> f.roots(RF) # real floating-point roots ([-1.414213562373095, 1.414213562373095, -1.500000000000000], [1, 1, 2]) >>> f.roots(CF) # complex floating-point roots ([-1.414213562373095, 1.414213562373095, 1.000000000000000*I, -1.000000000000000*I, -1.500000000000000], [1, 1, 1, 1, 2])
- inv_series(n)[source]¶
Reciprocal of this polynomial viewed as a power series, truncated to length n.
>>> ZZx([1,2,3]).inv_series(10) 1 - 2*x + x^2 + 4*x^3 - 11*x^4 + 10*x^5 + 13*x^6 - 56*x^7 + 73*x^8 + 22*x^9 >>> ZZx([2,3,4]).inv_series(5) Traceback (most recent call last): ... FlintDomainError: f.inv_series(n) is not an element of {Ring of polynomials over Integer ring (fmpz)} for {f = 2 + 3*x + 4*x^2}, {n = 5} >>> QQx([2,3,4]).inv_series(5) (1/2) + (-3/4)*x + (1/8)*x^2 + (21/16)*x^3 + (-71/32)*x^4
- log_series(n)[source]¶
Logarithm of this polynomial viewed as a power series, truncated to length n.
>>> QQx([1,1]).log_series(8) x + (-1/2)*x^2 + (1/3)*x^3 + (-1/4)*x^4 + (1/5)*x^5 + (-1/6)*x^6 + (1/7)*x^7 >>> RRx([2,1]).log_series(3) [0.693147180559945 +/- 4.12e-16] + 0.5000000000000000*x - 0.1250000000000000*x^2 >>> RRx([0,0]).log_series(3) Traceback (most recent call last): ... FlintDomainError: f.log_series(n) is not an element of {Ring of polynomials over Real numbers (arb, prec = 53)} for {f = 0}, {n = 3}
- exp_series(n)[source]¶
Exponential of this polynomial viewed as a power series, truncated to length n.
>>> QQx([0,1]).exp_series(8) 1 + x + (1/2)*x^2 + (1/6)*x^3 + (1/24)*x^4 + (1/120)*x^5 + (1/720)*x^6 + (1/5040)*x^7 >>> QQx([1,1]).exp_series(2) Traceback (most recent call last): ... FlintUnableError: failed to compute f.exp_series(n) in {Ring of polynomials over Rational field (fmpq)} for {f = 1 + x}, {n = 2} >>> RRx([1,1]).exp_series(2) [2.718281828459045 +/- 5.41e-16] + [2.718281828459045 +/- 5.41e-16]*x >>> RRx([2,3]).log_series(3).exp_series(3) [2.000000000000000 +/- 6.97e-16] + [3.00000000000000 +/- 1.61e-15]*x + [+/- 1.49e-15]*x^2
- class flint_ctypes.DirichletGroup_dirichlet_char(q, **kwargs)[source]¶
Group of Dirichlet characters of given modulus.
>>> G = DirichletGroup(10) >>> G.q 10 >>> len(G) 4 >>> [G(i) for i in [1,3,7,9]] [chi_10(1, .), chi_10(3, .), chi_10(7, .), chi_10(9, .)]
- class flint_ctypes.Mat(element_domain, nrows=None, ncols=None)[source]¶
Parent class for matrix domains.
There are two kinds of matrix domains:
Mat(R), the set of matrices of any size over the domain R.
Mat(R, n, m), the set of n x m matrices over the domain R. If R is a ring and n = m, then this is also a ring.
While Mat(R) may be more convenient, e.g. for representing linear transformations of arbitrary dimension under a single parent, fixed-shape matrix domains have advantages such as allowing automatic conversion from scalars to scalar matrices of the right size.
>>> Mat(ZZ) Matrices (any shape) over Integer ring (fmpz) >>> Mat(ZZ, 2) Ring of 2 x 2 matrices over Integer ring (fmpz) >>> Mat(ZZ, 2, 3) Space of 2 x 3 matrices over Integer ring (fmpz) >>> Mat(ZZ)([[1, 2, 3], [4, 5, 6]]) [[1, 2, 3], [4, 5, 6]] >>> Mat(ZZ, 2, 2)(5) [[5, 0], [0, 5]]
- class flint_ctypes.gr_mat(*args, **kwargs)[source]¶
- __init__(*args, **kwargs)[source]¶
>>> ZZ(QQ(1)) 1 >>> ZZ(QQ(1) / 3) Traceback (most recent call last): ... FlintDomainError: 1/3 is not defined in Integer ring (fmpz)
- hessenberg(algorithm=None)[source]¶
Return this matrix reduced to upper Hessenberg form:
>>> B = Mat(QQ, 3, 3)([[4, 2, 3], [-1, 5, -3], [-4, 1, 2]]); >>> B.hessenberg() [[4, 14, 3], [-1, -7, -3], [0, 37, 14]]
Options: - algorithm:
None(default),"gauss"or"householder"
- eigenvalues(domain=None)[source]¶
Computes the eigenvalues in the coefficient ring of this matrix, returning a tuple (
eigenvalues,multiplicities). If the ring is not algebraically closed, the sum of multiplicities can be smaller than the dimension of the matrix. Ifdomainis given, returns eigenvalues in that ring instead.>>> Mat(ZZ)([[1,2],[3,4]]).eigenvalues() ([], []) >>> Mat(ZZ)([[1,2],[3,-4]]).eigenvalues() ([2, -5], [1, 1]) >>> Mat(ZZ)([[1,2],[3,4]]).eigenvalues(domain=QQbar) ([Root a = 5.37228 of a^2-5*a-2, Root a = -0.372281 of a^2-5*a-2], [1, 1]) >>> Mat(ZZ)([[1,2],[3,4]]).eigenvalues(domain=RR) ([[-0.3722813232690143 +/- 3.01e-17], [5.372281323269014 +/- 3.31e-16]], [1, 1])
The matrix must be square:
>>> Mat(ZZ)([[1,2,3],[4,5,6]]).eigenvalues() Traceback (most recent call last): ... ValueError
- diagonalization()[source]¶
Matrix diagonalization: returns (D, L, R) where D is a vector of eigenvalues, LAR = diag(D) and LR = 1.
>>> A = Mat(QQ)([[1,2],[-1,4]]) >>> D, L, R = A.diagonalization() >>> L*A*R [[3, 0], [0, 2]] >>> D [3, 2] >>> L*R [[1, 0], [0, 1]]
>>> A = Mat(CC)([[1,2],[-1,4]]) >>> D, L, R = A.diagonalization() >>> D [([2.00000000000000 +/- 1.86e-15] + [+/- 1.86e-15]*I), ([3.00000000000000 +/- 2.90e-15] + [+/- 1.86e-15]*I)] >>> L*A*R [[([2.00000000000 +/- 1.10e-12] + [+/- 1.08e-12]*I), ([+/- 1.44e-12] + [+/- 1.42e-12]*I)], [([+/- 9.76e-13] + [+/- 9.63e-13]*I), ([3.00000000000 +/- 1.27e-12] + [+/- 1.25e-12]*I)]] >>> L*R [[([1.00000000000 +/- 3.26e-13] + [+/- 3.20e-13]*I), ([+/- 3.72e-13] + [+/- 3.67e-13]*I)], [([+/- 2.77e-13] + [+/- 2.73e-13]*I), ([1.00000000000 +/- 3.17e-13] + [+/- 3.13e-13]*I)]]
>>> A = Mat(CF)([[1,2],[-1,4]]) >>> D, L, R = A.diagonalization() >>> D [2.000000000000000, 3.000000000000000] >>> L*A*R [[2.000000000000000, -8.275113827716402e-16], [0, 3.000000000000000]] >>> L*R [[1.000000000000000, -8.275113803054639e-17], [0, 1.000000000000000]]
- class flint_ctypes.gr_vec(*args, **kwargs)[source]¶
- flint_ctypes.PolynomialRing¶
alias of
flint_ctypes.PolynomialRing_gr_poly
- flint_ctypes.PowerSeriesRing¶
- flint_ctypes.PowerSeriesModRing¶
- flint_ctypes.ModularGroup¶
alias of
flint_ctypes.ModularGroup_psl2z
- flint_ctypes.DirichletGroup¶
- flint_ctypes.SymmetricGroup¶
alias of
flint_ctypes.SymmetricGroup_perm