Transcendental functions

^

The expression x{^}n is powering. If the exponent is an integer, then exact operations are performed using binary (left-shift) powering techniques. In particular, in this case x cannot be a vector or matrix unless it is a square matrix (invertible if the exponent is negative). If x is a p-adic number, its precision will increase if v_p(n) > 0. Powering a binary quadratic form (types t_QFI and t_QFR) returns a reduced representative of the class, provided the input is reduced. In particular, x{^}1 is identical to x.

PARI is able to rewrite the multiplication x * x of two identical objects as x^2, or sqr(x). Here, identical means the operands are two different labels referencing the same chunk of memory; no equality test is performed. This is no longer true when more than two arguments are involved.

If the exponent is not of type integer, this is treated as a transcendental function (see Section [Label: se:trans]), and in particular has the effect of componentwise powering on vector or matrices.

As an exception, if the exponent is a rational number p/q and x an integer modulo a prime or a p-adic number, return a solution y of y^q = x^p if it exists. Currently, q must not have large prime factors. Beware that

  ? Mod(7,19)^(1/2)
  %1 = Mod(11, 19) /* is any square root */
  ? sqrt(Mod(7,19))
  %2 = Mod(8, 19)  /* is the smallest square root */
  ? Mod(7,19)^(3/5)
  %3 = Mod(1, 19)
  ? %3^(5/3)
  %4 = Mod(1, 19)  /* Mod(7,19) is just another cubic root */

If the exponent is a negative integer, an inverse must be computed. For non-invertible t_INTMOD, this will fail and implicitly exhibit a non trivial factor of the modulus:

  ? Mod(4,6)^(-1)
    ***   at top-level: Mod(4,6)^(-1)
    ***                         ^-----
    *** _^_: impossible inverse modulo: Mod(2, 6).

(Here, a factor 2 is obtained directly. In general, take the gcd of the representative and the modulus.) This is most useful when performing complicated operations modulo an integer N whose factorization is unknown. Either the computation succeeds and all is well, or a factor d is discovered and the computation may be restarted modulo d or N/d.

For non-invertible t_POLMOD, this will fail without exhibiting a factor.

  ? Mod(x^2, x^3-x)^(-1)
    ***   at top-level: Mod(x^2,x^3-x)^(-1)
    ***                               ^-----
    *** _^_: non-invertible polynomial in RgXQ_inv.
  ? a = Mod(3,4)*y^3 + Mod(1,4); b = y^6+y^5+y^4+y^3+y^2+y+1;
  ? Mod(a, b)^(-1);
    ***   at top-level: Mod(a,b)^(-1)
    ***                         ^-----
    *** _^_: impossible inverse modulo: Mod(0, 4).

In fact the latter polynomial is invertible, but the algorithm used (subresultant) assumes the base ring is a domain. If it is not the case, as here for Z/4Z, a result will be correct but chances are an error will occur first. In this specific case, one should work with 2-adics. In general, one can try the following approach

  ? inversemod(a, b) =
  { my(m);
    m = polsylvestermatrix(polrecip(a), polrecip(b));
    m = matinverseimage(m, matid(#m)[,1]);
    Polrev( vecextract(m, Str("..", poldegree(b))), variable(b) )
  }
  ? inversemod(a,b)
  %2 = Mod(2,4)*y^5 + Mod(3,4)*y^3 + Mod(1,4)*y^2 + Mod(3,4)*y + Mod(2,4)

This is not guaranteed to work either since it must invert pivots. See Section [Label: se:linear_algebra].

The library syntax is GEN gpow(GEN x, GEN n, long prec) for x{^}n.


Catalan

Catalan's constant G = sum_{n >= 0}((-1)^n)/((2n+1)^2) = 0.91596.... Note that Catalan is one of the few reserved names which cannot be used for user variables.

The library syntax is GEN mpcatalan(long prec).


Euler

Euler's constant gamma = 0.57721.... Note that Euler is one of the few reserved names which cannot be used for user variables.

The library syntax is GEN mpeuler(long prec).


I

The complex number sqrt{-1}.

The library syntax is GEN gen_I().


Pi

The constant Pi (3.14159...). Note that Pi is one of the few reserved names which cannot be used for user variables.

The library syntax is GEN mppi(long prec).


abs(x)

Absolute value of x (modulus if x is complex). Rational functions are not allowed. Contrary to most transcendental functions, an exact argument is not converted to a real number before applying abs and an exact result is returned if possible.

  ? abs(-1)
  %1 = 1
  ? abs(3/7 + 4/7*I)
  %2 = 5/7
  ? abs(1 + I)
  %3 = 1.414213562373095048801688724

If x is a polynomial, returns -x if the leading coefficient is real and negative else returns x. For a power series, the constant coefficient is considered instead.

The library syntax is GEN gabs(GEN x, long prec).


acos(x)

Principal branch of {cos}^{-1}(x) = -i log (x + isqrt{1-x^2}). In particular, {Re(acos}(x)) belongs to [0,Pi] and if x belongs to R and |x| > 1, then {acos}(x) is complex. The branch cut is in two pieces: ]- oo ,-1] , continuous with quadrant II, and [1,+ oo [, continuous with quadrant IV. We have {acos}(x) = Pi/2 - {asin}(x) for all x.

The library syntax is GEN gacos(GEN x, long prec).


acosh(x)

Principal branch of {cosh}^{-1}(x) = 2 log(sqrt{(x+1)/2} + sqrt{(x-1)/2}). In particular, {Re}({acosh}(x)) >= 0 and {In}({acosh}(x)) belongs to ]-Pi,Pi]0; if x belongs to R and x < 1, then {acosh}(x) is complex.

