Code coverage tests

This page documents the degree to which the PARI/GP source code is tested by our public test suite, distributed with the source distribution in directory src/test/. This is measured by the gcov utility; we then process gcov output using the lcov frond-end.

We test a few variants depending on Configure flags on the pari.math.u-bordeaux.fr machine (x86_64 architecture), and agregate them in the final report:

The target is 90% coverage for all mathematical modules (given that branches depending on DEBUGLEVEL or DEBUGMEM are not covered). This script is run to produce the results below.

LCOV - code coverage report
Current view: top level - basemath - trans2.c (source / functions) Hit Total Coverage
Test: PARI/GP v2.8.0 lcov report (development 19352-1b11b25) Lines: 958 1014 94.5 %
Date: 2016-08-25 06:11:27 Functions: 62 63 98.4 %
Legend: Lines: hit not hit

          Line data    Source code
       1             : /* Copyright (C) 2000  The PARI group.
       2             : 
       3             : This file is part of the PARI/GP package.
       4             : 
       5             : PARI/GP is free software; you can redistribute it and/or modify it under the
       6             : terms of the GNU General Public License as published by the Free Software
       7             : Foundation. It is distributed in the hope that it will be useful, but WITHOUT
       8             : ANY WARRANTY WHATSOEVER.
       9             : 
      10             : Check the License for details. You should have received a copy of it, along
      11             : with the package; see the file 'COPYING'. If not, write to the Free Software
      12             : Foundation, Inc., 51 Franklin Street, Fifth Floor, Boston, MA 02110-1301 USA. */
      13             : 
      14             : /********************************************************************/
      15             : /**                                                                **/
      16             : /**                   TRANSCENDENTAL FUNCTIONS                     **/
      17             : /**                          (part 2)                              **/
      18             : /**                                                                **/
      19             : /********************************************************************/
      20             : #include "pari.h"
      21             : #include "paripriv.h"
      22             : 
      23             : GEN
      24      184290 : trans_fix_arg(long *prec, GEN *s0, GEN *sig, GEN *tau, pari_sp *av, GEN *res)
      25             : {
      26             :   GEN s, p1;
      27             :   long l;
      28      184290 :   if (typ(*s0)==t_COMPLEX && gequal0(gel(*s0,2))) *s0 = gel(*s0,1);
      29      184290 :   s = *s0;
      30      184290 :   l = precision(s); if (!l) l = *prec;
      31      184290 :   if (l < LOWDEFAULTPREC) l = LOWDEFAULTPREC;
      32      184290 :   *res = cgetc(l); *av = avma;
      33      184290 :   if (typ(s) == t_COMPLEX)
      34             :   { /* s = sig + i t */
      35      175231 :     s = cxtofp(s, l+EXTRAPRECWORD);
      36      175231 :     *sig = gel(s,1);
      37      175231 :     *tau = gel(s,2);
      38             :   }
      39             :   else /* real number */
      40             :   {
      41        9059 :     *sig = s = gtofp(s, l+EXTRAPRECWORD);
      42        9059 :     *tau = gen_0;
      43        9059 :     p1 = trunc2nr(s, 0);
      44        9059 :     if (!signe(subri(s,p1))) *s0 = p1;
      45             :   }
      46      184290 :   *prec = l; return s;
      47             : }
      48             : 
      49             : /********************************************************************/
      50             : /**                                                                **/
      51             : /**                          ARCTANGENT                            **/
      52             : /**                                                                **/
      53             : /********************************************************************/
      54             : static GEN
      55      394299 : mpatan(GEN x)
      56             : {
      57      394299 :   long l, l1, l2, n, m, i, lp, e, s, sx = signe(x);
      58             :   pari_sp av0, av;
      59             :   double alpha, beta, delta;
      60             :   GEN y, p1, p2, p3, p4, p5, unr;
      61             :   int inv;
      62             : 
      63      394299 :   if (!sx) return real_0_bit(expo(x));
      64      394236 :   l = lp = realprec(x);
      65      394236 :   if (absrnz_equal1(x)) { /* |x| = 1 */
      66        3206 :     y = Pi2n(-2, l+EXTRAPRECWORD); if (sx < 0) setsigne(y,-1);
      67        3206 :     return y;
      68             :   }
      69      391030 :   if (l > AGM_ATAN_LIMIT)
      70             :   {
      71         175 :     av = avma; y = logagmcx(mkcomplex(gen_1, x), l);
      72         175 :     return gerepileuptoleaf(av, gel(y,2));
      73             :   }
      74             : 
      75      390855 :   e = expo(x); inv = (e >= 0); /* = (|x| > 1 ) */
      76      390855 :   if (e > 0) lp += nbits2extraprec(e);
      77             : 
      78      390855 :   y = cgetr(lp); av0 = avma;
      79      390855 :   p1 = rtor(x, l+EXTRAPRECWORD); setabssign(p1); /* p1 = |x| */
      80      390855 :   if (inv) p1 = invr(p1);
      81      390855 :   e = expo(p1);
      82      390855 :   if (e < -100)
      83        2255 :     alpha = 1.65149612947 - e; /* log_2(Pi) - e */
      84             :   else
      85      388600 :     alpha = log2(M_PI / atan(rtodbl(p1)));
      86      390855 :   beta = (double)(prec2nbits(l)>>1);
      87      390855 :   delta = 1 + beta - alpha/2;
      88      390855 :   if (delta <= 0) { n = 1; m = 0; }
      89             :   else
      90             :   {
      91      388646 :     double fi = alpha-2;
      92      388646 :     if (delta >= fi*fi)
      93             :     {
      94      388386 :       double t = 1 + sqrt(delta);
      95      388386 :       n = (long)t;
      96      388386 :       m = (long)(t - fi);
      97             :     }
      98             :     else
      99             :     {
     100         260 :       n = (long)(1+beta/fi);
     101         260 :       m = 0;
     102             :     }
     103             :   }
     104      390855 :   l2 = l + nbits2extraprec(m);
     105      390855 :   p2 = rtor(p1, l2); av = avma;
     106     3878878 :   for (i=1; i<=m; i++)
     107             :   {
     108     3488023 :     p5 = addsr(1, sqrr(p2)); setprec(p5,l2);
     109     3488023 :     p5 = addsr(1, sqrtr_abs(p5)); setprec(p5,l2);
     110     3488023 :     affrr(divrr(p2,p5), p2); avma = av;
     111             :   }
     112      390855 :   p3 = sqrr(p2); l1 = minss(LOWDEFAULTPREC+EXTRAPRECWORD, l2); /* l1 increases to l2 */;
     113      390855 :   unr = real_1(l2); setprec(unr,l1);
     114      390855 :   p4 = cgetr(l2); setprec(p4,l1);
     115      390855 :   affrr(divru(unr,2*n+1), p4);
     116      390855 :   s = 0; e = expo(p3); av = avma;
     117     4334139 :   for (i = n; i > 1; i--) /* n >= 1. i = 1 done outside for efficiency */
     118             :   {
     119     3943284 :     setprec(p3,l1); p5 = mulrr(p4,p3);
     120     3943284 :     l1 += dvmdsBIL(s - e, &s); if (l1 > l2) l1 = l2;
     121     3943284 :     setprec(unr,l1); p5 = subrr(divru(unr,2*i-1), p5);
     122     3943284 :     setprec(p4,l1); affrr(p5,p4); avma = av;
     123             :   }
     124      390855 :   setprec(p3, l2); p5 = mulrr(p4,p3); /* i = 1 */
     125      390855 :   setprec(unr,l2); p4 = subrr(unr, p5);
     126             : 
     127      390855 :   p4 = mulrr(p2,p4); shiftr_inplace(p4, m);
     128      390855 :   if (inv) p4 = subrr(Pi2n(-1, lp), p4);
     129      390855 :   if (sx < 0) togglesign(p4);
     130      390855 :   affrr_fixlg(p4,y); avma = av0; return y;
     131             : }
     132             : 
     133             : GEN
     134       24568 : gatan(GEN x, long prec)
     135             : {
     136             :   pari_sp av;
     137             :   GEN a, y;
     138             : 
     139       24568 :   switch(typ(x))
     140             :   {
     141       16420 :     case t_REAL: return mpatan(x);
     142             :     case t_COMPLEX: /* atan(x) = -i atanh(ix) */
     143        7630 :       if (ismpzero(gel(x,2))) return gatan(gel(x,1), prec);
     144        7399 :       av = avma; return gerepilecopy(av, mulcxmI(gatanh(mulcxI(x),prec)));
     145             :     default:
     146         518 :       av = avma; if (!(y = toser_i(x))) break;
     147          28 :       if (valp(y) < 0) pari_err_DOMAIN("atan","valuation", "<", gen_0, x);
     148          21 :       if (lg(y)==2) return gerepilecopy(av, y);
     149             :       /* lg(y) > 2 */
     150          14 :       a = integser(gdiv(derivser(y), gaddsg(1,gsqr(y))));
     151          14 :       if (!valp(y)) a = gadd(a, gatan(gel(y,2),prec));
     152          14 :       return gerepileupto(av, a);
     153             :   }
     154         490 :   return trans_eval("atan",gatan,x,prec);
     155             : }
     156             : /********************************************************************/
     157             : /**                                                                **/
     158             : /**                             ARCSINE                            **/
     159             : /**                                                                **/
     160             : /********************************************************************/
     161             : /* |x| < 1, x != 0 */
     162             : static GEN
     163          98 : mpasin(GEN x) {
     164          98 :   pari_sp av = avma;
     165          98 :   GEN z, a = sqrtr(subsr(1, sqrr(x)));
     166          98 :   if (realprec(x) > AGM_ATAN_LIMIT)
     167             :   {
     168           7 :     z = logagmcx(mkcomplex(a,x), realprec(x));
     169           7 :     z = gel(z,2);
     170             :   }
     171             :   else
     172          91 :     z = mpatan(divrr(x, a));
     173          98 :   return gerepileuptoleaf(av, z);
     174             : }
     175             : 
     176             : static GEN mpacosh(GEN x);
     177             : GEN
     178        8155 : gasin(GEN x, long prec)
     179             : {
     180             :   long sx;
     181             :   pari_sp av;
     182             :   GEN a, y, p1;
     183             : 
     184        8155 :   switch(typ(x))
     185             :   {
     186         469 :     case t_REAL: sx = signe(x);
     187         469 :       if (!sx) return real_0_bit(expo(x));
     188         462 :       if (absrnz_equal1(x)) { /* |x| = 1 */
     189          28 :         if (sx > 0) return Pi2n(-1, realprec(x)); /* 1 */
     190          14 :         y = Pi2n(-1, realprec(x)); setsigne(y, -1); return y; /* -1 */
     191             :       }
     192         434 :       if (expo(x) < 0) return mpasin(x);
     193         336 :       y = cgetg(3,t_COMPLEX);
     194         336 :       gel(y,1) = Pi2n(-1, realprec(x));
     195         336 :       gel(y,2) = mpacosh(x);
     196         336 :       if (sx < 0) togglesign(gel(y,1)); else togglesign(gel(y,2));
     197         336 :       return y;
     198             : 
     199             :     case t_COMPLEX: /* asin(z) = -i asinh(iz) */
     200        7623 :       if (ismpzero(gel(x,2))) return gasin(gel(x,1), prec);
     201        7392 :       av = avma;
     202        7392 :       return gerepilecopy(av, mulcxmI(gasinh(mulcxI(x), prec)));
     203             :     default:
     204          63 :       av = avma; if (!(y = toser_i(x))) break;
     205          42 :       if (gequal0(y)) return gerepilecopy(av, y);
     206             :       /* lg(y) > 2*/
     207          35 :       if (valp(y) < 0) pari_err_DOMAIN("asin","valuation", "<", gen_0, x);
     208          28 :       p1 = gsubsg(1,gsqr(y));
     209          28 :       if (gequal0(p1))
     210             :       {
     211          21 :         GEN t = Pi2n(-1,prec);
     212          21 :         if (gsigne(gel(y,2)) < 0) setsigne(t, -1);
     213          21 :         return gerepileupto(av, scalarser(t, varn(y), valp(p1)>>1));
     214             :       }
     215           7 :       p1 = gdiv(derivser(y), gsqrt(p1,prec));
     216           7 :       a = integser(p1);
     217           7 :       if (!valp(y)) a = gadd(a, gasin(gel(y,2),prec));
     218           7 :       return gerepileupto(av, a);
     219             :   }
     220          21 :   return trans_eval("asin",gasin,x,prec);
     221             : }
     222             : /********************************************************************/
     223             : /**                                                                **/
     224             : /**                             ARCCOSINE                          **/
     225             : /**                                                                **/
     226             : /********************************************************************/
     227             : static GEN
     228          14 : acos0(long e) { return Pi2n(-1, nbits2prec(e<0? -e: 1)); }
     229             : 
     230             : /* |x| < 1, x != 0 */
     231             : static GEN
     232         105 : mpacos(GEN x)
     233             : {
     234         105 :   pari_sp av = avma;
     235         105 :   GEN z, a = sqrtr(subsr(1, sqrr(x)));
     236         105 :   if (realprec(x) > AGM_ATAN_LIMIT)
     237             :   {
     238          14 :     z = logagmcx(mkcomplex(x,a), realprec(x));
     239          14 :     z = gel(z,2);
     240             :   }
     241             :   else {
     242          91 :     z = mpatan(divrr(a, x));
     243          91 :     if (signe(x) < 0) z = addrr(mppi(realprec(z)), z);
     244             :   }
     245         105 :   return gerepileuptoleaf(av, z);
     246             : }
     247             : 
     248             : GEN
     249        7917 : gacos(GEN x, long prec)
     250             : {
     251             :   long sx;
     252             :   pari_sp av;
     253             :   GEN a, y, p1;
     254             : 
     255        7917 :   switch(typ(x))
     256             :   {
     257         245 :     case t_REAL: sx = signe(x);
     258         245 :       if (!sx) return acos0(expo(x));
     259         238 :       if (absrnz_equal1(x)) /* |x| = 1 */
     260          14 :         return sx > 0? real_0_bit( -(bit_prec(x)>>1) ) : mppi(realprec(x));
     261         224 :       if (expo(x) < 0) return mpacos(x);
     262             : 
     263         168 :       y = cgetg(3,t_COMPLEX); p1 = mpacosh(x);
     264         168 :       if (sx < 0) { gel(y,1) = mppi(realprec(x)); togglesign(p1); }
     265          84 :       else gel(y,1) = gen_0;
     266         168 :       gel(y,2) = p1; return y;
     267             : 
     268             :     case t_COMPLEX:
     269        7623 :       if (ismpzero(gel(x,2))) return gacos(gel(x,1), prec);
     270        7392 :       av = avma;
     271        7392 :       p1 = gadd(x, mulcxI(gsqrt(gsubsg(1,gsqr(x)), prec)));
     272        7392 :       y = glog(p1,prec); /* log(x + I*sqrt(1-x^2)) */
     273        7392 :       return gerepilecopy(av, mulcxmI(y));
     274             :     default:
     275          49 :       av = avma; if (!(y = toser_i(x))) break;
     276          35 :       if (valp(y) < 0) pari_err_DOMAIN("acos","valuation", "<", gen_0, x);
     277          28 :       if (lg(y) > 2)
     278             :       {
     279          21 :         p1 = gsubsg(1,gsqr(y));
     280          21 :         if (gequal0(p1)) { avma = av; return zeroser(varn(y), valp(p1)>>1); }
     281           7 :         p1 = integser(gdiv(gneg(derivser(y)), gsqrt(p1,prec)));
     282             :         /*y(t) = 1+O(t)*/
     283           7 :         if (gequal1(gel(y,2)) && !valp(y)) return gerepileupto(av, p1);
     284             :       }
     285           7 :       else p1 = y;
     286          14 :       a = (lg(y)==2 || valp(y))? Pi2n(-1, prec): gacos(gel(y,2),prec);
     287          14 :       return gerepileupto(av, gadd(a,p1));
     288             :   }
     289          14 :   return trans_eval("acos",gacos,x,prec);
     290             : }
     291             : /********************************************************************/
     292             : /**                                                                **/
     293             : /**                            ARGUMENT                            **/
     294             : /**                                                                **/
     295             : /********************************************************************/
     296             : 
     297             : /* we know that x and y are not both 0 */
     298             : static GEN
     299      377817 : mparg(GEN x, GEN y)
     300             : {
     301      377817 :   long prec, sx = signe(x), sy = signe(y);
     302             :   GEN z;
     303             : 
     304      377817 :   if (!sy)
     305             :   {
     306         113 :     if (sx > 0) return real_0_bit(expo(y) - expo(x));
     307          56 :     return mppi(realprec(x));
     308             :   }
     309      377704 :   prec = realprec(y); if (prec < realprec(x)) prec = realprec(x);
     310      377704 :   if (!sx)
     311             :   {
     312           7 :     z = Pi2n(-1, prec); if (sy < 0) setsigne(z,-1);
     313           7 :     return z;
     314             :   }
     315             : 
     316      377697 :   if (expo(x)-expo(y) > -2)
     317             :   {
     318      327654 :     z = mpatan(divrr(y,x)); if (sx > 0) return z;
     319       66576 :     return addrr_sign(z, signe(z), mppi(prec), sy);
     320             :   }
     321       50043 :   z = mpatan(divrr(x,y));
     322       50043 :   return addrr_sign(z, -signe(z), Pi2n(-1, prec), sy);
     323             : }
     324             : 
     325             : static GEN
     326      755634 : rfix(GEN x,long prec)
     327             : {
     328      755634 :   switch(typ(x))
     329             :   {
     330        2754 :     case t_INT: return itor(x, prec);
     331        2177 :     case t_FRAC: return fractor(x, prec);
     332      750703 :     case t_REAL: break;
     333           0 :     default: pari_err_TYPE("rfix (conversion to t_REAL)",x);
     334             :   }
     335      750703 :   return x;
     336             : }
     337             : 
     338             : static GEN
     339      377817 : cxarg(GEN x, GEN y, long prec)
     340             : {
     341      377817 :   pari_sp av = avma;
     342      377817 :   x = rfix(x,prec);
     343      377817 :   y = rfix(y,prec); return gerepileuptoleaf(av, mparg(x,y));
     344             : }
     345             : 
     346             : GEN
     347      378951 : garg(GEN x, long prec)
     348             : {
     349      378951 :   if (gequal0(x)) pari_err_DOMAIN("arg", "argument", "=", gen_0, x);
     350      378951 :   switch(typ(x))
     351             :   {
     352        1134 :     case t_REAL: prec = realprec(x); /* fall through */
     353        1134 :     case t_INT: case t_FRAC: return (gsigne(x)>0)? real_0(prec): mppi(prec);
     354      377817 :     case t_COMPLEX: return cxarg(gel(x,1),gel(x,2),prec);
     355             :   }
     356           0 :   return trans_eval("arg",garg,x,prec);
     357             : }
     358             : 
     359             : /********************************************************************/
     360             : /**                                                                **/
     361             : /**                      HYPERBOLIC COSINE                         **/
     362             : /**                                                                **/
     363             : /********************************************************************/
     364             : 
     365             : static GEN
     366        3675 : mpcosh(GEN x)
     367             : {
     368             :   pari_sp av;
     369             :   GEN z;
     370             : 
     371        3675 :   if (!signe(x)) { /* 1 + x */
     372           0 :     long e = expo(x);
     373           0 :     return e >= 0? real_0_bit(e): real_1_bit(-e);
     374             :   }
     375        3675 :   av = avma;
     376        3675 :   z = mpexp(x); z = addrr(z, invr(z)); shiftr_inplace(z, -1);
     377        3675 :   return gerepileuptoleaf(av, z);
     378             : }
     379             : 
     380             : GEN
     381        3738 : gcosh(GEN x, long prec)
     382             : {
     383             :   pari_sp av;
     384             :   GEN y, p1;
     385             : 
     386        3738 :   switch(typ(x))
     387             :   {
     388        3675 :     case t_REAL: return mpcosh(x);
     389             :     case t_COMPLEX:
     390           7 :       if (isintzero(gel(x,1))) return gcos(gel(x,2),prec);
     391             :       /* fall through */
     392             :     case t_PADIC:
     393           7 :       av = avma; p1 = gexp(x,prec); p1 = gadd(p1, ginv(p1));
     394           7 :       return gerepileupto(av, gmul2n(p1,-1));
     395             :     default:
     396          49 :       av = avma; if (!(y = toser_i(x))) break;
     397          28 :       if (gequal0(y) && valp(y) == 0) return gerepilecopy(av, y);
     398          28 :       p1 = gexp(y,prec); p1 = gadd(p1, ginv(p1));
     399          28 :       return gerepileupto(av, gmul2n(p1,-1));
     400             :   }
     401          21 :   return trans_eval("cosh",gcosh,x,prec);
     402             : }
     403             : /********************************************************************/
     404             : /**                                                                **/
     405             : /**                       HYPERBOLIC SINE                          **/
     406             : /**                                                                **/
     407             : /********************************************************************/
     408             : 
     409             : static GEN
     410          28 : mpsinh(GEN x)
     411             : {
     412             :   pari_sp av;
     413          28 :   long ex = expo(x), lx;
     414             :   GEN z, res;
     415             : 
     416          28 :   if (!signe(x)) return real_0_bit(ex);
     417          28 :   lx = realprec(x); res = cgetr(lx); av = avma;
     418          28 :   if (ex < 1 - BITS_IN_LONG) x = rtor(x, lx + nbits2extraprec(-ex)-1);
     419          28 :   z = mpexp(x); z = subrr(z, invr(z)); shiftr_inplace(z, -1);
     420          28 :   affrr(z, res); avma = av; return res;
     421             : }
     422             : 
     423             : GEN
     424        1260 : gsinh(GEN x, long prec)
     425             : {
     426             :   pari_sp av;
     427             :   GEN y, p1;
     428             : 
     429        1260 :   switch(typ(x))
     430             :   {
     431          28 :     case t_REAL: return mpsinh(x);
     432             :     case t_COMPLEX:
     433           7 :       if (isintzero(gel(x,1))) retmkcomplex(gen_0, gsin(gel(x,2),prec));
     434             :       /* fall through */
     435             :     case t_PADIC:
     436           7 :       av = avma; p1 = gexp(x,prec); p1 = gsub(p1, ginv(p1));
     437           7 :       return gerepileupto(av, gmul2n(p1,-1));
     438             :     default:
     439        1218 :       av = avma; if (!(y = toser_i(x))) break;
     440        1190 :       if (gequal0(y) && valp(y) == 0) return gerepilecopy(av, y);
     441        1190 :       p1 = gexp(y, prec); p1 = gsub(p1, ginv(p1));
     442        1190 :       return gerepileupto(av, gmul2n(p1,-1));
     443             :   }
     444          28 :   return trans_eval("sinh",gsinh,x,prec);
     445             : }
     446             : /********************************************************************/
     447             : /**                                                                **/
     448             : /**                      HYPERBOLIC TANGENT                        **/
     449             : /**                                                                **/
     450             : /********************************************************************/
     451             : 
     452             : static GEN
     453       76951 : mptanh(GEN x)
     454             : {
     455       76951 :   long lx, s = signe(x);
     456             :   GEN y;
     457             : 
     458       76951 :   if (!s) return real_0_bit(expo(x));
     459       76951 :   lx = realprec(x);
     460       76951 :   if (abscmprr(x, stor(prec2nbits(lx), LOWDEFAULTPREC)) >= 0) {
     461       24792 :     y = real_1(lx);
     462             :   } else {
     463       52159 :     pari_sp av = avma;
     464       52159 :     long ex = expo(x);
     465             :     GEN t;
     466       52159 :     if (ex < 1 - BITS_IN_LONG) x = rtor(x, lx + nbits2extraprec(-ex)-1);
     467       52159 :     t = exp1r_abs(gmul2n(x,1)); /* exp(|2x|) - 1 */
     468       52159 :     y = gerepileuptoleaf(av, divrr(t, addsr(2,t)));
     469             :   }
     470       76951 :   if (s < 0) togglesign(y); /* tanh is odd */
     471       76951 :   return y;
     472             : }
     473             : 
     474             : GEN
     475       77035 : gtanh(GEN x, long prec)
     476             : {
     477             :   pari_sp av;
     478             :   GEN y, t;
     479             : 
     480       77035 :   switch(typ(x))
     481             :   {
     482       76951 :     case t_REAL: return mptanh(x);
     483             :     case t_COMPLEX:
     484          14 :       if (isintzero(gel(x,1))) retmkcomplex(gen_0, gtan(gel(x,2),prec));
     485             :       /* fall through */
     486             :     case t_PADIC:
     487          14 :       av = avma;
     488          14 :       t = gexp(gmul2n(x,1),prec);
     489          14 :       t = gdivsg(-2, gaddgs(t,1));
     490          14 :       return gerepileupto(av, gaddsg(1,t));
     491             :     default:
     492          63 :       av = avma; if (!(y = toser_i(x))) break;
     493          28 :       if (gequal0(y)) return gerepilecopy(av, y);
     494          14 :       t = gexp(gmul2n(y, 1),prec);
     495          14 :       t = gdivsg(-2, gaddgs(t,1));
     496          14 :       return gerepileupto(av, gaddsg(1,t));
     497             :   }
     498          35 :   return trans_eval("tanh",gtanh,x,prec);
     499             : }
     500             : 
     501             : static GEN
     502           0 : mpcotanh(GEN x)
     503             : {
     504           0 :   long lx, s = signe(x);
     505             :   GEN y;
     506             : 
     507           0 :   if (!s) pari_err_DOMAIN("cotan", "argument", "=", gen_0, x);
     508             : 
     509           0 :   lx = realprec(x);
     510           0 :   if (abscmprr(x, stor(prec2nbits(lx), LOWDEFAULTPREC)) >= 0) {
     511           0 :     y = real_1(lx);
     512             :   } else {
     513           0 :     pari_sp av = avma;
     514           0 :     long ex = expo(x);
     515             :     GEN t;
     516           0 :     if (ex < 1 - BITS_IN_LONG) x = rtor(x, lx + nbits2extraprec(-ex)-1);
     517           0 :     t = exp1r_abs(gmul2n(x,1)); /* exp(|2x|) - 1 */
     518           0 :     y = gerepileuptoleaf(av, divrr(addsr(2,t), t));
     519             :   }
     520           0 :   if (s < 0) togglesign(y); /* cotanh is odd */
     521           0 :   return y;
     522             : }
     523             : 
     524             : GEN
     525          42 : gcotanh(GEN x, long prec)
     526             : {
     527             :   pari_sp av;
     528             :   GEN y, t;
     529             : 
     530          42 :   switch(typ(x))
     531             :   {
     532           0 :     case t_REAL: return mpcotanh(x);
     533             :     case t_COMPLEX:
     534           0 :       if (isintzero(gel(x,1))) retmkcomplex(gen_0, gcotan(gel(x,2),prec));
     535             :       /* fall through */
     536             :     case t_PADIC:
     537           7 :       av = avma;
     538           7 :       t = gexpm1(gmul2n(x,1),prec);
     539           7 :       return gerepileupto(av, gaddsg(1, gdivsg(2,t)));
     540             :     default:
     541          35 :       av = avma; if (!(y = toser_i(x))) break;
     542          28 :       if (gequal0(y)) return gerepilecopy(av, y);
     543          14 :       t = gexpm1(gmul2n(y,1),prec);
     544          14 :       return gerepileupto(av, gaddsg(1, gdivsg(2,t)));
     545             :   }
     546           7 :   return trans_eval("cotanh",gcotanh,x,prec);
     547             : }
     548             : 
     549             : /********************************************************************/
     550             : /**                                                                **/
     551             : /**                     AREA HYPERBOLIC SINE                       **/
     552             : /**                                                                **/
     553             : /********************************************************************/
     554             : 
     555             : /* x != 0 */
     556             : static GEN
     557         469 : mpasinh(GEN x)
     558             : {
     559             :   GEN z, res;
     560             :   pari_sp av;
     561         469 :   long lx = realprec(x), ex = expo(x);
     562             : 
     563         469 :   res = cgetr(lx); av = avma;
     564         469 :   if (ex < 1 - BITS_IN_LONG) x = rtor(x, lx + nbits2extraprec(-ex)-1);
     565         469 :   z = logr_abs( addrr_sign(x,1, sqrtr_abs( addrs(sqrr(x), 1) ), 1) );
     566         469 :   if (signe(x) < 0) togglesign(z);
     567         469 :   affrr(z, res); avma = av; return res;
     568             : }
     569             : 
     570             : GEN
     571       15652 : gasinh(GEN x, long prec)
     572             : {
     573             :   pari_sp av;
     574             :   GEN a, y, p1;
     575             : 
     576       15652 :   switch(typ(x))
     577             :   {
     578             :     case t_REAL:
     579         476 :       if (!signe(x)) return rcopy(x);
     580         469 :       return mpasinh(x);
     581             : 
     582             :     case t_COMPLEX:
     583       15015 :       if (ismpzero(gel(x,2))) return gasinh(gel(x,1), prec);
     584       14560 :       av = avma;
     585       14560 :       if (ismpzero(gel(x,1))) /* avoid cancellation */
     586         224 :         return gerepilecopy(av, mulcxI(gasin(gel(x,2), prec)));
     587       14336 :       p1 = gadd(x, gsqrt(gaddsg(1,gsqr(x)), prec));
     588       14336 :       y = glog(p1,prec); /* log (x + sqrt(1+x^2)) */
     589       14336 :       return gerepileupto(av, y);
     590             :     default:
     591         161 :       av = avma; if (!(y = toser_i(x))) break;
     592         133 :       if (gequal0(y)) return gerepilecopy(av, y);
     593         126 :       if (valp(y) < 0) pari_err_DOMAIN("asinh","valuation", "<", gen_0, x);
     594         119 :       p1 = gaddsg(1,gsqr(y));
     595         119 :       if (gequal0(p1))
     596             :       {
     597          14 :         GEN t = PiI2n(-1,prec);
     598          14 :         if ( gsigne(imag_i(gel(y,2))) < 0 ) setsigne(gel(t,2), -1);
     599          14 :         return gerepileupto(av, scalarser(t, varn(y), valp(p1)>>1));
     600             :       }
     601         105 :       p1 = gdiv(derivser(y), gsqrt(p1,prec));
     602         105 :       a = integser(p1);
     603         105 :       if (!valp(y)) a = gadd(a, gasinh(gel(y,2),prec));
     604         105 :       return gerepileupto(av, a);
     605             :   }
     606          28 :   return trans_eval("asinh",gasinh,x,prec);
     607             : }
     608             : /********************************************************************/
     609             : /**                                                                **/
     610             : /**                     AREA HYPERBOLIC COSINE                     **/
     611             : /**                                                                **/
     612             : /********************************************************************/
     613             : 
     614             : /* |x| >= 1, return ach(|x|) */
     615             : static GEN
     616         700 : mpacosh(GEN x)
     617             : {
     618         700 :   pari_sp av = avma;
     619             :   GEN z;
     620         700 :   if (absrnz_equal1(x)) return real_0_bit(- bit_prec(x) >> 1);
     621         693 :   z = logr_abs( addrr_sign(x, 1, sqrtr( subrs(sqrr(x), 1) ), 1) );
     622         693 :   return gerepileuptoleaf(av, z);
     623             : }
     624             : 
     625             : GEN
     626        7959 : gacosh(GEN x, long prec)
     627             : {
     628             :   pari_sp av;
     629             :   GEN y, p1;
     630             : 
     631        7959 :   switch(typ(x))
     632             :   {
     633             :     case t_REAL: {
     634         259 :       long s = signe(x), e = expo(x);
     635             :       GEN a, b;
     636         259 :       if (s > 0 && e >= 0) return mpacosh(x);
     637             :       /* x < 1 */
     638         147 :       y = cgetg(3,t_COMPLEX); a = gen_0;
     639         147 :       if (s == 0) b = acos0(e);
     640         140 :       else if (e < 0) b = mpacos(x); /* -1 < x < 1 */
     641             :       else {
     642          91 :         if (!absrnz_equal1(x)) a = mpacosh(x);
     643          91 :         b = mppi(realprec(x));
     644             :       }
     645         147 :       gel(y,1) = a;
     646         147 :       gel(y,2) = b; return y;
     647             :     }
     648             :     case t_COMPLEX:
     649        7623 :       if (ismpzero(gel(x,2))) return gacosh(gel(x,1), prec);
     650        7392 :       av = avma;
     651        7392 :       p1 = gadd(x, gsqrt(gaddsg(-1,gsqr(x)), prec));
     652        7392 :       y = glog(p1,prec); /* log(x + sqrt(x^2-1)) */
     653        7392 :       if (signe(real_i(y)) < 0) y = gneg(y);
     654        7392 :       return gerepileupto(av, y);
     655             :     default: {
     656             :       GEN a;
     657             :       long v;
     658          77 :       av = avma; if (!(y = toser_i(x))) break;
     659          49 :       v = valp(y);
     660          49 :       if (v < 0) pari_err_DOMAIN("acosh","valuation", "<", gen_0, x);
     661          42 :       if (gequal0(y))
     662             :       {
     663           7 :         if (!v) return gerepilecopy(av, y);
     664           7 :         return gerepileupto(av, gadd(y, PiI2n(-1, prec)));
     665             :       }
     666          35 :       p1 = gsubgs(gsqr(y),1);
     667          35 :       if (gequal0(p1)) { avma = av; return zeroser(varn(y), valp(p1)>>1); }
     668          21 :       p1 = gdiv(derivser(y), gsqrt(p1,prec));
     669          21 :       a = integser(p1);
     670          21 :       if (v)
     671           7 :         p1 = PiI2n(-1, prec); /* I Pi/2 */
     672             :       else
     673             :       {
     674          14 :         p1 = gel(y,2); if (gequal1(p1)) return gerepileupto(av,a);
     675           7 :         p1 = gacosh(p1, prec);
     676             :       }
     677          14 :       return gerepileupto(av, gadd(p1,a));
     678             :     }
     679             :   }
     680          28 :   return trans_eval("acosh",gacosh,x,prec);
     681             : }
     682             : /********************************************************************/
     683             : /**                                                                **/
     684             : /**                     AREA HYPERBOLIC TANGENT                    **/
     685             : /**                                                                **/
     686             : /********************************************************************/
     687             : 
     688             : /* |x| < 1, x != 0 */
     689             : static GEN
     690          98 : mpatanh(GEN x)
     691             : {
     692          98 :   pari_sp av = avma;
     693          98 :   long ex = expo(x);
     694             :   GEN z;
     695          98 :   if (ex < 1 - BITS_IN_LONG) x = rtor(x, realprec(x) + nbits2extraprec(-ex)-1);
     696          98 :   z = invr( subsr(1,x) ); shiftr_inplace(z, 1); /* 2/(1-x)*/
     697          98 :   z = logr_abs( addrs(z,-1) );
     698          98 :   shiftr_inplace(z, -1); return gerepileuptoleaf(av, z);
     699             : }
     700             : 
     701             : GEN
     702       15540 : gatanh(GEN x, long prec)
     703             : {
     704             :   long sx;
     705             :   pari_sp av;
     706             :   GEN a, y, z;
     707             : 
     708       15540 :   switch(typ(x))
     709             :   {
     710             :     case t_REAL:
     711         469 :       sx = signe(x);
     712         469 :       if (!sx) return real_0_bit(expo(x));
     713         462 :       if (expo(x) < 0) return mpatanh(x);
     714             : 
     715         364 :       y = cgetg(3,t_COMPLEX);
     716         364 :       av = avma;
     717         364 :       z = subrs(x,1);
     718         364 :       if (!signe(z)) pari_err_DOMAIN("atanh", "argument", "=", gen_1, x);
     719         350 :       z = invr(z); shiftr_inplace(z, 1); /* 2/(x-1)*/
     720         350 :       z = addrs(z,1);
     721         350 :       if (!signe(z)) pari_err_DOMAIN("atanh", "argument", "=", gen_m1, x);
     722         336 :       z = logr_abs(z);
     723         336 :       shiftr_inplace(z, -1); /* (1/2)log((1+x)/(x-1)) */
     724         336 :       gel(y,1) = gerepileuptoleaf(av, z);
     725         336 :       gel(y,2) = Pi2n(-1, realprec(x));
     726         336 :       if (sx > 0) togglesign(gel(y,2));
     727         336 :       return y;
     728             : 
     729             :     case t_COMPLEX: /* 2/(1-z) - 1 = (1+z) / (1-z) */
     730       15022 :       if (ismpzero(gel(x,2))) return gatanh(gel(x,1), prec);
     731       14567 :       av = avma; z = glog( gaddgs(gdivsg(2,gsubsg(1,x)),-1), prec );
     732       14567 :       return gerepileupto(av, gmul2n(z,-1));
     733             : 
     734             :     default:
     735          49 :       av = avma; if (!(y = toser_i(x))) break;
     736          28 :       if (valp(y) < 0) pari_err_DOMAIN("atanh","valuation", "<", gen_0, x);
     737          21 :       z = gdiv(derivser(y), gsubsg(1,gsqr(y)));
     738          14 :       a = integser(z);
     739          14 :       if (!valp(y)) a = gadd(a, gatanh(gel(y,2),prec));
     740          14 :       return gerepileupto(av, a);
     741             :   }
     742          21 :   return trans_eval("atanh",gatanh,x,prec);
     743             : }
     744             : /********************************************************************/
     745             : /**                                                                **/
     746             : /**               CACHE BERNOULLI NUMBERS B_2k                     **/
     747             : /**                                                                **/
     748             : /********************************************************************/
     749             : static GEN
     750     1346477 : bern(GEN B, long pr)
     751             : {
     752     1346477 :   if (typ(B) != t_REAL) return fractor(B, pr);
     753      704476 :   if (realprec(B) < pr) return rtor(B,pr);
     754           0 :   return B;
     755             : }
     756             : static const long BERN_MINNB = 5;
     757             : /* need B[2..2*nb] at least prec accuracy. If prec = 0, compute exactly */
     758             : void
     759      182312 : mpbern(long nb, long prec)
     760             : {
     761      182312 :   const pari_sp av = avma;
     762      182312 :   long n, pr, n_is_small = 1, lbern = 0;
     763             :   GEN B;
     764             :   pari_timer T;
     765             : 
     766             :   /* pr = accuracy for computation, prec = required accuracy for result */
     767      182312 :   if (prec)
     768             :   {
     769      182235 :     pr = prec;
     770      182235 :     incrprec(pr);
     771             :   }
     772             :   else
     773          77 :     pr = prec = LONG_MAX; /* oo */
     774      182312 :   if (nb < BERN_MINNB) nb = BERN_MINNB;
     775      182312 :   if (bernzone)
     776             :   { /* don't recompute known Bernoulli */
     777             :     long i, min, max;
     778      182159 :     lbern = lg(bernzone);
     779      182159 :     if (lbern-1 < nb)
     780         447 :     { min = lbern-1; max = nb; }
     781             :     else
     782      181712 :     { min = nb; max = lbern-1; }
     783             :     /* skip B_2, ..., B_{2*MINNB}, always included as t_FRAC */
     784     7870845 :     for (n = BERN_MINNB+1; n <= min; n++)
     785             :     {
     786     7688797 :       GEN c = gel(bernzone,n);
     787             :       /* also stop if prec = 0 (compute exactly) */
     788     7688797 :       if (typ(c) == t_REAL && realprec(c) < prec) break;
     789             :     }
     790             :     /* B[1..n-1] are OK */
     791      363863 :     if (n > nb) return;
     792         455 :     B = cgetg_block(max+1, t_VEC);
     793         455 :     for (i = 1; i < n; i++) gel(B,i) = gel(bernzone,i);
     794             :     /* keep B[nb+1..max] */
     795         455 :     for (i = nb+1; i <= max; i++) gel(B,i) = gel(bernzone,i);
     796             :   }
     797             :   else
     798             :   {
     799         153 :     B = cgetg_block(nb+1, t_VEC);
     800         153 :     gel(B,1) = gclone(mkfrac(gen_1, utoipos(6)));
     801         153 :     gel(B,2) = gclone(mkfrac(gen_m1, utoipos(30)));
     802         153 :     gel(B,3) = gclone(mkfrac(gen_1, utoipos(42)));
     803         153 :     gel(B,4) = gel(B,2);
     804         153 :     gel(B,5) = gclone(mkfrac(utoipos(5), utoipos(66)));
     805         153 :     n = BERN_MINNB+1;
     806             :   }
     807         608 :   avma = av;
     808         608 :   if (DEBUGLEVEL) {
     809           0 :     err_printf("caching Bernoulli numbers 2 to 2*%ld, prec = %ld\n",
     810             :                nb, prec == LONG_MAX? 0: prec);
     811           0 :     timer_start(&T);
     812             :   }
     813             : 
     814             :   /* B_{2n} = (2n-1) / (4n+2) -
     815             :    * sum_{a = 1}^{n-1} (2n)...(2n+2-2a) / (2...(2a-1)2a) B_{2a} */
     816         608 :   n_is_small = 1;
     817       15458 :   for (; n <= nb; n++, avma = av)
     818             :   { /* compute and store B[n] = B_{2n} */
     819             :     GEN S;
     820       14850 :     if (n < lbern)
     821             :     {
     822        5941 :       GEN b = gel(bernzone,n);
     823        5941 :       if (typ(b)!=t_REAL || realprec(b)>=prec) { gel(B,n) = b; continue; }
     824             :     }
     825             :     /* Not cached, must compute */
     826             :     /* huge accuracy ? May as well compute exactly */
     827       18261 :     if (n_is_small && (prec == LONG_MAX ||
     828        3411 :                        2*n * log((double)2*n) < prec2nbits_mul(prec, LOG2)))
     829        2892 :       S = bernfrac_using_zeta(2*n);
     830             :     else
     831             :     {
     832             : #ifdef LONG_IS_64BIT
     833        9978 :       const ulong mul_overflow = 3037000500UL;
     834             : #else
     835        1980 :       const ulong mul_overflow = 46341UL;
     836             : #endif
     837       11958 :       ulong u = 8, v = 5, a = n-1, b = 2*n-3;
     838       11958 :       n_is_small = 0;
     839       11958 :       S = bern(gel(B,a), pr); /* B_2a */
     840             :       for (;;)
     841             :       { /* b = 2a-1, u = 2v-2, 2a + v = 2n+3 */
     842     1346477 :         if (a == 1) { S = mulri(S, muluu(u,v)); break; } /* a=b=1, v=2n+1, u=4n */
     843             :         /* beware overflow */
     844     1334519 :         S = (v <= mul_overflow)? mulru(S, u*v): mulri(S, muluu(u,v));
     845     1334519 :         S = (a <= mul_overflow)? divru(S, a*b): divri(S, muluu(a,b));
     846     1334519 :         u += 4; v += 2; a--; b -= 2;
     847     1334519 :         S = addrr(bern(gel(B,a), pr), S);
     848     1334519 :         if ((a & 127) == 0) S = gerepileuptoleaf(av, S);
     849     1334519 :       }
     850       11958 :       S = divru(subsr(2*n, S), 2*n+1);
     851       11958 :       shiftr_inplace(S, -2*n);
     852       11958 :       if (realprec(S) != prec) S = rtor(S, prec);
     853             :     }
     854       14850 :     gel(B,n) = gclone(S); /* S = B_2n */
     855             :   }
     856         608 :   if (DEBUGLEVEL) timer_printf(&T, "Bernoulli");
     857         608 :   swap(B, bernzone);
     858         608 :   if (B)
     859             :   { /* kill old non-reused values */
     860       36615 :     for (n = lbern-1; n; n--)
     861             :     {
     862       36160 :       if (gel(B,n) != gel(bernzone,n)) gunclone(gel(B,n));
     863             :     }
     864         455 :     killblock(B);
     865             :   }
     866         608 :   avma = av;
     867             : }
     868             : 
     869             : GEN
     870         588 : bernfrac(long n)
     871             : {
     872             :   long k;
     873         588 :   if (n < 0) pari_err_DOMAIN("bernfrac", "index", "<", gen_0, stoi(n));
     874         581 :   if (n == 0) return gen_1;
     875         497 :   if (n == 1) return mkfrac(gen_m1,gen_2);
     876         420 :   if (odd(n)) return gen_0;
     877         329 :   k = n >> 1;
     878         329 :   if (!bernzone && k <= BERN_MINNB) mpbern(BERN_MINNB, 0);
     879         329 :   if (bernzone && k < lg(bernzone))
     880             :   {
     881         119 :     GEN B = gel(bernzone, k), C;
     882         119 :     if (typ(B) != t_REAL) return B;
     883           0 :     C = bernfrac_using_zeta(n);
     884           0 :     gel(bernzone, k) = gclone(C);
     885           0 :     gunclone(B); return C;
     886             :   }
     887         210 :   return bernfrac_using_zeta(n);
     888             : }
     889             : 
     890             : /* mpbern as exact fractions */
     891             : static GEN
     892          21 : bernvec_old(long nb)
     893             : {
     894             :   long n, i;
     895             :   GEN y;
     896             : 
     897          21 :   if (nb < 0) return cgetg(1, t_VEC);
     898           3 :   if (nb > 46340 && BITS_IN_LONG == 32) pari_err_IMPL( "bernvec for n > 46340");
     899             : 
     900          21 :   y = cgetg(nb+2, t_VEC); gel(y,1) = gen_1;
     901         126 :   for (n = 1; n <= nb; n++)
     902             :   { /* compute y[n+1] = B_{2n} */
     903         105 :     pari_sp av = avma;
     904         105 :     GEN b = gmul2n(utoineg(2*n - 1), -1); /* 1 + (2n+1)B_1 = -(2n-1) /2 */
     905         105 :     GEN c = gen_1;
     906         105 :     ulong u1 = 2*n + 1, u2 = n, d1 = 1, d2 = 1;
     907             : 
     908         336 :     for (i = 1; i < n; i++)
     909             :     {
     910         231 :       c = diviiexact(muliu(c, u1*u2), utoipos(d1*d2));/*= binomial(2n+1, 2*i) */
     911         231 :       b = gadd(b, gmul(c, gel(y,i+1)));
     912         231 :       u1 -= 2; u2--; d1++; d2 += 2;
     913             :     }
     914         105 :     gel(y,n+1) = gerepileupto(av, gdivgs(b, -(1+2*n)));
     915             :   }
     916          21 :   return y;
     917             : }
     918             : GEN
     919          28 : bernvec(long nb)
     920             : {
     921          28 :   long i, l = nb+2;
     922          28 :   GEN y = cgetg(l, t_VEC);
     923          28 :   if (nb < 20) return bernvec_old(nb);
     924           7 :   for (i = 1; i < l; i++) gel(y,i) = bernfrac((i-1) << 1);
     925           7 :   return y;
     926             : }
     927             : 
     928             : /* x := pol_x(v); B_k(x) = \sum_{i=0}^k binomial(k, i) B_i x^{k-i} */
     929             : static GEN
     930          77 : bernpol_i(long k, long v)
     931             : {
     932             :   GEN B, C;
     933             :   long i;
     934          77 :   if (v < 0) v = 0;
     935          77 :   if (k < 0) pari_err_DOMAIN("bernpol", "index", "<", gen_0, stoi(k));
     936          70 :   mpbern(k >> 1, 0); /* cache B_2, ..., B_2[k/2] */
     937          70 :   C = vecbinome(k);
     938          70 :   B = cgetg(k + 3, t_POL);
     939          70 :   for (i = 0; i <= k; ++i) gel(B, k-i+2) = gmul(gel(C,i+1), bernfrac(i));
     940          70 :   B[1] = evalsigne(1) | evalvarn(v);
     941          70 :   return B;
     942             : }
     943             : GEN
     944          49 : bernpol(long k, long v) {
     945          49 :   pari_sp av = avma;
     946          49 :   return gerepileupto(av, bernpol_i(k, v));
     947             : }
     948             : /* x := pol_x(v); return 1^e + ... + x^e = x^e + (B_{e+1}(x) - B_{e+1})/(e+1) */
     949             : static GEN
     950          42 : faulhaber(long e, long v)
     951             : {
     952             :   GEN B;
     953          42 :   if (e == 0) return pol_x(v);
     954          28 :   B = RgX_integ(bernpol_i(e, v)); /* (B_{e+1}(x) - B_{e+1}) / (e+1) */
     955          28 :   gel(B,e+2) = gaddgs(gel(B,e+2), 1); /* add x^e, in place */
     956          28 :   return B;
     957             : }
     958             : /* sum_v T(v), T a polynomial expression in v */
     959             : GEN
     960          49 : sumformal(GEN T, long v)
     961             : {
     962          49 :   pari_sp av = avma, av2;
     963             :   long i, t, d;
     964             :   GEN R;
     965             : 
     966          49 :   T = simplify_shallow(T);
     967          49 :   t = typ(T);
     968          49 :   if (is_scalar_t(t))
     969          14 :     return gerepileupto(av, monomialcopy(T, 1, v < 0? 0: v));
     970          35 :   if (t != t_POL) pari_err_TYPE("sumformal [not a t_POL]", T);
     971          28 :   if (v < 0) v = varn(T);
     972          28 :   av2 = avma;
     973          28 :   R = gen_0;
     974          28 :   d = poldegree(T,v);
     975          91 :   for (i = d; i >= 0; i--)
     976             :   {
     977          63 :     GEN c = polcoeff0(T, i, v);
     978          63 :     if (gequal0(c)) continue;
     979          42 :     R = gadd(R, gmul(c, faulhaber(i, v)));
     980          42 :     if (gc_needed(av2,3))
     981             :     {
     982           0 :       if(DEBUGMEM>1) pari_warn(warnmem,"sumformal, i = %ld/%ld", i,d);
     983           0 :       R = gerepileupto(av2, R);
     984             :     }
     985             :   }
     986          28 :   return gerepileupto(av, R);
     987             : }
     988             : 
     989             : /********************************************************************/
     990             : /**                                                                **/
     991             : /**                         EULER'S GAMMA                          **/
     992             : /**                                                                **/
     993             : /********************************************************************/
     994             : 
     995             : /* x / (i*(i+1)) */
     996             : GEN
     997    28865077 : divrunu(GEN x, ulong i)
     998             : {
     999    28865077 :   if (i <= LOWMASK) /* i(i+1) < 2^BITS_IN_LONG*/
    1000    28865077 :     return divru(x, i*(i+1));
    1001             :   else
    1002           0 :     return divru(divru(x, i), i+1);
    1003             : }
    1004             : /* x / (i*(i+1)) */
    1005             : GEN
    1006      782439 : divgunu(GEN x, ulong i)
    1007             : {
    1008             : #ifdef LONG_IS_64BIT
    1009      673968 :   if (i < 3037000500L) /* i(i+1) < 2^63 */
    1010             : #else
    1011      108471 :   if (i < 46341L) /* i(i+1) < 2^31 */
    1012             : #endif
    1013      782439 :     return gdivgs(x, i*(i+1));
    1014             :   else
    1015           0 :     return gdivgs(gdivgs(x, i), i+1);
    1016             : }
    1017             : 
    1018             : /* arg(s+it) */
    1019             : double
    1020      177349 : darg(double s, double t)
    1021             : {
    1022             :   double x;
    1023      177349 :   if (!t) return (s>0)? 0.: M_PI;
    1024      171367 :   if (!s) return (t>0)? M_PI/2: -M_PI/2;
    1025      171360 :   x = atan(t/s);
    1026      171360 :   return (s>0)? x
    1027      171360 :               : ((t>0)? x+M_PI : x-M_PI);
    1028             : }
    1029             : 
    1030             : void
    1031      177349 : dcxlog(double s, double t, double *a, double *b)
    1032             : {
    1033      177349 :   *a = log(s*s + t*t) / 2; /* log |s| = Re(log(s)) */
    1034      177349 :   *b = darg(s,t);          /* Im(log(s)) */
    1035      177349 : }
    1036             : 
    1037             : double
    1038       12782 : dabs(double s, double t) { return sqrt( s*s + t*t ); }
    1039             : double
    1040          14 : dnorm(double s, double t) { return s*s + t*t; }
    1041             : 
    1042             : #if 0
    1043             : /* x, z t_REAL. Compute unique x in ]-z,z] congruent to x mod 2z */
    1044             : static GEN
    1045             : red_mod_2z(GEN x, GEN z)
    1046             : {
    1047             :   GEN Z = gmul2n(z, 1), d = subrr(z, x);
    1048             :   /* require little accuracy */
    1049             :   if (!signe(d)) return x;
    1050             :   setprec(d, nbits2prec(expo(d) - expo(Z)));
    1051             :   return addrr(mulir(floorr(divrr(d, Z)), Z), x);
    1052             : }
    1053             : #endif
    1054             : 
    1055             : /* lngamma(1+z) = -Euler*z + sum_{i > 1} zeta(i)/i (-z)^i
    1056             :  * at relative precision prec, |z| < 1 is small */
    1057             : static GEN
    1058        2510 : lngamma1(GEN z, long prec)
    1059             : { /* sum_{i > l} |z|^(i-1) = |z|^l / (1-|z|) < 2^-B
    1060             :    * for l > (B+1) / |log2(|z|)| */
    1061        2510 :   long i, l = ceil((bit_accuracy(prec) + 1) / - dbllog2(z));
    1062        2510 :   GEN zet, me = mpeuler(prec), s = gen_0;
    1063        2510 :   setsigne(me, -1); /* -Euler */
    1064        2510 :   if (l <= 1) return gmul(me, z);
    1065        2335 :   zet = veczeta(gen_1, gen_2, l-1, prec); /* z[i] = zeta(i+1) */
    1066       24085 :   for (i = l; i > 1; i--)
    1067             :   {
    1068       21750 :     GEN c = divru(gel(zet,i-1), i);
    1069       21750 :     if (odd(i)) setsigne(c, -1);
    1070       21750 :     s = gadd(gmul(s,z), c);
    1071             :   }
    1072        2335 :   return gmul(z, gadd(gmul(s,z), me));
    1073             : }
    1074             : 
    1075             : static GEN
    1076      176856 : cxgamma(GEN s0, int dolog, long prec)
    1077             : {
    1078             :   GEN s, u, a, y, res, tes, sig, tau, invn2, p1, nnx, pi, pi2, sqrtpi2;
    1079             :   long i, lim, nn, esig, et;
    1080             :   pari_sp av, av2;
    1081      176856 :   int funeq = 0;
    1082             :   pari_timer T;
    1083             : 
    1084      176856 :   if (DEBUGLEVEL>5) timer_start(&T);
    1085      176856 :   s = trans_fix_arg(&prec,&s0,&sig,&tau,&av,&res);
    1086             : 
    1087      176856 :   esig = expo(sig);
    1088      176856 :   et = signe(tau)? expo(tau): 0;
    1089      176856 :   if ((signe(sig) <= 0 || esig < -1) && et <= 16)
    1090             :   { /* s <--> 1-s */
    1091        3220 :     funeq = 1; s = gsubsg(1, s); sig = real_i(s);
    1092             :   }
    1093             : 
    1094             :   /* find "optimal" parameters [lim, nn] */
    1095      176856 :   if (esig > 300 || et > 300)
    1096          35 :   { /* |s| is HUGE ! Play safe and avoid inf / NaN */
    1097             :     GEN S, iS, l2, la, u;
    1098             :     double logla, l;
    1099             : 
    1100          35 :     S = gprec_w(s,LOWDEFAULTPREC);
    1101             :     /* l2 ~ |lngamma(s))|^2 */
    1102          35 :     l2 = gnorm(gmul(S, glog(S, LOWDEFAULTPREC)));
    1103          35 :     l = (prec2nbits_mul(prec, LOG2) - rtodbl(glog(l2,LOWDEFAULTPREC))/2) / 2.;
    1104          35 :     if (l < 0) l = 0.;
    1105             : 
    1106          35 :     iS = imag_i(S);
    1107          35 :     if (et > 0 && l > 0)
    1108          21 :     {
    1109          21 :       GEN t = gmul(iS, dbltor(M_PI / l)), logt = glog(t,LOWDEFAULTPREC);
    1110          21 :       la = gmul(t, logt);
    1111          21 :       if      (gcmpgs(la, 3) < 0)   { logla = log(3.); la = stoi(3); }
    1112          14 :       else if (gcmpgs(la, 150) > 0) { logla = rtodbl(logt); la = t; }
    1113           7 :       else logla = rtodbl(mplog(la));
    1114             :     }
    1115             :     else
    1116             :     {
    1117          14 :       logla = log(3.); la = stoi(3);
    1118             :     }
    1119          35 :     lim = (long)ceil(l / (1.+ logla));
    1120          35 :     if (lim == 0) lim = 1;
    1121             : 
    1122          35 :     u = gmul(la, dbltor((lim-0.5)/M_PI));
    1123          35 :     l2 = gsub(gsqr(u), gsqr(iS));
    1124          35 :     if (signe(l2) > 0)
    1125             :     {
    1126          14 :       l2 = gsub(gsqrt(l2,3), sig);
    1127          14 :       if (signe(l2) > 0) nn = itos( gceil(l2) ); else nn = 1;
    1128             :     }
    1129             :     else
    1130          21 :       nn = 1;
    1131             :   }
    1132             :   else
    1133             :   { /* |s| is moderate. Use floats  */
    1134      176821 :     double ssig = rtodbl(sig);
    1135      176821 :     double st = typ(s) == t_REAL? 0.0: rtodbl(imag_i(s));
    1136             :     double la, l,l2,u,v, rlogs, ilogs;
    1137             : 
    1138      176821 :     if (fabs(ssig-1) + fabs(st) < 0.0001)
    1139             :     { /* s ~ 1: loggamma(1+u) ~ - Euler * u, cancellation */
    1140        2076 :       if (funeq) /* s0 ~ 0: use lngamma(s0)+log(s0) = lngamma(s0+1) */
    1141             :       {
    1142          14 :         if (dolog)
    1143           7 :           y = gsub(lngamma1(s0,prec), glog(s0,prec));
    1144             :         else
    1145           7 :           y = gdiv(gexp(lngamma1(s0,prec), prec), s0);
    1146             :       }
    1147             :       else
    1148             :       {
    1149        4138 :         if (isint1(s0)) { avma = av; return dolog? real_0(prec): real_1(prec); }
    1150        1271 :         y = lngamma1(gsubgs(s0,1),prec);
    1151        1271 :         if (!dolog) y = gexp(y,prec);
    1152             :       }
    1153        1285 :       avma = av; return affc_fixlg(y, res);
    1154             :     }
    1155      174745 :     dcxlog(ssig,st, &rlogs,&ilogs);
    1156             :     /* Re (s - 1/2) log(s) */
    1157      174745 :     u = (ssig - 0.5)*rlogs - st * ilogs;
    1158             :     /* Im (s - 1/2) log(s) */
    1159      174745 :     v = (ssig - 0.5)*ilogs + st * rlogs;
    1160             :     /* l2 = | (s - 1/2) log(s) - s + log(2Pi)/2 |^2 ~ |lngamma(s))|^2 */
    1161      174745 :     u = u - ssig + log(2.*M_PI)/2;
    1162      174745 :     v = v - st;
    1163      174745 :     l2 = u*u + v*v;
    1164      174745 :     if (l2 < 0.000001) l2 = 0.000001;
    1165      174745 :     l = (prec2nbits_mul(prec, LOG2) - log(l2)/2) / 2.;
    1166      174745 :     if (l < 0) l = 0.;
    1167             : 
    1168      174745 :     la = 3.; /* FIXME: heuristic... */
    1169      174745 :     if (st > 1 && l > 0)
    1170             :     {
    1171       63322 :       double t = st * M_PI / l;
    1172       63322 :       la = t * log(t);
    1173       63322 :       if (la < 3) la = 3.;
    1174       63322 :       if (la > 150) la = t;
    1175             :     }
    1176      174745 :     lim = (long)ceil(l / (1.+ log(la)));
    1177      174745 :     if (lim == 0) lim = 1;
    1178             : 
    1179      174745 :     u = (lim-0.5) * la / M_PI;
    1180      174745 :     l2 = u*u - st*st;
    1181      174745 :     if (l2 > 0)
    1182             :     {
    1183      160648 :       nn = (long)ceil(sqrt(l2) - ssig);
    1184      160648 :       if (nn < 1) nn = 1;
    1185             :     }
    1186             :     else
    1187       14097 :       nn = 1;
    1188      174745 :     if (DEBUGLEVEL>5) err_printf("lim, nn: [%ld, %ld], la = %lf\n",lim,nn,la);
    1189             :   }
    1190      174780 :   incrprec(prec);
    1191             : 
    1192      174780 :   av2 = avma;
    1193      174780 :   y = s;
    1194      174780 :   if (typ(s0) == t_INT)
    1195             :   {
    1196         810 :     if (signe(s0) <= 0)
    1197           0 :       pari_err_DOMAIN("gamma","argument", "=",
    1198             :                        strtoGENstr("non-positive integer"), s0);
    1199         810 :     if (is_bigint(s0)) {
    1200           7 :       for (i=1; i < nn; i++)
    1201             :       {
    1202           0 :         y = mulri(y, addis(s0, i));
    1203           0 :         if (gc_needed(av2,3))
    1204             :         {
    1205           0 :           if(DEBUGMEM>1) pari_warn(warnmem,"gamma");
    1206           0 :           y = gerepileuptoleaf(av2, y);
    1207             :         }
    1208             :       }
    1209             :     } else {
    1210         803 :       ulong ss = itou(s0);
    1211        7313 :       for (i=1; i < nn; i++)
    1212             :       {
    1213        6510 :         y = mulru(y, ss + i);
    1214        6510 :         if (gc_needed(av2,3))
    1215             :         {
    1216           0 :           if(DEBUGMEM>1) pari_warn(warnmem,"gamma");
    1217           0 :           y = gerepileuptoleaf(av2, y);
    1218             :         }
    1219             :       }
    1220             :     }
    1221         810 :     if (dolog) y = logr_abs(y);
    1222             :   }
    1223      173970 :   else if (!dolog || typ(s) == t_REAL)
    1224             :   { /* Compute lngamma mod 2 I Pi */
    1225     6411558 :     for (i=1; i < nn; i++)
    1226             :     {
    1227     6237637 :       y = gmul(y, gaddgs(s,i));
    1228     6237637 :       if (gc_needed(av2,3))
    1229             :       {
    1230           0 :         if(DEBUGMEM>1) pari_warn(warnmem,"gamma");
    1231           0 :         y = gerepileupto(av2, y);
    1232             :       }
    1233             :     }
    1234      173921 :     if (dolog) y = logr_abs(y);
    1235             :   }
    1236             :   else
    1237             :   { /* dolog && complex s: be careful with imaginary part */
    1238          49 :     y = glog(y, prec);
    1239         175 :     for (i=1; i < nn; i++)
    1240             :     {
    1241         126 :       y = gadd(y, glog(gaddgs(s,i), prec));
    1242         126 :       if (gc_needed(av2,3))
    1243             :       {
    1244           0 :         if(DEBUGMEM>1) pari_warn(warnmem,"gamma");
    1245           0 :         y = gerepileupto(av2, y);
    1246             :       }
    1247             :     }
    1248             :   }
    1249      174780 :   if (DEBUGLEVEL>5) timer_printf(&T,"product from 0 to N-1");
    1250             : 
    1251      174780 :   nnx = gaddgs(s, nn);
    1252      174780 :   a = ginv(nnx); invn2 = gsqr(a);
    1253      174780 :   av2 = avma;
    1254      174780 :   mpbern(lim,prec);
    1255      174780 :   tes = divrunu(bernreal(2*lim,prec), 2*lim-1); /* B2l / (2l-1) 2l*/
    1256      174780 :   if (DEBUGLEVEL>5) timer_printf(&T,"Bernoullis");
    1257     7823867 :   for (i = 2*lim-2; i > 1; i -= 2)
    1258             :   {
    1259     7649087 :     u = divrunu(bernreal(i,prec), i-1); /* Bi / i(i-1) */
    1260     7649087 :     tes = gadd(u, gmul(invn2,tes));
    1261     7649087 :     if (gc_needed(av2,3))
    1262             :     {
    1263           0 :       if(DEBUGMEM>1) pari_warn(warnmem,"gamma");
    1264           0 :       tes = gerepileupto(av2, tes);
    1265             :     }
    1266             :   }
    1267      174780 :   if (DEBUGLEVEL>5) timer_printf(&T,"Bernoulli sum");
    1268             : 
    1269      174780 :   p1 = gsub(gmul(gsub(nnx, ghalf), glog(nnx,prec)), nnx);
    1270      174780 :   p1 = gadd(p1, gmul(tes, a));
    1271             : 
    1272      174780 :   pi = mppi(prec); pi2 = shiftr(pi, 1); sqrtpi2 = sqrtr(pi2);
    1273             : 
    1274      174780 :   if (dolog)
    1275             :   {
    1276        3658 :     if (funeq)
    1277             :     { /* (recall that s = 1 - s0) */
    1278             : 
    1279             :       /* We compute log(sin(Pi s0)) so that it has branch cuts along
    1280             :       * (-oo, 0] and [1, oo).  To do this in a numerically stable way
    1281             :       * we must compute the log first then mangle its imaginary part.
    1282             :       * The rounding operation below is stable because we're rounding
    1283             :       * a number which is already within 1/4 of an integer. */
    1284             : 
    1285             :       /* z = log( sin(Pi s0) / (sqrt(2Pi)/2) ) */
    1286           7 :       GEN z = glog(gdiv(gsin(gmul(pi,s0),prec), shiftr(sqrtpi2,-1)), prec);
    1287             :       /* b = (2 Re(s) - 1) / 4 */
    1288           7 :       GEN b = shiftr(subrs(shiftr(sig, 1), 1), -2);
    1289           7 :       y = gsub(y, z);
    1290           7 :       if (gsigne(imag_i(s)) > 0) togglesign(b);
    1291             :       /* z = 2Pi round( Im(z)/2Pi - b ) */
    1292           7 :       z = gmul(roundr(gsub(gdiv(imag_i(z), pi2), b)), pi2);
    1293           7 :       if (signe(z)) { /* y += I*z */
    1294           7 :         if (typ(y) == t_COMPLEX)
    1295           7 :           gel(y,2) = gadd(gel(y,2), z);
    1296             :         else
    1297           0 :           y = gadd(y, mkcomplex(gen_0, z));
    1298             :       }
    1299           7 :       p1 = gneg(p1);
    1300             :     }
    1301             :     else /* y --> sqrt(2Pi) / y */
    1302        3651 :       y = gsub(logr_abs(sqrtpi2), y);
    1303        3658 :     y = gadd(p1, y);
    1304             :   }
    1305             :   else
    1306             :   {
    1307      171122 :     if (funeq)
    1308             :     { /* y --> y Pi/(sin(Pi s) * sqrt(2Pi)) = y sqrt(Pi/2)/sin(Pi s) */
    1309        3199 :       y = gdiv(gmul(shiftr(sqrtpi2,-1),y), gsin(gmul(pi,s0), prec));
    1310             :       /* don't use s above: sin(pi s0) = sin(pi s) and the former is
    1311             :        * more accurate, esp. if s0 ~ 0 */
    1312        3199 :       p1 = gneg(p1);
    1313             :     }
    1314             :     else /* y --> sqrt(2Pi) / y */
    1315      167923 :       y = gdiv(sqrtpi2, y);
    1316      171122 :     y = gmul(gexp(p1, prec), y);
    1317             :   }
    1318      174780 :   avma = av; return affc_fixlg(y, res);
    1319             : }
    1320             : 
    1321             : /* Gamma((m+1) / 2) */
    1322             : static GEN
    1323        5782 : gammahs(long m, long prec)
    1324             : {
    1325        5782 :   GEN y = cgetr(prec), z;
    1326        5782 :   pari_sp av = avma;
    1327        5782 :   long ma = labs(m);
    1328             : 
    1329        5782 :   if (ma > 200 + 50*(prec-2)) /* heuristic */
    1330             :   {
    1331           7 :     z = stor(m + 1, prec); shiftr_inplace(z, -1);
    1332           7 :     affrr(cxgamma(z,0,prec), y);
    1333           7 :     avma = av; return y;
    1334             :   }
    1335        5775 :   z = sqrtr( mppi(prec) );
    1336        5775 :   if (m)
    1337             :   {
    1338        2583 :     GEN p1 = mulu_interval(ma/2 + 1, ma);
    1339        2583 :     long v = vali(p1);
    1340        2583 :     p1 = shifti(p1, -v); v -= ma;
    1341        2583 :     if (m >= 0) z = mulri(z,p1);
    1342             :     else
    1343             :     {
    1344          63 :       z = divri(z,p1); v = -v;
    1345          63 :       if ((m&3) == 2) setsigne(z,-1);
    1346             :     }
    1347        2583 :     shiftr_inplace(z, v);
    1348             :   }
    1349        5775 :   affrr(z, y); avma = av; return y;
    1350             : }
    1351             : GEN
    1352          21 : ggammah(GEN x, long prec)
    1353             : {
    1354          21 :   switch(typ(x))
    1355             :   {
    1356             :     case t_INT:
    1357             :     {
    1358          14 :       long k = itos(x);
    1359          14 :       if (labs(k) > 962353) pari_err_OVERFLOW("gammah");
    1360          14 :       return gammahs(k<<1, prec);
    1361             :     }
    1362             :     case t_REAL: case t_COMPLEX: case t_PADIC: case t_SER: {
    1363           7 :       pari_sp av = avma;
    1364           7 :       return gerepileupto(av, ggamma(gadd(x,ghalf), prec));
    1365             :     }
    1366             :   }
    1367           0 :   return trans_eval("gammah",ggammah,x,prec);
    1368             : }
    1369             : 
    1370             : /* find n such that n+v_p(n!)>=k p^2/(p-1)^2 */
    1371             : static long
    1372          28 : nboft(long k, long p)
    1373             : {
    1374          28 :   pari_sp av = avma;
    1375             :   long s, n;
    1376             : 
    1377          28 :   if (k <= 0) return 0;
    1378          28 :   k = itou( gceil(gdiv(mului(k, sqru(p)), sqru(p-1))) );
    1379          28 :   avma = av;
    1380          28 :   for (s=0, n=0; n+s < k; n++, s += u_lval(n, p));
    1381          28 :   return n;
    1382             : }
    1383             : 
    1384             : /* Using Dwork's expansion, compute \Gamma(px+1)=-\Gamma(px) with x a unit.
    1385             :  * See p-Adic Gamma Functions and Dwork Cohomology, Maurizio Boyarsky
    1386             :  * Transactions of the AMS, Vol. 257, No. 2. (Feb., 1980), pp. 359-369.
    1387             :  * Inspired by a GP script by Fernando Rodriguez-Villegas */
    1388             : static GEN
    1389          28 : gadw(GEN x, long p)
    1390             : {
    1391          28 :   pari_sp ltop = avma;
    1392          28 :   GEN s, t, u = cgetg(p+1, t_VEC);
    1393          28 :   long j, k, kp, n = nboft(precp(x)+valp(x)+1, p);
    1394             : 
    1395          28 :   t = s = gaddsg(1, zeropadic(gel(x,2), n));
    1396          28 :   gel(u, 1) = s;
    1397          28 :   gel(u, 2) = s;
    1398          42 :   for (j = 2; j < p; ++j)
    1399          14 :     gel(u, j+1) = gdivgs(gel(u, j), j);
    1400         378 :   for (k = 1, kp = p; k < n; ++k, kp += p) /* kp = k*p */
    1401             :   {
    1402             :     GEN c;
    1403         350 :     gel(u, 1) = gdivgs(gadd(gel(u, 1), gel(u, p)), kp);
    1404         812 :     for (j = 1; j < p; ++j)
    1405         462 :       gel(u, j+1) = gdivgs(gadd(gel(u, j), gel(u, j+1)), kp + j);
    1406             : 
    1407         350 :     t = gmul(t, gaddgs(x, k-1));
    1408         350 :     c = leafcopy(gel(u,1)); setvalp(c, valp(c) + k); /* c = u[1] * p^k */
    1409         350 :     s = gadd(s, gmul(c, t));
    1410         350 :     if ((k&0xFL)==0) gerepileall(ltop, 3, &u,&s,&t);
    1411             :   }
    1412          28 :   return gneg(s);
    1413             : }
    1414             : 
    1415             : /*Use Dwork expansion*/
    1416             : /*This is a O(p*e*log(pe)) algorithm, should be used when p small
    1417             :  * If p==2 this is a O(pe) algorithm. */
    1418             : static GEN
    1419          28 : Qp_gamma_Dwork(GEN x, long p)
    1420             : {
    1421          28 :   pari_sp ltop = avma;
    1422          28 :   long k = padic_to_Fl(x, p);
    1423             :   GEN p1;
    1424             :   long j;
    1425          28 :   long px = precp(x);
    1426          28 :   if (p==2 && px)
    1427             :   {
    1428          14 :     x = shallowcopy(x);
    1429          14 :     setprecp(x, px+1);
    1430          14 :     gel(x,3) = shifti(gel(x,3),1);
    1431             :   }
    1432          28 :   if (k)
    1433             :   {
    1434          21 :     GEN x_k = gsubgs(x,k);
    1435          21 :     x = gdivgs(x_k, p);
    1436          21 :     p1 = gadw(x, p); if (!odd(k)) p1 = gneg(p1);
    1437          21 :     for (j = 1; j < k; ++j) p1 = gmul(p1, gaddgs(x_k, j));
    1438             :   }
    1439             :   else
    1440           7 :     p1 = gneg(gadw(gdivgs(x, p), p));
    1441          28 :   return gerepileupto(ltop, p1);
    1442             : }
    1443             : 
    1444             : /* Compute Qp_gamma using the definition. This is a O(x*M(log(pe))) algorithm.
    1445             :  * This should be used if x is very small. */
    1446             : static GEN
    1447          49 : Qp_gamma_Morita(long n, GEN p, long e)
    1448             : {
    1449          49 :   pari_sp ltop=avma;
    1450          49 :   GEN p2 = gaddsg((n&1)?-1:1, zeropadic(p, e));
    1451             :   long i;
    1452          49 :   long pp=is_bigint(p)? 0: itos(p);
    1453         154 :   for (i = 2; i < n; i++)
    1454         105 :     if (!pp || i%pp)
    1455             :     {
    1456          63 :       p2 = gmulgs(p2, i);
    1457          63 :       if ((i&0xFL) == 0xFL)
    1458           0 :         p2 = gerepileupto(ltop, p2);
    1459             :     }
    1460          49 :   return gerepileupto(ltop, p2);
    1461             : }
    1462             : 
    1463             : /* x\in\N: Gamma(-x)=(-1)^(1+x+x\p)*Gamma(1+x) */
    1464             : static GEN
    1465          28 : Qp_gamma_neg_Morita(long n, GEN p, long e)
    1466             : {
    1467          28 :   GEN g = ginv(Qp_gamma_Morita(n+1, p, e));
    1468          28 :   return ((n^sdivsi(n,p)) & 1)? g: gneg(g);
    1469             : }
    1470             : 
    1471             : /* p-adic Gamma function for x a p-adic integer */
    1472             : /* If n < p*e : use Morita's definition.
    1473             :  * Else : use Dwork's expansion.
    1474             :  * If both n and p are big : itos(p) will fail.
    1475             :  * TODO: handle p=2 better (Qp_gamma_Dwork is slow for p=2). */
    1476             : GEN
    1477          77 : Qp_gamma(GEN x)
    1478             : {
    1479          77 :   GEN n, m, N, p = gel(x,2);
    1480          77 :   long s, e = precp(x);
    1481          77 :   if (absequaliu(p, 2) && e == 2) e = 1;
    1482          77 :   if (valp(x) < 0) pari_err_DOMAIN("gamma","v_p(x)", "<", gen_0, x);
    1483          77 :   n = gtrunc(x);
    1484          77 :   m = gtrunc(gneg(x));
    1485          77 :   N = cmpii(n,m)<=0?n:m;
    1486          77 :   s = itos_or_0(N);
    1487          77 :   if (s && cmpsi(s, muliu(p,e)) < 0) /* s < p*e */
    1488          49 :     return (N == n) ? Qp_gamma_Morita(s,p,e): Qp_gamma_neg_Morita(s,p,e);
    1489          28 :   return Qp_gamma_Dwork(x, itos(p));
    1490             : }
    1491             : 
    1492             : /* gamma(1+x) - 1, |x| < 1 is "small" */
    1493             : GEN
    1494        1211 : ggamma1m1(GEN x, long prec) { return gexpm1(lngamma1(x, prec), prec); }
    1495             : 
    1496             : /* lngamma(y) with 0 constant term, using (lngamma y)' = y' psi(y) */
    1497             : static GEN
    1498        9709 : serlngamma0(GEN y, long prec)
    1499             : {
    1500             :   GEN t;
    1501        9709 :   if (valp(y)) pari_err_DOMAIN("lngamma","valuation", "!=", gen_0, y);
    1502        9702 :   t = derivser(y);
    1503             :   /* don't compute psi if y'=0 */
    1504        9702 :   if (signe(t)) t = gmul(t, gpsi(y,prec));
    1505        9702 :   return integser(t);
    1506             : }
    1507             : 
    1508             : GEN
    1509      191317 : ggamma(GEN x, long prec)
    1510             : {
    1511             :   pari_sp av;
    1512             :   GEN y, z;
    1513             : 
    1514      191317 :   switch(typ(x))
    1515             :   {
    1516             :     case t_INT:
    1517        3913 :       if (signe(x) <= 0)
    1518           0 :         pari_err_DOMAIN("gamma","argument", "=",
    1519             :                          strtoGENstr("non-positive integer"), x);
    1520        3913 :       if (abscmpiu(x,481177) > 0) pari_err_OVERFLOW("gamma");
    1521        3913 :       return mpfactr(itos(x) - 1, prec);
    1522             : 
    1523             :     case t_REAL: case t_COMPLEX:
    1524      171388 :       return cxgamma(x, 0, prec);
    1525             : 
    1526             :     case t_FRAC:
    1527             :     {
    1528        6272 :       GEN a = gel(x,1), b = gel(x,2), c;
    1529             :       long m;
    1530        6272 :       if (absequaliu(b,2))
    1531             :       {
    1532        5768 :         if (is_bigint(a) || labs(m = itos(a)) > 962354)
    1533             :         {
    1534           0 :           pari_err_OVERFLOW("gamma");
    1535           0 :           return NULL; /* not reached */
    1536             :         }
    1537        5768 :         return gammahs(m-1, prec);
    1538             :       }
    1539         504 :       av = avma; c = subii(a,b);
    1540         504 :       if (expi(c) - expi(b) < -10)
    1541             :       {
    1542           7 :         x = mkfrac(c,b);
    1543           7 :         if (lgefint(b) >= prec) x = fractor(x,prec);
    1544           7 :         y = mpexp(lngamma1(x, prec));
    1545             :       }
    1546             :       else
    1547             :       {
    1548         497 :         x = fractor(x, prec);
    1549         497 :         y = cxgamma(x, 0, prec);
    1550             :       }
    1551         504 :       return gerepileupto(av, y);
    1552             :     }
    1553             : 
    1554          70 :     case t_PADIC: return Qp_gamma(x);
    1555             :     default:
    1556        9674 :       av = avma; if (!(y = toser_i(x))) break;
    1557        9674 :       if (lg(y) == 2) pari_err_DOMAIN("gamma", "argument", "=", gen_0,y);
    1558             :       /* exp(lngamma) */
    1559        9667 :       if (valp(y) > 0 || lg(y) == 2)
    1560         112 :         z = gdiv(gexp(glngamma(gaddgs(y,1),prec),prec),y);
    1561             :       else
    1562             :       {
    1563        9555 :         GEN Y = y, y0 = simplify_shallow(gel(y,2));
    1564        9555 :         z = NULL;
    1565        9555 :         if (isint(y0, &y0))
    1566             :         { /* fun eq. avoids log singularity of lngamma at negative ints */
    1567        5152 :           if (signe(y0) < 0) { Y = gsubsg(1, y); y0 = subsi(1, y0); }
    1568        5152 :           if (abscmpiu(y0, 50) < 0) z = mpfact(itos(y0)-1); /* more precise */
    1569             :         }
    1570        9555 :         if (!z) z = ggamma(y0,prec);
    1571        9555 :         z = gmul(z, gexp(serlngamma0(Y,prec),prec));
    1572        9555 :         if (Y != y)
    1573             :         {
    1574           7 :           GEN pi = mppi(prec);
    1575           7 :           z = gdiv(mpodd(y0)? pi: negr(pi),
    1576             :                    gmul(z, gsin(gmul(pi,serchop0(y)), prec)));
    1577             :         }
    1578             :       }
    1579        9667 :       return gerepileupto(av, z);
    1580             :   }
    1581           0 :   return trans_eval("gamma",ggamma,x,prec);
    1582             : }
    1583             : 
    1584             : GEN
    1585       12314 : mpfactr(long n, long prec)
    1586             : {
    1587       12314 :   GEN f = cgetr(prec);
    1588       12314 :   pari_sp av = avma;
    1589             : 
    1590       12314 :   if (n+1 > 350 + 70*(prec-2)) /* heuristic */
    1591          21 :     affrr(cxgamma(stor(n+1, prec), 0, prec), f);
    1592             :   else
    1593       12293 :     affir(mpfact(n), f);
    1594       12314 :   avma = av; return f;
    1595             : }
    1596             : 
    1597             : GEN
    1598        5253 : glngamma(GEN x, long prec)
    1599             : {
    1600        5253 :   pari_sp av = avma;
    1601             :   GEN y, y0, t;
    1602             : 
    1603        5253 :   switch(typ(x))
    1604             :   {
    1605             :     case t_INT:
    1606         154 :       if (signe(x) <= 0)
    1607           0 :         pari_err_DOMAIN("lngamma","argument", "=",
    1608             :                          strtoGENstr("non-positive integer"), x);
    1609         154 :       if (abscmpiu(x,200 + 50*(prec-2)) > 0) /* heuristic */
    1610          19 :         return cxgamma(x, 1, prec);
    1611         135 :       return gerepileuptoleaf(av, logr_abs( itor(mpfact(itos(x) - 1), prec) ));
    1612             :     case t_FRAC:
    1613             :     {
    1614         602 :       GEN a = gel(x,1), b = gel(x,2), c = subii(a,b);
    1615         602 :       long e = expi(b) - expi(c);
    1616         602 :       if (e > 10)
    1617             :       {
    1618           7 :         x = mkfrac(c,b);
    1619           7 :         if (lgefint(b) >= prec) x = fractor(x,prec + nbits2nlong(e));
    1620           7 :         y = lngamma1(x, prec);
    1621             :       }
    1622             :       else
    1623             :       {
    1624         595 :         x = fractor(x, e > 1? prec+EXTRAPRECWORD: prec);
    1625         595 :         y = cxgamma(x, 1, prec);
    1626             :       }
    1627         602 :       return gerepileupto(av, y);
    1628             :     }
    1629             : 
    1630             :     case t_REAL: case t_COMPLEX:
    1631        4329 :       return cxgamma(x, 1, prec);
    1632             : 
    1633             :     default:
    1634         161 :       if (!(y = toser_i(x))) break;
    1635         161 :       if (lg(y) == 2) pari_err_DOMAIN("lngamma", "argument", "=", gen_0,y);
    1636         154 :       t = serlngamma0(y,prec);
    1637         140 :       y0 = simplify_shallow(gel(y,2));
    1638             :       /* no constant term if y0 = 1 or 2 */
    1639         140 :       if (!isint(y0,&y0) || signe(y0) <= 0 || abscmpiu(y0,2) > 2)
    1640           7 :         t = gadd(t, glngamma(y0,prec));
    1641         140 :       return gerepileupto(av, t);
    1642             : 
    1643           7 :     case t_PADIC: return gerepileupto(av, Qp_log(Qp_gamma(x)));
    1644             :   }
    1645           0 :   return trans_eval("lngamma",glngamma,x,prec);
    1646             : }
    1647             : /********************************************************************/
    1648             : /**                                                                **/
    1649             : /**                  PSI(x) = GAMMA'(x)/GAMMA(x)                   **/
    1650             : /**                                                                **/
    1651             : /********************************************************************/
    1652             : 
    1653             : GEN
    1654          21 : cxpsi(GEN s0, long prec)
    1655             : {
    1656             :   pari_sp av, av2;
    1657             :   GEN sum,z,a,res,tes,in2,sig,tau,s,unr;
    1658             :   long lim,nn,k;
    1659          21 :   const long la = 3;
    1660          21 :   int funeq = 0;
    1661             :   pari_timer T;
    1662             : 
    1663          21 :   if (DEBUGLEVEL>2) timer_start(&T);
    1664          21 :   s = trans_fix_arg(&prec,&s0,&sig,&tau,&av,&res);
    1665          21 :   if (signe(sig) <= 0) { funeq = 1; s = gsub(gen_1, s); sig = real_i(s); }
    1666          21 :   if (typ(s0) == t_INT && signe(s0) <= 0)
    1667           0 :     pari_err_DOMAIN("psi","argument", "=",
    1668             :                     strtoGENstr("non-positive integer"), s0);
    1669             : 
    1670          21 :   if (expo(sig) > 300 || (typ(s) == t_COMPLEX && gexpo(gel(s,2)) > 300))
    1671           7 :   { /* |s| is HUGE. Play safe */
    1672           7 :     GEN L, S = gprec_w(s,LOWDEFAULTPREC), rS = real_i(S), iS = imag_i(S);
    1673             :     double l;
    1674             : 
    1675           7 :     l = rtodbl( gnorm(glog(S, 3)) );
    1676           7 :     l = log(l) / 2.;
    1677           7 :     lim = 2 + (long)ceil((prec2nbits_mul(prec, LOG2) - l) / (2*(1+log((double)la))));
    1678           7 :     if (lim < 2) lim = 2;
    1679             : 
    1680           7 :     l = (2*lim-1)*la / (2.*M_PI);
    1681           7 :     L = gsub(dbltor(l*l), gsqr(iS));
    1682           7 :     if (signe(L) < 0) L = gen_0;
    1683             : 
    1684           7 :     L = gsub(gsqrt(L, 3), rS);
    1685           7 :     if (signe(L) > 0) nn = (long)ceil(rtodbl(L)); else nn = 1;
    1686           7 :     if (DEBUGLEVEL>2) err_printf("lim, nn: [%ld, %ld]\n",lim,nn);
    1687             :   }
    1688             :   else
    1689             :   {
    1690          14 :     double ssig = rtodbl(sig);
    1691          14 :     double st = typ(s) == t_REAL? 0.0: rtodbl(imag_i(s));
    1692             :     double l;
    1693             :     {
    1694             :       double rlog, ilog; /* log (s - Euler) */
    1695          14 :       dcxlog(ssig - 0.57721566, st, &rlog,&ilog);
    1696          14 :       l = dnorm(rlog,ilog);
    1697             :     }
    1698          14 :     if (l < 0.000001) l = 0.000001;
    1699          14 :     l = log(l) / 2.;
    1700          14 :     lim = 2 + (long)ceil((prec2nbits_mul(prec, LOG2) - l) / (2*(1+log((double)la))));
    1701          14 :     if (lim < 2) lim = 2;
    1702             : 
    1703          14 :     l = (2*lim-1)*la / (2.*M_PI);
    1704          14 :     l = l*l - st*st;
    1705          14 :     if (l < 0.) l = 0.;
    1706          14 :     nn = (long)ceil( sqrt(l) - ssig );
    1707          14 :     if (nn < 1) nn = 1;
    1708          14 :     if (DEBUGLEVEL>2) err_printf("lim, nn: [%ld, %ld]\n",lim,nn);
    1709             :   }
    1710          21 :   incrprec(prec); unr = real_1(prec); /* one extra word of precision */
    1711             : 
    1712          21 :   a = gdiv(unr, gaddgs(s, nn)); /* 1 / (s+n) */
    1713          21 :   av2 = avma; sum = gmul2n(a,-1);
    1714         959 :   for (k = 0; k < nn; k++)
    1715             :   {
    1716         938 :     sum = gadd(sum, gdiv(unr, gaddgs(s, k)));
    1717         938 :     if ((k & 127) == 0) sum = gerepileupto(av2, sum);
    1718             :   }
    1719          21 :   z = gsub(glog(gaddgs(s, nn), prec), sum);
    1720          21 :   if (DEBUGLEVEL>2) timer_printf(&T,"sum from 0 to N-1");
    1721             : 
    1722          21 :   in2 = gsqr(a);
    1723          21 :   mpbern(lim,prec);
    1724          21 :   av2 = avma; tes = divru(bernreal(2*lim, prec), 2*lim);
    1725        1155 :   for (k=2*lim-2; k>=2; k-=2)
    1726             :   {
    1727        1134 :     tes = gadd(gmul(in2,tes), divru(bernreal(k, prec), k));
    1728        1134 :     if ((k & 255) == 0) tes = gerepileupto(av2, tes);
    1729             :   }
    1730          21 :   if (DEBUGLEVEL>2) timer_printf(&T,"Bernoulli sum");
    1731          21 :   z = gsub(z, gmul(in2,tes));
    1732          21 :   if (funeq)
    1733             :   {
    1734           7 :     GEN pi = mppi(prec);
    1735           7 :     z = gadd(z, gmul(pi, gcotan(gmul(pi,s), prec)));
    1736             :   }
    1737          21 :   avma = av; return affc_fixlg(z, res);
    1738             : }
    1739             : 
    1740             : /* n > 0; return psi(1+x) + O(x^n), x = pol_x(v) */
    1741             : static GEN
    1742        4718 : serpsi1(long n, long v, long prec)
    1743             : {
    1744        4718 :   long i, l = n+3;
    1745        4718 :   GEN z, g, s = cgetg(l, t_SER);
    1746        4718 :   s[1] = evalsigne(1)|evalvalp(0)|evalvarn(v);
    1747        4718 :   g = mpeuler(prec); setsigne(g, -1);
    1748        4718 :   z = veczeta(gen_1, gen_2, n, prec); /* zeta(2..n) */
    1749        4718 :   gel(s,2) = g;
    1750       20356 :   for (i = 2; i < l-1; i++)
    1751             :   {
    1752       15638 :     GEN c = gel(z,i-1); /* zeta(i) */
    1753       15638 :     if (odd(i)) setsigne(c, -1);
    1754       15638 :     gel(s,i+1) = c;
    1755             :   }
    1756        4718 :   return s;
    1757             : }
    1758             : /* T an RgX, return T(X + z0) + O(X^L) */
    1759             : static GEN
    1760      693076 : tr(GEN T, GEN z0, long L)
    1761             : {
    1762      693076 :   GEN s = RgX_to_ser(RgX_translate(T, z0), L+3);
    1763      693076 :   setvarn(s, 0); return s;
    1764             : }
    1765             : /* z0 a complex number with Re(z0) > 1/2; return psi(z0+x) + O(x^L)
    1766             :  * using Luke's rational approximation for psi(x) */
    1767             : static GEN
    1768        3178 : serpsiz0(GEN z0, long L, long v, long prec)
    1769             : {
    1770             :   pari_sp av;
    1771             :   GEN A,A1,A2, B,B1,B2, Q;
    1772             :   long n;
    1773        3178 :   n = gprecision(z0); if (n) prec = n;
    1774        3178 :   z0 = gtofp(z0, prec + EXTRAPRECWORD);
    1775             :   /* Start from n = 3; in Luke's notation, A2 := A_{n-2}, A1 := A_{n-1},
    1776             :    * A := A_n. Same for B */
    1777        3178 :   av = avma;
    1778        3178 :   A2= gdivgs(mkpoln(2, gen_1, utoipos(6)), 2);
    1779        3178 :   B2 = scalarpol_shallow(utoipos(4), 0);
    1780        3178 :   A1= gdivgs(mkpoln(3, gen_1, utoipos(82), utoipos(96)), 6);
    1781        3178 :   B1 = mkpoln(2, utoipos(8), utoipos(28));
    1782        3178 :   A = gdivgs(mkpoln(4, gen_1, utoipos(387), utoipos(2906), utoipos(1920)), 12);
    1783        3178 :   B = mkpoln(3, utoipos(14), utoipos(204), utoipos(310));
    1784        3178 :   A2= tr(A2,z0, L);
    1785        3178 :   B2= tr(B2,z0, L);
    1786        3178 :   A1= tr(A1,z0, L);
    1787        3178 :   B1= tr(B1,z0, L);
    1788        3178 :   A = tr(A, z0, L);
    1789        3178 :   B = tr(B, z0, L); Q = gdiv(A, B);
    1790             :   /* work with z0+x as a variable */
    1791      223610 :   for (n = 4;; n++)
    1792             :   {
    1793      223610 :     GEN Q0 = Q, a, b, r, c3,c2,c1,c0 = muluu(2*n-3, n+1);
    1794      223610 :     GEN u = subiu(muluu(n, 7*n-9), 6);
    1795      223610 :     GEN t = addiu(muluu(n, 7*n-19), 4);
    1796             :     /* c1=(2*n-1)*(3*(n-1)*z+7*n^2-9*n-6);
    1797             :      * c2=(2*n-3)*(z-n-1)*(-3*(n-1)*z+7*n^2-19*n+4);
    1798             :      * c3=(2*n-1)*(n-3)*(z-n)*(z-(n+1))*(z+(n-4)); */
    1799      223610 :     c1 = deg1pol_shallow(muluu(3*(n-1),2*n-1), muliu(u,2*n-1), 0);
    1800      223610 :     c2 = ZX_mul(deg1pol_shallow(utoipos(2*n-3), negi(muluu(2*n-3,n+1)), 0),
    1801      223610 :                 deg1pol_shallow(utoineg(3*(n-1)), t, 0));
    1802      223610 :     r = mkvec3(utoipos(n), utoipos(n+1), stoi(4-n));
    1803      223610 :     c3 = ZX_Z_mul(roots_to_pol(r,0), muluu(2*n-1,n-3));
    1804      223610 :     c1 = tr(c1, z0, L+3);
    1805      223610 :     c2 = tr(c2, z0, L+3);
    1806      223610 :     c3 = tr(c3, z0, L+3);
    1807             : 
    1808             :     /* A_{n+1}, B_{n+1} */
    1809      223610 :     a = gdiv(gadd(gadd(gmul(c1,A),gmul(c2,A1)),gmul(c3,A2)), c0);
    1810      223610 :     b = gdiv(gadd(gadd(gmul(c1,B),gmul(c2,B1)),gmul(c3,B2)), c0);
    1811      223610 :     Q = gdiv(a,b);
    1812      223610 :     if (gexpo(gsub(Q,Q0)) < -prec2nbits(prec)) break;
    1813      220432 :     A2 = A1; A1 = A; A = a;
    1814      220432 :     B2 = B1; B1 = B; B = b;
    1815      220432 :     if (gc_needed(av,1))
    1816             :     {
    1817           0 :       if(DEBUGMEM>1) pari_warn(warnmem,"serpsiz0, n = %ld", n);
    1818           0 :       gerepileall(av, 7, &A,&A1,&A2, &B,&B1,&B2, &Q);
    1819             :     }
    1820      220432 :   }
    1821        3178 :   Q = gmul(Q, gmul2n(gsubsg(1, ginv(tr(pol_x(v),z0, L))), 1));
    1822        3178 :   setvarn(Q, v);
    1823        3178 :   return gadd(negr(mpeuler(prec)), Q);
    1824             : }
    1825             : /* sum (-1)^k*H(m,k)x^k + O(x^L); L > 0;
    1826             :  * H(m,k) = (-1)^{k * \delta_{m > 0}} sum_{1<=i<m} 1/i^(k+1) */
    1827             : static GEN
    1828         301 : Hseries(long m, long L, long v, long prec)
    1829             : {
    1830         301 :   long i, k, bit, l = L+3, M = m < 0? 1-m: m;
    1831         301 :   pari_sp av = avma;
    1832         301 :   GEN H = cgetg(l, t_SER);
    1833         301 :   H[1] = evalsigne(1)|evalvarn(v)|evalvalp(0);
    1834         301 :   prec++;
    1835         301 :   bit = -prec2nbits(prec);
    1836         301 :   for(k = 2; k < l; k++) gel(H,k) = gen_1; /* i=1 */
    1837         322 :   for (i = 2; i < M; i++)
    1838             :   {
    1839          21 :     GEN ik = invr(utor(i, prec));
    1840         154 :     for (k = 2; k < l; k++)
    1841             :     {
    1842         133 :       if (k > 2) { ik = divru(ik, i); if (expo(ik) < bit) break; }
    1843         133 :       gel(H,k) = gadd(gel(H,k), ik);
    1844             :     }
    1845          21 :     if (gc_needed(av,3))
    1846             :     {
    1847           0 :       if(DEBUGMEM>1) pari_warn(warnmem,"Hseries, i = %ld/%ld", i,M);
    1848           0 :       H = gerepilecopy(av, H);
    1849             :     }
    1850             :   }
    1851         301 :   if (m > 0)
    1852         287 :     for (k = 3; k < l; k+=2) togglesign_safe(&gel(H,k));
    1853         301 :   return H;
    1854             : }
    1855             : 
    1856             : static GEN
    1857        7910 : serpsi(GEN y, long prec)
    1858             : {
    1859        7910 :   GEN Q = NULL, z0, Y = y, Y2;
    1860        7910 :   long L = lg(y)-2, v  = varn(y), vy = valp(y);
    1861             : 
    1862        7910 :   if (!L) pari_err_DOMAIN("psi", "argument", "=", gen_0,y);
    1863        7903 :   if (vy < 0) pari_err_DOMAIN("psi", "series valuation", "<", gen_0,y);
    1864        7903 :   if (vy)
    1865          14 :     z0 = gen_0;
    1866             :   else
    1867             :   {
    1868        7889 :     z0 = simplify_shallow(gel(y,2));
    1869        7889 :     (void)isint(z0, &z0);
    1870             :   }
    1871        7903 :   if (typ(z0) == t_INT && !is_bigint(z0))
    1872             :   {
    1873        4725 :     long m = itos(z0);
    1874        4725 :     if (abscmpiu(muluu(prec2nbits(prec),L), labs(m)) > 0)
    1875             :     { /* psi(m+x) = psi(1+x) + sum_{1 <= i < m} 1/(i+x) for m > 0
    1876             :                     psi(1+x) - sum_{0 <= i < -m} 1/(i+x) for m <= 0 */
    1877        4725 :       GEN H = NULL;
    1878        4725 :       if (m <= 0) L--; /* lose series accuracy due to 1/x term */
    1879        4725 :       if (L)
    1880             :       {
    1881        4718 :         Q = serpsi1(L, v, prec);
    1882        4718 :         if (m && m != 1) { H = Hseries(m, L, v, prec); Q = gadd(Q, H); }
    1883        4718 :         if (m <= 0) Q = gsub(Q, ginv(pol_x(v)));
    1884             :       }
    1885             :       else
    1886             :       {
    1887           7 :         Q = scalarser(gen_m1, v, 1);
    1888           7 :         setvalp(Q,-1);
    1889             :       }
    1890             :     }
    1891             :   }
    1892        7903 :   if (!Q)
    1893             :   { /* use psi(1-y)=psi(y)+Pi*cotan(Pi*y) ? */
    1894        3178 :     if (gcmp(real_i(z0),ghalf) < 0) { z0 = gsubsg(1,z0); Y = gsubsg(1,y); }
    1895        3178 :     Q = serpsiz0(z0, L, v, prec);
    1896             :   }
    1897        7903 :   Y2 = serchop0(Y); if (signe(Y2)) Q = gsubst(Q, v, Y2);
    1898             :   /* psi(z0 + Y2) = psi(Y) */
    1899        7903 :   if (Y != y)
    1900             :   { /* psi(y) = psi(Y) + Pi cotan(Pi Y) */
    1901          49 :     GEN pi = mppi(prec);
    1902          49 :     if (typ(z0) == t_INT) Y = Y2; /* in this case cotan(Pi*Y2) = cotan(Pi*Y) */
    1903          49 :     Q = gadd(Q, gmul(pi, gcotan(gmul(pi,Y), prec)));
    1904             :   }
    1905        7903 :   return Q;
    1906             : }
    1907             : 
    1908             : GEN
    1909        7945 : gpsi(GEN x, long prec)
    1910             : {
    1911             :   pari_sp av;
    1912             :   GEN y;
    1913        7945 :   switch(typ(x))
    1914             :   {
    1915          21 :     case t_REAL: case t_COMPLEX: return cxpsi(x,prec);
    1916             :     default:
    1917        7924 :       av = avma; if (!(y = toser_i(x))) break;
    1918        7910 :       return gerepileupto(av, serpsi(y,prec));
    1919             :   }
    1920          14 :   return trans_eval("psi",gpsi,x,prec);
    1921             : }

Generated by: LCOV version 1.11