The library syntax is GEN gacosh(GEN x, long prec).


agm(x,y)

Arithmetic-geometric mean of x and y. In the case of complex or negative numbers, the optimal AGM is returned (the largest in absolute value over all choices of the signs of the square roots). p-adic or power series arguments are also allowed. Note that a p-adic agm exists only if x/y is congruent to 1 modulo p (modulo 16 for p = 2). x and y cannot both be vectors or matrices.

The library syntax is GEN agm(GEN x, GEN y, long prec).


arg(x)

Argument of the complex number x, such that -Pi < {arg}(x) <= Pi.

The library syntax is GEN garg(GEN x, long prec).


asin(x)

Principal branch of {sin}^{-1}(x) = -i log(ix + sqrt{1 - x^2}). In particular, {Re(asin}(x)) belongs to [-Pi/2,Pi/2] and if x belongs to R and |x| > 1 then {asin}(x) is complex. The branch cut is in two pieces: ]- oo ,-1], continuous with quadrant II, and [1,+ oo [ continuous with quadrant IV. The function satisfies i {asin}(x) = {asinh}(ix).

The library syntax is GEN gasin(GEN x, long prec).


asinh(x)

Principal branch of {sinh}^{-1}(x) = log(x + sqrt{1+x^2}). In particular {Im(asinh}(x)) belongs to [-Pi/2,Pi/2]. The branch cut is in two pieces: [-i oo ,-i], continuous with quadrant III and [i,+i oo [ continuous with quadrant I.

The library syntax is GEN gasinh(GEN x, long prec).


atan(x)

Principal branch of {tan}^{-1}(x) = log ((1+ix)/(1-ix)) / 2i. In particular the real part of {atan}(x)) belongs to ]-Pi/2,Pi/2[. The branch cut is in two pieces: ]-i oo ,-i[, continuous with quadrant IV, and ]i,+i oo [ continuous with quadrant II. The function satisfies i {atan}(x) = -i{atanh}(ix) for all x != ± i.

The library syntax is GEN gatan(GEN x, long prec).


atanh(x)

Principal branch of {tanh}^{-1}(x) = log ((1+x)/(1-x)) / 2. In particular the imaginary part of {atanh}(x) belongs to [-Pi/2,Pi/2]; if x belongs to R and |x| > 1 then {atanh}(x) is complex.

The library syntax is GEN gatanh(GEN x, long prec).


bernfrac(x)

Bernoulli number B_x, where B_0 = 1, B_1 = -1/2, B_2 = 1/6,..., expressed as a rational number. The argument x should be of type integer.

The library syntax is GEN bernfrac(long x).


bernpol(n, {v = 'x})

Bernoulli polynomial B_n in variable v.

  ? bernpol(1)
  %1 = x - 1/2
  ? bernpol(3)
  %2 = x^3 - 3/2*x^2 + 1/2*x

The library syntax is GEN bernpol(long n, long v = -1), where v is a variable number.


bernreal(x)

Bernoulli number B_x, as bernfrac, but B_x is returned as a real number (with the current precision).

The library syntax is GEN bernreal(long x, long prec).


bernvec(x)

Creates a vector containing, as rational numbers, the Bernoulli numbers B_0, B_2,..., B_{2x}. This routine is obsolete. Use bernfrac instead each time you need a Bernoulli number in exact form.

Note. This routine is implemented using repeated independent calls to bernfrac, which is faster than the standard recursion in exact arithmetic. It is only kept for backward compatibility: it is not faster than individual calls to bernfrac, its output uses a lot of memory space, and coping with the index shift is awkward.

The library syntax is GEN bernvec(long x).


besselh1(nu,x)

H^1-Bessel function of index nu and argument x.

The library syntax is GEN hbessel1(GEN nu, GEN x, long prec).


besselh2(nu,x)

H^2-Bessel function of index nu and argument x.

The library syntax is GEN hbessel2(GEN nu, GEN x, long prec).


besseli(nu,x)

I-Bessel function of index nu and argument x. If x converts to a power series, the initial factor (x/2)^nu/Gamma(nu+1) is omitted (since it cannot be represented in PARI when nu is not integral).

The library syntax is GEN ibessel(GEN nu, GEN x, long prec).


besselj(nu,x)

J-Bessel function of index nu and argument x. If x converts to a power series, the initial factor (x/2)^nu/Gamma(nu+1) is omitted (since it cannot be represented in PARI when nu is not integral).

The library syntax is GEN jbessel(GEN nu, GEN x, long prec).


besseljh(n,x)

J-Bessel function of half integral index. More precisely, besseljh(n,x) computes J_{n+1/2}(x) where n must be of type integer, and x is any element of C. In the present version 2.7.0, this function is not very accurate when x is small.

The library syntax is GEN jbesselh(GEN n, GEN x, long prec).


besselk(nu,x)

K-Bessel function of index nu and argument x.

The library syntax is GEN kbessel(GEN nu, GEN x, long prec).


besseln(nu,x)

N-Bessel function of index nu and argument x.

The library syntax is GEN nbessel(GEN nu, GEN x, long prec).


cos(x)

Cosine of x.

The library syntax is GEN gcos(GEN x, long prec).


cosh(x)

Hyperbolic cosine of x.

The library syntax is GEN gcosh(GEN x, long prec).


cotan(x)

Cotangent of x.

The library syntax is GEN gcotan(GEN x, long prec).


dilog(x)

Principal branch of the dilogarithm of x, i.e. analytic continuation of the power series log_2(x) = sum_{n >= 1}x^n/n^2.

The library syntax is GEN dilog(GEN x, long prec).


eint1(x,{n})

Exponential integral int_x^ oo (e^{-t})/(t)dt = incgam(0, x), where the latter expression extends the function definition from real x > 0 to all complex x != 0.

If n is present, we must have x > 0; the function returns the n-dimensional vector [eint1(x),...,eint1(nx)]. Contrary to other transcendental functions, and to the default case (n omitted), the values are correct up to a bounded absolute, rather than relative, error 10^-n, where n is precision(x) if x is a t_REAL and defaults to realprecision otherwise. (In the most important application, to the computation of L-functions via approximate functional equations, those values appear as weights in long sums and small individual relative errors are less useful than controlling the absolute error.) This is faster than repeatedly calling eint1(i * x), but less precise.

The library syntax is GEN veceint1(GEN x, GEN n = NULL, long prec). Also available is GEN eint1(GEN x, long prec).


erfc(x)

Complementary error function, analytic continuation of (2/sqrtPi)int_x^ oo e^{-t^2}dt = incgam(1/2,x^2)/sqrtPi, where the latter expression extends the function definition from real x to all complex x != 0.

The library syntax is GEN gerfc(GEN x, long prec).


eta(z,{flag = 0})

Variants of Dedekind's eta function. If flag = 0, return prod_{n = 1}^ oo (1-q^n), where q depends on x in the following way:

* q = e^{2iPi x} if x is a complex number (which must then have positive imaginary part); notice that the factor q^{1/24} is missing!

* q = x if x is a t_PADIC, or can be converted to a power series (which must then have positive valuation).

If flag is non-zero, x is converted to a complex number and we return the true eta function, q^{1/24}prod_{n = 1}^ oo (1-q^n), where q = e^{2iPi x}.

The library syntax is GEN eta0(GEN z, long flag, long prec).

Also available is GEN trueeta(GEN x, long prec) (flag = 1).


exp(x)

Exponential of x. p-adic arguments with positive valuation are accepted.

The library syntax is GEN gexp(GEN x, long prec). For a t_PADIC x, the function GEN Qp_exp(GEN x) is also available.


expm1(x)

Return exp(x)-1, computed in a way that is also accurate when the real part of x is near 0. Only accept real or complex arguments. A naive direct computation would suffer from catastrophic cancellation; PARI's direct computation of exp(x) alleviates this well known problem at the expense of computing exp(x) to a higher accuracy when x is small. Using expm1 is recommanded instead:

  ? default(realprecision, 10000); x = 1e-100;
  ? a = expm1(x);
  time = 4 ms.
  ? b = exp(x)-1;
  time = 28 ms.
  ? default(realprecision, 10040); x = 1e-100;
  ? c = expm1(x);  \\ reference point
  ? abs(a-c)/c  \\ relative error in expm1(x)
  %7 = 0.E-10017
  ? abs(b-c)/c  \\ relative error in exp(x)-1
  %8 = 1.7907031188259675794 E-9919
As the example above shows, when x is near 0, expm1 is both faster and more accurate than exp(x)-1.

The library syntax is GEN gexpm1(GEN x, long prec).


gamma(s)

For s a complex number, evaluates Euler's gamma function Gamma(s) = int_0^ oo t^{s-1}exp(-t)dt. Error if s is a non-positive integer, where Gamma has a pole.

For s a t_PADIC, evaluates the Morita gamma function at s, that is the unique continuous p-adic function on the p-adic integers extending Gamma_p(k) = (-1)^k prod_{j < k}'j, where the prime means that p does not divide j.

  ? gamma(1/4 + O(5^10))
  %1= 1 + 4*5 + 3*5^4 + 5^6 + 5^7 + 4*5^9 + O(5^10)
  ? algdep(%,4)
  %2 = x^4 + 4*x^2 + 5

The library syntax is GEN ggamma(GEN s, long prec). For a t_PADIC x, the function GEN Qp_gamma(GEN x) is also available.


gammah(x)

Gamma function evaluated at the argument x+1/2.

The library syntax is GEN ggammah(GEN x, long prec).


hyperu(a,b,x)

U-confluent hypergeometric function with parameters a and b. The parameters a and b can be complex but the present implementation requires x to be positive.

The library syntax is GEN hyperu(GEN a, GEN b, GEN x, long prec).


incgam(s,x,{g})

Incomplete gamma function int_x^ oo e^{-t}t^{s-1}dt, extended by analytic continuation to all complex x, s not both 0. The relative error is bounded in terms of the precision of s (the accuracy of x is ignored when determining the output precision). When g is given, assume that g = Gamma(s). For small |x|, this will speed up the computation.

The library syntax is GEN incgam0(GEN s, GEN x, GEN g = NULL, long prec). Also available is GEN incgam(GEN s, GEN x, long prec).


incgamc(s,x)

Complementary incomplete gamma function. The arguments x and s are complex numbers such that s is not a pole of Gamma and |x|/(|s|+1) is not much larger than 1 (otherwise the convergence is very slow). The result returned is int_0^x e^{-t}t^{s-1}dt.

The library syntax is GEN incgamc(GEN s, GEN x, long prec).


lambertw(y)

Lambert W function, solution of the implicit equation xe^x = y, for y > 0.

The library syntax is GEN glambertW(GEN y, long prec).


lngamma(x)

Principal branch of the logarithm of the gamma function of x. This function is analytic on the complex plane with non-positive integers removed, and can have much larger arguments than gamma itself.

For x a power series such that x(0) is not a pole of gamma, compute the Taylor expansion. (PARI only knows about regular power series and can't include logarithmic terms.)

  ? lngamma(1+x+O(x^2))
  %1 = -0.57721566490153286060651209008240243104*x + O(x^2)
  ? lngamma(x+O(x^2))
   ***   at top-level: lngamma(x+O(x^2))
   ***                 ^-----------------
   *** lngamma: domain error in lngamma: valuation != 0
  ? lngamma(-1+x+O(x^2))
   *** lngamma: Warning: normalizing a series with 0 leading term.
   ***   at top-level: lngamma(-1+x+O(x^2))
   ***                 ^--------------------
   *** lngamma: domain error in intformal: residue(series, pole) != 0

The library syntax is GEN glngamma(GEN x, long prec).


log(x)

Principal branch of the natural logarithm of x belongs to C^*, i.e. such that {Im(log}(x)) belongs to ]-Pi,Pi]. The branch cut lies along the negative real axis, continuous with quadrant 2, i.e. such that lim_{b\to 0^+} log (a+bi) = log a for a belongs to R^*. The result is complex (with imaginary part equal to Pi) if x belongs to R and x < 0. In general, the algorithm uses the formula log(x) ~ (Pi)/(2{agm}(1, 4/s)) - m log 2, if s = x 2^m is large enough. (The result is exact to B bits provided s > 2^{B/2}.) At low accuracies, the series expansion near 1 is used.

p-adic arguments are also accepted for x, with the convention that log(p) = 0. Hence in particular exp(log(x))/x is not in general equal to 1 but to a (p-1)-th root of unity (or ±1 if p = 2) times a power of p.

The library syntax is GEN glog(GEN x, long prec). For a t_PADIC x, the function GEN Qp_log(GEN x) is also available.


polylog(m,x,{flag = 0})

One of the different polylogarithms, depending on flag:

If flag = 0 or is omitted: m-th polylogarithm of x, i.e. analytic continuation of the power series {Li}_m(x) = sum_{n >= 1}x^n/n^m (x < 1). Uses the functional equation linking the values at x and 1/x to restrict to the case |x| <= 1, then the power series when |x|^2 <= 1/2, and the power series expansion in log(x) otherwise.

Using flag, computes a modified m-th polylogarithm of x. We use Zagier's notations; let Re_m denote Re or Im depending on whether m is odd or even:

If flag = 1: compute ~ D_m(x), defined for |x| <= 1 by Re_m(sum_{k = 0}^{m-1} ((-log|x|)^k)/(k!){Li}_{m-k}(x) +((-log|x|)^{m-1})/(m!)log|1-x|).

If flag = 2: compute D_m(x), defined for |x| <= 1 by Re_m(sum_{k = 0}^{m-1}((-log|x|)^k)/(k!){Li}_{m-k}(x) -(1)/(2)((-log|x|)^m)/(m!)).

If flag = 3: compute P_m(x), defined for |x| <= 1 by Re_m(sum_{k = 0}^{m-1}(2^kB_k)/(k!)(log|x|)^k{Li}_{m-k}(x) -(2^{m-1}B_m)/(m!)(log|x|)^m).

These three functions satisfy the functional equation f_m(1/x) = (-1)^{m-1}f_m(x).

The library syntax is GEN polylog0(long m, GEN x, long flag, long prec). Also available is GEN gpolylog(long m, GEN x, long prec) (flag = 0).


psi(x)

The psi-function of x, i.e. the logarithmic derivative Gamma'(x)/Gamma(x).

The library syntax is GEN gpsi(GEN x, long prec).


sin(x)

Sine of x.

The library syntax is GEN gsin(GEN x, long prec).


sinh(x)

Hyperbolic sine of x.

The library syntax is GEN gsinh(GEN x, long prec).


sqr(x)

Square of x. This operation is not completely straightforward, i.e. identical to x * x, since it can usually be computed more efficiently (roughly one-half of the elementary multiplications can be saved). Also, squaring a 2-adic number increases its precision. For example,

  ? (1 + O(2^4))^2
  %1 = 1 + O(2^5)
  ? (1 + O(2^4)) * (1 + O(2^4))
  %2 = 1 + O(2^4)

Note that this function is also called whenever one multiplies two objects which are known to be identical, e.g. they are the value of the same variable, or we are computing a power.

  ? x = (1 + O(2^4)); x * x
  %3 = 1 + O(2^5)
  ? (1 + O(2^4))^4
  %4 = 1 + O(2^6)

(note the difference between %2 and %3 above).

The library syntax is GEN gsqr(GEN x).


sqrt(x)

Principal branch of the square root of x, defined as sqrt{x} = exp(log x / 2). In particular, we have {Arg}({sqrt}(x)) belongs to ]-Pi/2, Pi/2], and if x belongs to R and x < 0, then the result is complex with positive imaginary part.

Intmod a prime p, t_PADIC and t_FFELT are allowed as arguments. In the first 2 cases (t_INTMOD, t_PADIC), the square root (if it exists) which is returned is the one whose first p-adic digit is in the interval [0,p/2]. For other arguments, the result is undefined.

The library syntax is GEN gsqrt(GEN x, long prec). For a t_PADIC x, the function GEN Qp_sqrt(GEN x) is also available.


sqrtn(x,n,{&z})

Principal branch of the nth root of x, i.e. such that {Arg}({sqrt}(x)) belongs to ]-Pi/n, Pi/n]. Intmod a prime and p-adics are allowed as arguments.

If z is present, it is set to a suitable root of unity allowing to recover all the other roots. If it was not possible, z is set to zero. In the case this argument is present and no square root exist, 0 is returned instead or raising an error.

  ? sqrtn(Mod(2,7), 2)
  %1 = Mod(4, 7)
  ? sqrtn(Mod(2,7), 2, &z); z
  %2 = Mod(6, 7)
  ? sqrtn(Mod(2,7), 3)
    ***   at top-level: sqrtn(Mod(2,7),3)
    ***                 ^-----------------
    *** sqrtn: nth-root does not exist in gsqrtn.
  ? sqrtn(Mod(2,7), 3,  &z)
  %2 = 0
  ? z
  %3 = 0

The following script computes all roots in all possible cases:

  sqrtnall(x,n)=
  { my(V,r,z,r2);
    r = sqrtn(x,n, &z);
    if (!z, error("Impossible case in sqrtn"));
    if (type(x) == "t_INTMOD" || type(x)=="t_PADIC",
      r2 = r*z; n = 1;
      while (r2!=r, r2*=z;n++));
    V = vector(n); V[1] = r;
    for(i=2, n, V[i] = V[i-1]*z);
    V
  }
  addhelp(sqrtnall,"sqrtnall(x,n):compute the vector of nth-roots of x");

The library syntax is GEN gsqrtn(GEN x, GEN n, GEN *z = NULL, long prec). If x is a t_PADIC, the function GEN Qp_sqrt(GEN x, GEN n, GEN *z) is also available.


tan(x)

Tangent of x.

The library syntax is GEN gtan(GEN x, long prec).


tanh(x)

Hyperbolic tangent of x.

The library syntax is GEN gtanh(GEN x, long prec).


teichmuller(x)

Teichmüller character of the p-adic number x, i.e. the unique (p-1)-th root of unity congruent to x / p^{v_p(x)} modulo p.

The library syntax is GEN teich(GEN x).


theta(q,z)

Jacobi sine theta-function theta_1(z, q) = 2q^{1/4} sum_{n >= 0} (-1)^n q^{n(n+1)} sin((2n+1)z).

The library syntax is GEN theta(GEN q, GEN z, long prec).


thetanullk(q,k)

k-th derivative at z = 0 of theta(q,z).

The library syntax is GEN thetanullk(GEN q, long k, long prec).

GEN vecthetanullk(GEN q, long k, long prec) returns the vector of all (d^itheta)/(dz^i)(q,0) for all odd i = 1, 3,..., 2k-1. GEN vecthetanullk_tau(GEN tau, long k, long prec) returns vecthetanullk_tau at q = exp(2iPi tau).


weber(x,{flag = 0})

One of Weber's three f functions. If flag = 0, returns f(x) = exp(-iPi/24).eta((x+1)/2)/eta(x) {such that} j = (f^{24}-16)^3/f^{24}, where j is the elliptic j-invariant (see the function ellj). If flag = 1, returns f_1(x) = eta(x/2)/eta(x) {such that} j = (f_1^{24}+16)^3/f_1^{24}. Finally, if flag = 2, returns f_2(x) = sqrt{2}eta(2x)/eta(x) {such that} j = (f_2^{24}+16)^3/f_2^{24}. Note the identities f^8 = f_1^8+f_2^8 and ff_1f_2 = sqrt2.

The library syntax is GEN weber0(GEN x, long flag, long prec). Also available are GEN weberf(GEN x, long prec), GEN weberf1(GEN x, long prec) and GEN weberf2(GEN x, long prec).


zeta(s)

For s a complex number, Riemann's zeta function zeta(s) = sum_{n >= 1}n^{-s}, computed using the Euler-Maclaurin summation formula, except when s is of type integer, in which case it is computed using Bernoulli numbers for s <= 0 or s > 0 and even, and using modular forms for s > 0 and odd.

For s a p-adic number, Kubota-Leopoldt zeta function at s, that is the unique continuous p-adic function on the p-adic integers that interpolates the values of (1 - p^{-k}) zeta(k) at negative integers k such that k = 1 (mod p-1) (resp. k is odd) if p is odd (resp. p = 2).

The library syntax is GEN gzeta(GEN s, long prec).