1 /* $Id: trans2.c 10281 2008-06-09 11:10:14Z kb $
2 
3 Copyright (C) 2000  The PARI group.
4 
5 This file is part of the PARI/GP package.
6 
7 PARI/GP is free software; you can redistribute it and/or modify it under the
8 terms of the GNU General Public License as published by the Free Software
9 Foundation. It is distributed in the hope that it will be useful, but WITHOUT
10 ANY WARRANTY WHATSOEVER.
11 
12 Check the License for details. You should have received a copy of it, along
13 with the package; see the file 'COPYING'. If not, write to the Free Software
14 Foundation, Inc., 59 Temple Place - Suite 330, Boston, MA 02111-1307, USA. */
15 
16 /********************************************************************/
17 /**                                                                **/
18 /**                   TRANSCENDENTAL FUNCTIONS                     **/
19 /**                          (part 2)                              **/
20 /**                                                                **/
21 /********************************************************************/
22 #include "pari.h"
23 #include "paripriv.h"
24 
25 /********************************************************************/
26 /**                                                                **/
27 /**                          ARCTANGENT                            **/
28 /**                                                                **/
29 /********************************************************************/
30 static GEN
mpatan(GEN x)31 mpatan(GEN x)
32 {
33   long l, l1, l2, n, m, i, lp, e, s, sx = signe(x);
34   pari_sp av0, av;
35   double alpha, beta, delta;
36   GEN y, p1, p2, p3, p4, p5, unr;
37   int inv;
38 
39   if (!sx) return real_0_bit(expo(x));
40   l = lp = lg(x);
41   if (absrnz_egal1(x)) { /* |x| = 1 */
42     y = Pi2n(-2, l+1); if (sx < 0) setsigne(y,-1);
43     return y;
44   }
45   if (l > AGM_ATAN_LIMIT)
46   {
47     av = avma;
48     return gerepileuptoleaf(av, (GEN)logagmcx(mkcomplex(gen_1, x), l)[2]);
49   }
50 
51   e = expo(x); inv = (e >= 0); /* = (|x| > 1 ) */
52   if (e > 0) lp += (e>>TWOPOTBITS_IN_LONG);
53 
54   y = cgetr(lp); av0 = avma;
55   p1 = cgetr(l+1); affrr(x,p1); setsigne(p1, 1); /* p1 = |x| */
56   if (inv) p1 = divsr(1, p1);
57   e = expo(p1);
58   if (e < -100)
59     alpha = 1.65149612947 - e; /* log_2(Pi) - e */
60   else
61     alpha = log2(PI / atan(rtodbl(p1)));
62   beta = (double)(bit_accuracy(l)>>1);
63   delta = 1 + beta - alpha/2;
64   if (delta <= 0) { n = 1; m = 0; }
65   else
66   {
67     double fi = alpha-2;
68 #if 0
69     const double gama = 1.; /* optimize this */
70     if (delta >= gama*fi*fi)
71     {
72       n = (long)(1+sqrt(gama*delta));
73       m = (long)(1+sqrt(delta/gama) - fi);
74     }
75 #else
76     if (delta >= fi*fi)
77     {
78       double t = 1 + sqrt(delta);
79       n = (long)t;
80       m = (long)(t - fi);
81     }
82 #endif
83     else
84     {
85       n = (long)(1+beta/fi);
86       m = 0;
87     }
88   }
89   l2 = l+1+(m>>TWOPOTBITS_IN_LONG);
90   p2 = cgetr(l2); affrr(p1,p2); av = avma;
91   for (i=1; i<=m; i++)
92   {
93     p5 = addsr(1, mulrr(p2,p2)); setlg(p5,l2);
94     p5 = addsr(1, sqrtr_abs(p5)); setlg(p5,l2);
95     affrr(divrr(p2,p5), p2); avma = av;
96   }
97   p3 = mulrr(p2,p2); l1 = 4;
98   unr = real_1(l2); setlg(unr,4);
99   p4 = cgetr(l2); setlg(p4,4);
100   affrr(divrs(unr,2*n+1), p4);
101   s = 0; e = expo(p3); av = avma;
102   for (i = n; i > 1; i--) /* n >= 1. i = 1 done outside for efficiency */
103   {
104     setlg(p3,l1); p5 = mulrr(p4,p3);
105     s -= e; l1 += (s>>TWOPOTBITS_IN_LONG);
106     s %= BITS_IN_LONG;
107     if (l1 > l2) l1 = l2;
108     setlg(unr,l1); p5 = subrr(divrs(unr,2*i-1), p5);
109     setlg(p4,l1); affrr(p5,p4); avma = av;
110   }
111   setlg(p3, l2); p5 = mulrr(p4,p3); /* i = 1 */
112   setlg(unr,l2); p4 = subrr(unr, p5);
113 
114   p4 = mulrr(p2,p4); setexpo(p4, expo(p4)+m);
115   if (inv) p4 = subrr(Pi2n(-1, lp), p4);
116   if (sx < 0) setsigne(p4,-signe(p4));
117   affr_fixlg(p4,y); avma = av0; return y;
118 }
119 
120 GEN
gatan(GEN x,long prec)121 gatan(GEN x, long prec)
122 {
123   pari_sp av;
124   GEN a, y;
125 
126   switch(typ(x))
127   {
128     case t_REAL:
129       return mpatan(x);
130 
131     case t_COMPLEX:
132       av = avma; return gerepilecopy(av, mulcxmI(gath(mulcxI(x),prec)));
133 
134     case t_INTMOD: case t_PADIC: pari_err(typeer,"gatan");
135 
136     default:
137       av = avma; if (!(y = toser_i(x))) break;
138       if (valp(y) < 0) pari_err(negexper,"gatan");
139       if (lg(y)==2) return gcopy(y);
140       /* lg(y) > 2 */
141       a = integ(gdiv(derivser(y), gaddsg(1,gsqr(y))), varn(y));
142       if (!valp(y)) a = gadd(a, gatan(gel(y,2),prec));
143       return gerepileupto(av, a);
144   }
145   return transc(gatan,x,prec);
146 }
147 /********************************************************************/
148 /**                                                                **/
149 /**                             ARCSINE                            **/
150 /**                                                                **/
151 /********************************************************************/
152 /* |x| < 1, x != 0 */
153 static GEN
mpasin(GEN x)154 mpasin(GEN x) {
155   pari_sp av = avma;
156   GEN z, a = sqrtr(subsr(1, mulrr(x,x)));
157   if (lg(x) > AGM_ATAN_LIMIT)
158     z = (GEN)logagmcx(mkcomplex(a,x), lg(x))[2];
159   else
160     z = mpatan(divrr(x, a));
161   return gerepileuptoleaf(av, z);
162 }
163 
164 static GEN mpach(GEN x);
165 GEN
gasin(GEN x,long prec)166 gasin(GEN x, long prec)
167 {
168   long sx;
169   pari_sp av;
170   GEN a, y, p1;
171 
172   switch(typ(x))
173   {
174     case t_REAL: sx = signe(x);
175       if (!sx) return real_0_bit(expo(x));
176       if (absrnz_egal1(x)) { /* |x| = 1 */
177         if (sx > 0) return Pi2n(-1, lg(x)); /* 1 */
178         y = Pi2n(-1, lg(x)); setsigne(y, -1); return y; /* -1 */
179       }
180       if (expo(x) < 0) return mpasin(x);
181       y = cgetg(3,t_COMPLEX);
182       gel(y,1) = Pi2n(-1, lg(x));
183       gel(y,2) = mpach(x);
184       if (sx < 0)
185       {
186         setsigne(y[1],-signe(y[1]));
187         setsigne(y[2],-signe(y[2]));
188       }
189       return y;
190 
191     case t_COMPLEX:
192       av = avma;
193       return gerepilecopy(av, mulcxmI(gash(mulcxI(x), prec)));
194 
195     case t_INTMOD: case t_PADIC: pari_err(typeer,"gasin");
196     default:
197       av = avma; if (!(y = toser_i(x))) break;
198       if (gcmp0(y)) return gcopy(y);
199       /* lg(y) > 2*/
200       if (valp(y) < 0) pari_err(negexper,"gasin");
201       p1 = gsubsg(1,gsqr(y));
202       if (gcmp0(p1))
203       {
204         GEN t = Pi2n(-1,prec);
205         if (gsigne(gel(y,2)) < 0) setsigne(t, -1);
206         return gerepileupto(av, scalarser(t, varn(y), valp(p1)>>1));
207       }
208       p1 = gdiv(derivser(y), gsqrt(p1,prec));
209       a = integ(p1,varn(y));
210       if (!valp(y)) a = gadd(a, gasin(gel(y,2),prec));
211       return gerepileupto(av, a);
212   }
213   return transc(gasin,x,prec);
214 }
215 /********************************************************************/
216 /**                                                                **/
217 /**                             ARCCOSINE                          **/
218 /**                                                                **/
219 /********************************************************************/
220 static GEN
acos0(long e)221 acos0(long e) {
222   long l = e >> TWOPOTBITS_IN_LONG; if (l >= 0) l = -1;
223   return Pi2n(-1, 2-l);
224 }
225 
226 /* |x| < 1, x != 0 */
227 static GEN
mpacos(GEN x)228 mpacos(GEN x)
229 {
230   pari_sp av = avma;
231   GEN z, a = sqrtr(subsr(1, mulrr(x,x)));
232   if (lg(x) > AGM_ATAN_LIMIT)
233     z = (GEN)logagmcx(mkcomplex(x,a), lg(x))[2];
234   else {
235     z = mpatan(divrr(a, x));
236     if (signe(x) < 0) z = addrr(mppi(lg(z)), z);
237   }
238   return gerepileuptoleaf(av, z);
239 }
240 
241 GEN
gacos(GEN x,long prec)242 gacos(GEN x, long prec)
243 {
244   long sx;
245   pari_sp av;
246   GEN a, y, p1;
247 
248   switch(typ(x))
249   {
250     case t_REAL: sx = signe(x);
251       if (!sx) return acos0(expo(x));
252       if (absrnz_egal1(x)) /* |x| = 1 */
253         return sx > 0? real_0_bit( -(bit_accuracy(lg(x))>>1) ) : mppi(lg(x));
254       if (expo(x) < 0) return mpacos(x);
255 
256       y = cgetg(3,t_COMPLEX); p1 = mpach(x);
257       if (sx < 0) gel(y,1) = mppi(lg(x));
258       else {
259 	gel(y,1) = gen_0;
260         setsigne(p1,-signe(p1));
261       }
262       gel(y,2) = p1; return y;
263 
264     case t_COMPLEX: av = avma;
265       return gerepilecopy(av, mulcxmI(gach(x,prec)));
266 
267     case t_INTMOD: case t_PADIC: pari_err(typeer,"gacos");
268     case t_SER:
269       av = avma; if (!(y = toser_i(x))) break;
270       if (valp(y) < 0) pari_err(negexper,"gacos");
271       if (lg(y) > 2)
272       {
273 	p1 = gsubsg(1,gsqr(y));
274 	if (gcmp0(p1)) return zeroser(varn(y), valp(p1)>>1);
275 	p1 = integ(gdiv(gneg(derivser(y)), gsqrt(p1,prec)), varn(y));
276 	if (gcmp1(gel(y,2)) && !valp(y)) /*y = 1+O(y^k), k>=1*/
277 	  return gerepileupto(av, p1);
278       }
279       else p1 = y;
280       a = (lg(y)==2 || valp(y))? Pi2n(-1, prec): gacos(gel(y,2),prec);
281       return gerepileupto(av, gadd(a,p1));
282   }
283   return transc(gacos,x,prec);
284 }
285 /********************************************************************/
286 /**                                                                **/
287 /**                            ARGUMENT                            **/
288 /**                                                                **/
289 /********************************************************************/
290 
291 /* we know that x and y are not both 0 */
292 static GEN
mparg(GEN x,GEN y)293 mparg(GEN x, GEN y)
294 {
295   long prec, sx = signe(x), sy = signe(y);
296   GEN z;
297 
298   if (!sy)
299   {
300     if (sx > 0) return real_0_bit(expo(y) - expo(x));
301     return mppi(lg(x));
302   }
303   prec = lg(y); if (prec < lg(x)) prec = lg(x);
304   if (!sx)
305   {
306     z = Pi2n(-1, prec); if (sy < 0) setsigne(z,-1);
307     return z;
308   }
309 
310   if (expo(x)-expo(y) > -2)
311   {
312     z = mpatan(divrr(y,x)); if (sx > 0) return z;
313     return addrr_sign(z, signe(z), mppi(prec), sy);
314   }
315   z = mpatan(divrr(x,y));
316   return addrr_sign(z, -signe(z), Pi2n(-1, prec), sy);
317 }
318 
319 static GEN
rfix(GEN x,long prec)320 rfix(GEN x,long prec)
321 {
322   switch(typ(x))
323   {
324     case t_INT: return itor(x, prec);
325     case t_FRAC: return rdivii(gel(x,1),gel(x,2), prec);
326     case t_REAL: break;
327     default: pari_err(typeer,"rfix (conversion to t_REAL)");
328   }
329   return x;
330 }
331 
332 static GEN
cxarg(GEN x,GEN y,long prec)333 cxarg(GEN x, GEN y, long prec)
334 {
335   pari_sp av = avma;
336   x = rfix(x,prec);
337   y = rfix(y,prec); return gerepileuptoleaf(av, mparg(x,y));
338 }
339 
340 GEN
garg(GEN x,long prec)341 garg(GEN x, long prec)
342 {
343   long tx = typ(x);
344   pari_sp av;
345 
346   if (gcmp0(x)) pari_err(talker,"zero argument in garg");
347   switch(tx)
348   {
349     case t_REAL: prec = lg(x); /* fall through */
350     case t_INT: case t_FRAC:
351       return (gsigne(x)>0)? real_0(prec): mppi(prec);
352 
353     case t_QUAD:
354       av = avma;
355       return gerepileuptoleaf(av, garg(quadtoc(x, prec), prec));
356 
357     case t_COMPLEX:
358       return cxarg(gel(x,1),gel(x,2),prec);
359 
360     case t_VEC: case t_COL: case t_MAT:
361       return transc(garg,x,prec);
362   }
363   pari_err(typeer,"garg");
364   return NULL; /* not reached */
365 }
366 
367 /********************************************************************/
368 /**                                                                **/
369 /**                      HYPERBOLIC COSINE                         **/
370 /**                                                                **/
371 /********************************************************************/
372 
373 static GEN
mpch(GEN x)374 mpch(GEN x)
375 {
376   pari_sp av;
377   GEN z;
378 
379   if (gcmp0(x)) { /* 1 + x */
380     long e = expo(x);
381     if (e > 0) return real_0_bit(e);
382     return real_1(3 + ((-e)>>TWOPOTBITS_IN_LONG));
383   }
384   av = avma;
385   z = mpexp(x); z = addrr(z, ginv(z)); setexpo(z, expo(z)-1);
386   return gerepileuptoleaf(av, z);
387 }
388 
389 GEN
gch(GEN x,long prec)390 gch(GEN x, long prec)
391 {
392   pari_sp av;
393   GEN y, p1;
394 
395   switch(typ(x))
396   {
397     case t_REAL: return mpch(x);
398     case t_COMPLEX: case t_PADIC:
399       av = avma; p1 = gexp(x,prec); p1 = gadd(p1, ginv(p1));
400       return gerepileupto(av, gmul2n(p1,-1));
401     case t_INTMOD: pari_err(typeer,"gch");
402     default:
403       av = avma; if (!(y = toser_i(x))) break;
404       if (gcmp0(y) && valp(y) == 0) return gcopy(y);
405       p1 = gexp(y,prec); p1 = gadd(p1, ginv(p1));
406       return gerepileupto(av, gmul2n(p1,-1));
407   }
408   return transc(gch,x,prec);
409 }
410 /********************************************************************/
411 /**                                                                **/
412 /**                       HYPERBOLIC SINE                          **/
413 /**                                                                **/
414 /********************************************************************/
415 
416 static GEN
mpsh(GEN x)417 mpsh(GEN x)
418 {
419   pari_sp av;
420   long ex = expo(x), lx;
421   GEN z, res;
422 
423   if (!signe(x)) return real_0_bit(ex);
424   lx = lg(x); res = cgetr(lx); av = avma;
425   if (ex < 1 - BITS_IN_LONG) x = rtor(x, lx + nbits2nlong(-ex)-1);
426   z = mpexp(x); z = addrr(z, divsr(-1,z)); setexpo(z, expo(z)-1);
427   affrr(z, res); avma = av; return res;
428 }
429 
430 GEN
gsh(GEN x,long prec)431 gsh(GEN x, long prec)
432 {
433   pari_sp av;
434   GEN y, p1;
435 
436   switch(typ(x))
437   {
438     case t_REAL: return mpsh(x);
439     case t_COMPLEX: case t_PADIC:
440       av = avma; p1 = gexp(x,prec); p1 = gsub(p1, ginv(p1));
441       return gerepileupto(av, gmul2n(p1,-1));
442     case t_INTMOD:
443     default:
444       av = avma; if (!(y = toser_i(x))) break;
445       if (gcmp0(y) && valp(y) == 0) return gcopy(y);
446       p1 = gexp(y, prec); p1 = gsub(p1, ginv(p1));
447       return gerepileupto(av, gmul2n(p1,-1));
448   }
449   return transc(gsh,x,prec);
450 }
451 /********************************************************************/
452 /**                                                                **/
453 /**                      HYPERBOLIC TANGENT                        **/
454 /**                                                                **/
455 /********************************************************************/
456 
457 static GEN
mpth(GEN x)458 mpth(GEN x)
459 {
460   long lx, s = signe(x);
461   GEN y;
462 
463   if (!s) return real_0_bit(expo(x));
464   lx = lg(x);
465   if (absr_cmp(x, stor(bit_accuracy(lx), 3)) >= 0) {
466     y = real_1(lx);
467   } else {
468     pari_sp av = avma;
469     long ex = expo(x);
470     GEN t;
471     if (ex < 1 - BITS_IN_LONG) x = rtor(x, lx + nbits2nlong(-ex)-1);
472     t = exp1r_abs(gmul2n(x,1)); /* exp(|2x|) - 1 */
473     y = gerepileuptoleaf(av, divrr(t, addsr(2,t)));
474   }
475   if (s < 0) togglesign(y); /* tanh is odd */
476   return y;
477 }
478 
479 GEN
gth(GEN x,long prec)480 gth(GEN x, long prec)
481 {
482   pari_sp av;
483   GEN y, t;
484 
485   switch(typ(x))
486   {
487     case t_REAL: return mpth(x);
488     case t_COMPLEX: case t_PADIC:
489       av = avma;
490       t = gexp(gmul2n(x,1),prec);
491       t = gdivsg(-2, gaddgs(t,1));
492       return gerepileupto(av, gaddsg(1,t));
493     case t_INTMOD: pari_err(typeer,"gth");
494     default:
495       av = avma; if (!(y = toser_i(x))) break;
496       if (gcmp0(y)) return gcopy(y);
497       t = gexp(gmul2n(y, 1),prec);
498       t = gdivsg(-2, gaddgs(t,1));
499       return gerepileupto(av, gaddsg(1,t));
500   }
501   return transc(gth,x,prec);
502 }
503 /********************************************************************/
504 /**                                                                **/
505 /**                     ARG-HYPERBOLIC SINE                        **/
506 /**                                                                **/
507 /********************************************************************/
508 
509 /* x != 0 */
510 static GEN
mpash(GEN x)511 mpash(GEN x)
512 {
513   GEN z, res;
514   pari_sp av;
515   long lx = lg(x), ex = expo(x);
516 
517   res = cgetr(lx); av = avma;
518   if (ex < 1 - BITS_IN_LONG) x = rtor(x, lx + nbits2nlong(-ex)-1);
519   z = logr_abs( addrr_sign(x,1, sqrtr( addrs(mulrr(x,x), 1) ), 1) );
520   if (signe(x) < 0) togglesign(z);
521   affrr(z, res); avma = av; return res;
522 }
523 
524 GEN
gash(GEN x,long prec)525 gash(GEN x, long prec)
526 {
527   long sx, sy, sz;
528   pari_sp av;
529   GEN a, y, p1;
530 
531   if (gcmp0(x)) return gcopy(x);
532   switch(typ(x))
533   {
534     case t_REAL:
535       return mpash(x);
536 
537     case t_COMPLEX: av = avma;
538       p1 = gadd(x, gsqrt(gaddsg(1,gsqr(x)), prec));
539       y = glog(p1,prec);
540       sz = (typ(y)==t_COMPLEX)? gsigne(gel(y,1)): gsigne(y);
541       if (typ(p1) == t_COMPLEX) {
542         sx = gsigne(gel(p1,1));
543         sy = gsigne(gel(p1,2));
544       } else {
545         sx = gsigne(p1);
546         sy = 0;
547       }
548       if (sx > 0 || (!sx && sy*sz<=0)) return gerepileupto(av, y);
549 
550       p1 = mppi(prec); if (sy<0) setsigne(p1,-1);
551       return gerepileupto(av, gadd(gneg_i(y), pureimag(p1)));
552     case t_INTMOD: case t_PADIC: pari_err(typeer,"gash");
553     default:
554       av = avma; if (!(y = toser_i(x))) break;
555       if (gcmp0(y)) return gcopy(y);
556       if (valp(y) < 0) pari_err(negexper,"gash");
557       p1 = gaddsg(1,gsqr(y));
558       if (gcmp0(p1))
559       {
560         GEN t = PiI2n(-1,prec);
561         if ( gsigne(imag_i(gel(y,2))) < 0 ) setsigne(gel(t,2), -1);
562         return gerepileupto(av, scalarser(t, varn(y), valp(p1)>>1));
563       }
564       p1 = gdiv(derivser(y), gsqrt(p1,prec));
565       a = integ(p1,varn(y));
566       if (!valp(y)) a = gadd(a, gash(gel(y,2),prec));
567       return gerepileupto(av, a);
568   }
569   return transc(gash,x,prec);
570 }
571 /********************************************************************/
572 /**                                                                **/
573 /**                     ARG-HYPERBOLIC COSINE                      **/
574 /**                                                                **/
575 /********************************************************************/
576 
577 /* |x| >= 1, return ach(|x|) */
578 static GEN
mpach(GEN x)579 mpach(GEN x)
580 {
581   pari_sp av = avma;
582   GEN z = logr_abs( addrr_sign(x, 1, sqrtr( subrs(mulrr(x,x), 1) ), 1) );
583   return gerepileuptoleaf(av, z);
584 }
585 
586 GEN
gach(GEN x,long prec)587 gach(GEN x, long prec)
588 {
589   pari_sp av;
590   GEN a, y, p1;
591   long v;
592 
593   switch(typ(x))
594   {
595     case t_REAL:
596       if (signe(x) == 0) { y=cgetimag(); gel(y,2) = acos0(expo(x)); return y; }
597       if (signe(x) > 0 && expo(x) >= 0) return mpach(x); /* x >= 1 */
598       /* -1 < x < 1 */
599       if (expo(x) < 0) { y = cgetimag(); gel(y,2) = mpacos(x); return y; }
600       /* x <= -1 */
601       if (absrnz_egal1(x)) { y = cgetimag(); gel(y,2) = mppi(lg(x)); return y; }
602       y = cgetg(3,t_COMPLEX);
603       av = avma; p1 = mpach(x);
604       setsigne(p1, -signe(p1));
605       gel(y,1) = p1;
606       gel(y,2) = mppi(lg(x)); return y;
607 
608     case t_COMPLEX:
609       av = avma;
610       p1 = gadd(x, gsqrt(gaddsg(-1,gsqr(x)), prec)); /* x + sqrt(x^2-1) */
611       y = glog(p1,prec);
612       if (typ(y) == t_COMPLEX && signe(y[2]) < 0) y = gneg(y);
613       return gerepileupto(av, y);
614 
615     case t_INTMOD: case t_PADIC: pari_err(typeer,"gach");
616 
617     default:
618       av = avma; if (!(y = toser_i(x))) break;
619       v = valp(y);
620       if (v < 0) pari_err(negexper,"gach");
621       if (gcmp0(y))
622       {
623         if (!v) return gcopy(y);
624         return gerepileupto(av, gadd(y, PiI2n(-1, prec)));
625       }
626       p1 = gsubgs(gsqr(y),1);
627       if (gcmp0(p1)) { avma = av; return zeroser(varn(y), valp(p1)>>1); }
628       p1 = gdiv(derivser(y), gsqrt(p1,prec));
629       a = integ(p1, varn(y));
630       if (v)
631         p1 = PiI2n(-1, prec); /* I Pi/2 */
632       else
633       {
634         p1 = gel(y,2); if (gcmp1(p1)) return gerepileupto(av,a);
635         p1 = gach(p1, prec);
636       }
637       return gerepileupto(av, gadd(p1,a));
638   }
639   return transc(gach,x,prec);
640 }
641 /********************************************************************/
642 /**                                                                **/
643 /**                     ARG-HYPERBOLIC TANGENT                     **/
644 /**                                                                **/
645 /********************************************************************/
646 
647 /* |x| < 1, x != 0 */
648 static GEN
mpath(GEN x)649 mpath(GEN x)
650 {
651   pari_sp av = avma;
652   long ex = expo(x);
653   GEN z;
654   if (ex < 1 - BITS_IN_LONG) x = rtor(x, lg(x) + nbits2nlong(-ex)-1);
655   z = logr_abs( addrs(divsr(2,subsr(1,x)), -1) );
656   setexpo(z, expo(z)-1); return gerepileuptoleaf(av, z);
657 }
658 
659 GEN
gath(GEN x,long prec)660 gath(GEN x, long prec)
661 {
662   pari_sp av;
663   GEN a, y, p1;
664 
665   switch(typ(x))
666   {
667     case t_REAL:
668       if (!signe(x)) return real_0_bit(expo(x));
669       if (expo(x) < 0) return mpath(x);
670 
671       y = cgetg(3,t_COMPLEX);
672       av = avma;
673       p1 = addrs(divsr(2,addsr(-1,x)),1);
674       if (!signe(p1)) pari_err(talker,"singular argument in atanh");
675       p1 = logr_abs(p1);
676       setexpo(p1, expo(p1)-1);
677       gel(y,1) = gerepileuptoleaf(av, p1);
678       gel(y,2) = Pi2n(-1, lg(x)); return y;
679 
680     case t_COMPLEX:
681       av = avma; p1 = glog( gaddgs(gdivsg(2,gsubsg(1,x)),-1), prec );
682       return gerepileupto(av, gmul2n(p1,-1));
683 
684     case t_INTMOD: case t_PADIC: pari_err(typeer,"gath");
685     default:
686       av = avma; if (!(y = toser_i(x))) break;
687       if (valp(y) < 0) pari_err(negexper,"gath");
688       p1 = gdiv(derivser(y), gsubsg(1,gsqr(y)));
689       a = integ(p1, varn(y));
690       if (!valp(y)) a = gadd(a, gath(gel(y,2),prec));
691       return gerepileupto(av, a);
692   }
693   return transc(gath,x,prec);
694 }
695 /********************************************************************/
696 /**                                                                **/
697 /**               CACHE BERNOULLI NUMBERS B_2k                     **/
698 /**                                                                **/
699 /********************************************************************/
700 /* is B_{2k} precomputed at precision >= prec ? */
701 int
OK_bern(long k,long prec)702 OK_bern(long k, long prec)
703 {
704   return (bernzone && bernzone[1] >= k && bernzone[2] >= prec);
705 }
706 
707 #define BERN(i)       (B + 3 + (i)*B[2])
708 #define set_bern(c0, i, B) STMT_START { \
709   *(BERN(i)) = c0; affrr(B, BERN(i)); } STMT_END
710 /* compute B_0,B_2,...,B_2*nb */
711 void
mpbern(long nb,long prec)712 mpbern(long nb, long prec)
713 {
714   long i, l, c0;
715   pari_sp av;
716   GEN B;
717   pari_timer T;
718 
719   prec++; /* compute one more word of accuracy than required */
720   if (OK_bern(nb, prec)) return;
721   if (nb < 0) nb = 0;
722   l = 3 + prec*(nb+1);
723   B = newbloc(l);
724   B[0] = evaltyp(t_STR) | evallg(l); /* dummy non-recursive type */
725   B[1] = nb;
726   B[2] = prec;
727   av = avma;
728 
729   c0 = evaltyp(t_REAL) | evallg(prec);
730   *(BERN(0)) = c0; affsr(1, BERN(0));
731   if (bernzone && bernzone[2] >= prec)
732   { /* don't recompute known Bernoulli */
733     for (i = 1; i <= bernzone[1]; i++) set_bern(c0, i, bern(i));
734   }
735   else i = 1;
736   if (DEBUGLEVEL) {
737     fprintferr("caching Bernoulli numbers 2*%ld to 2*%ld, prec = %ld\n",
738                i,nb,prec);
739     TIMERstart(&T);
740   }
741 
742   if (i == 1 && nb > 0)
743   {
744     set_bern(c0, 1, divrs(real_1(prec), 6)); /* B2 = 1/6 */
745     i = 2;
746   }
747   for (   ; i <= nb; i++, avma = av)
748   { /* i > 1 */
749     long n = 8, m = 5, d1 = i-1, d2 = 2*i-3;
750     GEN S = BERN(d1);
751 
752     for (;;)
753     {
754       S = divrs(mulrs(S, n*m), d1*d2);
755       if (d1 == 1) break;
756       n += 4; m += 2; d1--; d2 -= 2;
757       S = addrr(BERN(d1), S);
758       if ((d1 & 127) == 0) { set_bern(c0, i, S); S = BERN(i); avma = av; }
759     }
760     S = divrs(subsr(2*i, S), 2*i+1);
761     setexpo(S, expo(S) - 2*i);
762     set_bern(c0, i, S); /* S = B_2i */
763   }
764   if (DEBUGLEVEL) msgTIMER(&T, "Bernoulli");
765   if (bernzone) gunclone(bernzone);
766   avma = av; bernzone = B;
767 }
768 #undef BERN
769 
770 GEN
bernreal(long n,long prec)771 bernreal(long n, long prec)
772 {
773   GEN B;
774 
775   if (n==1) { B = stor(-1, prec); setexpo(B,-1); return B; }
776   if (n<0 || n&1) return gen_0;
777   n >>= 1; mpbern(n+1,prec); B=cgetr(prec);
778   affrr(bern(n),B); return B;
779 }
780 
781 #if 0
782 /* k > 0 */
783 static GEN
784 bernfracspec(long k)
785 {
786   ulong n, K = k+1;
787   pari_sp av, lim;
788   GEN s, c, N, b;
789 
790   c = N = utoipos(K); s = gen_1; b = gen_0;
791   av = avma; lim = stack_lim(av,2);
792   for (n=2; ; n++) /* n <= k+1 */
793   {
794     c = diviiexact(muliu(c,k+2-n), utoipos(n));
795     if (n & 1) setsigne(c, 1); else setsigne(c, -1);
796     /* c = (-1)^(n-1) binomial(k+1, n),  s = 1^k + ... + (n-1)^k */
797 
798     b = gadd(b, gdivgs(mulii(c,s), n));
799     if (n == K) return gerepileupto(av, b);
800 
801     gel(N,2) = n; s = addii(s, powiu(N,k));
802     if (low_stack(lim, stack_lim(av,2)))
803     {
804       if (DEBUGMEM>1) pari_warn(warnmem,"bernfrac");
805       gerepileall(av,3, &c,&b,&s);
806     }
807   }
808 }
809 #endif
810 
811 static GEN
B2(void)812 B2(void){ GEN z = cgetg(3, t_FRAC);
813   gel(z,1) = gen_1;
814   gel(z,2) = utoipos(6); return z;
815 }
816 static GEN
B4(void)817 B4(void) { GEN z = cgetg(3, t_FRAC);
818   gel(z,1) = gen_m1;
819   gel(z,2) = utoipos(30); return z;
820 }
821 
822 GEN
bernfrac(long k)823 bernfrac(long k)
824 {
825   if (k < 6) switch(k)
826   {
827     case 0: return gen_1;
828     case 1: return gneg(ghalf);
829     case 2: return B2();
830     case 4: return B4();
831     default: return gen_0;
832   }
833   if (k & 1) return gen_0;
834   return bernfrac_using_zeta(k);
835 }
836 
837 /* mpbern as exact fractions */
838 static GEN
bernvec_old(long nb)839 bernvec_old(long nb)
840 {
841   long n, i;
842   GEN y;
843 
844   if (nb < 0) return cgetg(1, t_VEC);
845   if (nb > 46340 && BITS_IN_LONG == 32) pari_err(impl, "bernvec for n > 46340");
846 
847   y = cgetg(nb+2, t_VEC); gel(y,1) = gen_1;
848   for (n = 1; n <= nb; n++)
849   { /* compute y[n+1] = B_{2n} */
850     pari_sp av = avma;
851     GEN b = gmul2n(utoineg(2*n - 1), -1); /* 1 + (2n+1)B_1 = -(2n-1) /2 */
852     GEN c = gen_1;
853     ulong u1 = 2*n + 1, u2 = n, d1 = 1, d2 = 1;
854 
855     for (i = 1; i < n; i++)
856     {
857       c = diviiexact(muliu(c, u1*u2), utoipos(d1*d2));/*= binomial(2n+1, 2*i) */
858       b = gadd(b, gmul(c, gel(y,i+1)));
859       u1 -= 2; u2--; d1++; d2 += 2;
860     }
861     gel(y,n+1) = gerepileupto(av, gdivgs(b, -(1+2*n)));
862   }
863   return y;
864 }
865 GEN
bernvec(long nb)866 bernvec(long nb)
867 {
868   GEN y = cgetg(nb+2, t_VEC), z = y + 1;
869   long i;
870   if (nb < 20) return bernvec_old(nb);
871   for (i = nb; i > 2; i--) gel(z,i) = bernfrac_using_zeta(i << 1);
872   gel(y,3) = B4();
873   gel(y,2) = B2();
874   gel(y,1) = gen_1; return y;
875 }
876 
877 /********************************************************************/
878 /**                                                                **/
879 /**                         EULER'S GAMMA                          **/
880 /**                                                                **/
881 /********************************************************************/
882 
883 /* x / (i*(i+1)) */
884 GEN
divrsns(GEN x,long i)885 divrsns(GEN x, long i)
886 {
887 #ifdef LONG_IS_64BIT
888   if (i < 3037000500) /* i(i+1) < 2^63 */
889 #else
890   if (i < 46341) /* i(i+1) < 2^31 */
891 #endif
892     return divrs(x, i*(i+1));
893   else
894     return divrs(divrs(x, i), i+1);
895 }
896 /* x / (i*(i+1)) */
897 GEN
divgsns(GEN x,long i)898 divgsns(GEN x, long i)
899 {
900 #ifdef LONG_IS_64BIT
901   if (i < 3037000500) /* i(i+1) < 2^63 */
902 #else
903   if (i < 46341) /* i(i+1) < 2^31 */
904 #endif
905     return gdivgs(x, i*(i+1));
906   else
907     return gdivgs(gdivgs(x, i), i+1);
908 }
909 
910 /* arg(s+it) */
911 double
darg(double s,double t)912 darg(double s, double t)
913 {
914   double x;
915   if (!t) return (s>0)? 0.: PI;
916   if (!s) return (t>0)? PI/2: -PI/2;
917   x = atan(t/s);
918   return (s>0)? x
919               : ((t>0)? x+PI : x-PI);
920 }
921 
922 void
dcxlog(double s,double t,double * a,double * b)923 dcxlog(double s, double t, double *a, double *b)
924 {
925   *a = log(s*s + t*t) / 2; /* log |s| = Re(log(s)) */
926   *b = darg(s,t);          /* Im(log(s)) */
927 }
928 
929 double
dabs(double s,double t)930 dabs(double s, double t) { return sqrt( s*s + t*t ); }
931 double
dnorm(double s,double t)932 dnorm(double s, double t) { return s*s + t*t; }
933 
934 GEN
trans_fix_arg(long * prec,GEN * s0,GEN * sig,pari_sp * av,GEN * res)935 trans_fix_arg(long *prec, GEN *s0, GEN *sig, pari_sp *av, GEN *res)
936 {
937   GEN s, p1;
938   long l;
939   if (typ(*s0)==t_COMPLEX && gcmp0(gel(*s0,2))) *s0 = gel(*s0,1);
940   s = *s0;
941   l = precision(s); if (!l) l = *prec;
942   if (l < 3) l = 3;
943 
944   if (typ(s) == t_COMPLEX)
945   { /* s = sig + i t */
946     *res = cgetc(l); *av = avma;
947     s = ctofp(s, l+1); *sig = gel(s,1);
948   }
949   else /* real number */
950   {
951     *res = cgetr(l); *av = avma;
952     *sig = s = gtofp(s, l+1);
953     p1 = floorr(s);
954     if (!signe(subri(s,p1))) *s0 = p1;
955   }
956   *prec = l; return s;
957 }
958 
959 #if 0
960 /* x, z t_REAL. Compute unique x in ]-z,z] congruent to x mod 2z */
961 static GEN
962 red_mod_2z(GEN x, GEN z)
963 {
964   GEN Z = gmul2n(z, 1), d = subrr(z, x);
965   /* require little accuracy */
966   if (!signe(d)) return x;
967   setlg(d, 3 + ((expo(d) - expo(Z)) >> TWOPOTBITS_IN_LONG));
968   return addrr(mulir(floorr(divrr(d, Z)), Z), x);
969 }
970 #endif
971 
972 /* update lg(z) before affrr(y, z)  [ to cater for precision loss ]*/
973 void
affr_fixlg(GEN y,GEN z)974 affr_fixlg(GEN y, GEN z) {
975   long ly = lg(y), lz = lg(z);
976   if (ly < lz)
977   {
978     setlg(z, ly);
979     stackdummy((pari_sp)(z + lz), (pari_sp)(z + ly));
980   }
981   /* lz <= ly */
982   affrr(y, z);
983 }
984 
985 static GEN
cxgamma(GEN s0,int dolog,long prec)986 cxgamma(GEN s0, int dolog, long prec)
987 {
988   GEN s, u, a, y, res, tes, sig, invn2, p1, nnx, pi, pi2, sqrtpi2;
989   long i, lim, nn;
990   pari_sp av, av2, avlim;
991   int funeq = 0;
992 
993   if (DEBUGLEVEL>5) (void)timer2();
994   s = trans_fix_arg(&prec,&s0,&sig,&av,&res);
995 
996   if ((signe(sig) <= 0 || expo(sig) < -1)
997     && (typ(s) == t_REAL || gexpo(gel(s,2)) <= 16))
998   { /* s <--> 1-s */
999     funeq = 1; s = gsub(gen_1, s); sig = real_i(s);
1000   }
1001 
1002   { /* find "optimal" parameters [lim, nn] */
1003     double ssig = rtodbl(sig);
1004     double st = rtodbl(imag_i(s));
1005     double la, l,l2,u,v, rlogs, ilogs;
1006 
1007     dcxlog(ssig,st, &rlogs,&ilogs);
1008     /* Re (s - 1/2) log(s) */
1009     u = (ssig - 0.5)*rlogs - st * ilogs;
1010     /* Im (s - 1/2) log(s) */
1011     v = (ssig - 0.5)*ilogs + st * rlogs;
1012     /* l2 = | (s - 1/2) log(s) - s + log(2Pi)/2 |^2 ~ |lngamma(s))|^2 */
1013     u = u - ssig + log(2.*PI)/2;
1014     v = v - st;
1015     l2 = u*u + v*v;
1016     if (l2 < 0.000001) l2 = 0.000001;
1017     l = (bit_accuracy_mul(prec, LOG2) - log(l2)/2) / 2.;
1018     if (l < 0) l = 0.;
1019 
1020     la = 3.; /* FIXME: heuristic... */
1021     if (st > 1 && l > 0)
1022     {
1023       double t = st * PI / l;
1024       la = t * log(t);
1025       if (la < 3) la = 3.;
1026       if (la > 150) la = t;
1027     }
1028     lim = (long)ceil(l / (1.+ log(la)));
1029     if (lim == 0) lim = 1;
1030 
1031     u = (lim-0.5) * la / PI;
1032     l2 = u*u - st*st;
1033     if (l2 > 0)
1034     {
1035       nn = (long)ceil(sqrt(l2) - ssig);
1036       if (nn < 1) nn = 1;
1037     }
1038     else
1039       nn = 1;
1040 #if 0
1041 #define pariK2 (1.1239968) /* 1/(1-(log(2)/(2*pi))) */
1042 #define pariK4 (17.079468445347/BITS_IN_LONG) /* 2*e*pi/BIL */
1043     {/* same: old method */
1044       long e = gexpo(s);
1045       double beta;
1046       if (e > 1000)
1047       {
1048         nn = 0;
1049         beta = log(pariK4 / (prec-2)) / LOG2 + e;
1050         if (beta > 1.) beta += log(beta)/LOG2;
1051         lim = (long)((bit_accuracy(prec)>>1)/beta + 1);
1052       }
1053       else
1054       {
1055         double alpha = sqrt( dnorm(ssig, st) );
1056         beta = bit_accuracy_mul(prec,LOG2/(2*PI)) - alpha;
1057         if (beta >= 0) nn = (long)(1+pariK2*beta); else nn = 0;
1058         if (nn)
1059           lim = (long)(1+PI*(alpha+nn));
1060         else
1061         {
1062           beta = log( pariK4 * alpha / (prec-2) ) / LOG2;
1063           if (beta > 1.) beta += log(beta)/LOG2;
1064           lim = (long)((bit_accuracy(prec)>>1)/beta + 1);
1065         }
1066       }
1067       nn++;
1068     }
1069 #endif
1070     if (DEBUGLEVEL>5) fprintferr("lim, nn: [%ld, %ld], la = %lf\n",lim,nn,la);
1071   }
1072   prec++;
1073 
1074   av2 = avma; avlim = stack_lim(av2,3);
1075   y = s;
1076   if (typ(s0) == t_INT)
1077   {
1078     if (signe(s0) <= 0) pari_err(talker,"non-positive integer argument in cxgamma");
1079     if (is_bigint(s0)) {
1080       for (i=1; i < nn; i++)
1081       {
1082         y = mulri(y, addis(s0, i));
1083         if (low_stack(avlim,stack_lim(av2,3)))
1084         {
1085           if(DEBUGMEM>1) pari_warn(warnmem,"gamma");
1086           y = gerepileuptoleaf(av2, y);
1087         }
1088       }
1089     } else {
1090       ulong ss = itou(s0);
1091       for (i=1; i < nn; i++)
1092       {
1093         y = mulru(y, ss + i);
1094         if (low_stack(avlim,stack_lim(av2,3)))
1095         {
1096           if(DEBUGMEM>1) pari_warn(warnmem,"gamma");
1097           y = gerepileuptoleaf(av2, y);
1098         }
1099       }
1100     }
1101     if (dolog) y = logr_abs(y);
1102   }
1103   else if (!dolog || typ(s) == t_REAL)
1104   { /* Compute lngamma mod 2 I Pi */
1105     for (i=1; i < nn; i++)
1106     {
1107       y = gmul(y, gaddgs(s,i));
1108       if (low_stack(avlim,stack_lim(av2,3)))
1109       {
1110         if(DEBUGMEM>1) pari_warn(warnmem,"gamma");
1111         y = gerepileupto(av2, y);
1112       }
1113     }
1114     if (dolog) y = logr_abs(y);
1115   }
1116   else
1117   { /* dolog && complex s: be careful with imaginary part */
1118     y = glog(y, prec);
1119     for (i=1; i < nn; i++)
1120     {
1121       y = gadd(y, glog(gaddgs(s,i), prec));
1122       if (low_stack(avlim,stack_lim(av2,3)))
1123       {
1124         if(DEBUGMEM>1) pari_warn(warnmem,"gamma");
1125         y = gerepileupto(av2, y);
1126       }
1127     }
1128   }
1129   if (DEBUGLEVEL>5) msgtimer("product from 0 to N-1");
1130 
1131   nnx = gaddgs(s, nn);
1132   a = ginv(nnx); invn2 = gsqr(a);
1133   tes = divrsns(bernreal(2*lim,prec), 2*lim-1); /* B2l / (2l-1) 2l*/
1134   if (DEBUGLEVEL>5) msgtimer("Bernoullis");
1135   for (i = 2*lim-2; i > 1; i -= 2)
1136   {
1137     u = divrsns(bernreal(i,prec), i-1); /* Bi / i(i-1) */
1138     tes = gadd(u, gmul(invn2,tes));
1139   }
1140   if (DEBUGLEVEL>5) msgtimer("Bernoulli sum");
1141 
1142   p1 = gsub(gmul(gsub(nnx, ghalf), glog(nnx,prec)), nnx);
1143   p1 = gadd(p1, gmul(tes, a));
1144 
1145   pi = mppi(prec); pi2 = shiftr(pi, 1); sqrtpi2 = sqrtr(pi2);
1146 
1147   if (dolog)
1148   {
1149     if (funeq)
1150     { /* 2 Pi ceil( (2Re(s) - 3)/4 ) */
1151       GEN z = mulri(pi2, ceilr(shiftr(subrs(shiftr(sig,1), 3), -2)));
1152       /* y --> y + log Pi - log sqrt(2Pi) - log sin(Pi s)
1153        *     = y - log( sin(Pi s) / (sqrt(2Pi)/2) ) */
1154       y = gsub(y, glog(gdiv(gsin(gmul(pi,s0),prec), shiftr(sqrtpi2,-1)), prec));
1155       if (signe(z)) {
1156         if (gsigne(imag_i(s)) < 0) setsigne(z, -signe(z));
1157         if (typ(y) == t_COMPLEX)
1158           gel(y,2) = gadd(gel(y,2), z);
1159         else
1160           y = gadd(y, pureimag(z));
1161       }
1162       p1 = gneg(p1);
1163     }
1164     else /* y --> sqrt(2Pi) / y */
1165       y = gsub(logr_abs(sqrtpi2), y);
1166     y = gadd(p1, y);
1167   }
1168   else
1169   {
1170     if (funeq)
1171     { /* y --> y Pi/(sin(Pi s) * sqrt(2Pi)) = y sqrt(Pi/2)/sin(Pi s) */
1172       y = gdiv(gmul(shiftr(sqrtpi2,-1),y), gsin(gmul(pi,s0), prec));
1173       /* don't use s above: sin(pi s0) = sin(pi s) and the former is
1174        * more accurate, esp. if s0 ~ 0 */
1175       p1 = gneg(p1);
1176     }
1177     else /* y --> sqrt(2Pi) / y */
1178       y = gdiv(sqrtpi2, y);
1179     y = gmul(gexp(p1, prec), y);
1180   }
1181   if (typ(y) == t_REAL) affr_fixlg(y, res);
1182   else
1183   {
1184     if (typ(res) == t_REAL) return gerepileupto(av, y);
1185     affr_fixlg(gel(y,1), gel(res,1));
1186     affr_fixlg(gel(y,2), gel(res,2));
1187   }
1188   avma = av; return res;
1189 }
1190 
1191 /* Gamma((m+1) / 2) */
1192 static GEN
gammahs(long m,long prec)1193 gammahs(long m, long prec)
1194 {
1195   GEN y = cgetr(prec), z;
1196   pari_sp av = avma;
1197   long ma = labs(m);
1198 
1199   if (ma > 200 + 50*(prec-2)) /* heuristic */
1200   {
1201     z = stor(m + 1, prec); setexpo(z, expo(z)-1);
1202     affrr(cxgamma(z,0,prec), y);
1203     avma = av; return y;
1204   }
1205   z = sqrtr( mppi(prec) );
1206   if (m)
1207   {
1208     GEN p1 = seq_umul(ma/2 + 1, ma);
1209     long v = vali(p1);
1210     p1 = shifti(p1, -v); v -= ma;
1211     if (m >= 0) z = mulri(z,p1);
1212     else
1213     {
1214       z = divri(z,p1); v = -v;
1215       if ((m&3) == 2) setsigne(z,-1);
1216     }
1217     setexpo(z, expo(z) + v);
1218   }
1219   affrr(z, y); avma = av; return y;
1220 }
1221 GEN
ggamd(GEN x,long prec)1222 ggamd(GEN x, long prec)
1223 {
1224   pari_sp av, tetpil;
1225 
1226   switch(typ(x))
1227   {
1228     case t_INT:
1229     {
1230       long k = itos(x);
1231       if (labs(k) > 962353) pari_err(talker, "argument too large in ggamd");
1232       return gammahs(k<<1, prec);
1233     }
1234     case t_REAL: case t_FRAC: case t_COMPLEX: case t_QUAD: case t_PADIC:
1235       av=avma; x = gadd(x,ghalf); tetpil=avma;
1236       return gerepile(av,tetpil,ggamma(x,prec));
1237 
1238     case t_INTMOD: pari_err(typeer,"ggamd");
1239     case t_SER: pari_err(impl,"gamd of a power series");
1240   }
1241   return transc(ggamd,x,prec);
1242 }
1243 
1244 /* find n such that n-v_p(n!)>=k */
nboft(long k,long p)1245 static long nboft(long k, long p)
1246 {
1247   long s,n;
1248   for (s=0,n=0; n-s<k; s += u_lval(++n, p));
1249   return n;
1250 }
1251 
1252 /*
1253  * Using Dwork's expansion, compute \Gamma(px+1)=-\Gamma(px) with x a
1254  * unit.
1255  * See p$-Adic Gamma Functions and Dwork Cohomology,
1256  * Maurizio Boyarsky
1257  * Transactions of the American Mathematical Society,
1258  * Vol. 257, No. 2. (Feb., 1980), pp. 359-369.
1259  * Inspired by a GP script by Fernando Rodriguez-Villegas
1260  */
1261 
1262 static GEN
gadw(GEN x,long p)1263 gadw(GEN x, long p)
1264 {
1265   pari_sp ltop=avma;
1266   GEN s, t;
1267   long j, k;
1268   long n = nboft(precp(x)+valp(x)+1,p);
1269   GEN  u = cgetg(p+1, t_VEC);
1270   s = gaddsg(1, zeropadic(gel(x,2), n));
1271   t = s;
1272   gel(u, 1) = s;
1273   for (j = 1; j < p; ++j)
1274     gel(u, j + 1) = gdivgs(gel(u, j), j);
1275   for (k = 1; k < n; ++k)
1276   {
1277     gel(u, 1) = gdivgs(gdivgs(gadd(gel(u, 1), gel(u, p)), k), p);
1278     for (j = 1; j < p; ++j)
1279       gel(u, j + 1) = gdivgs(gadd(gel(u, j), gel(u, j + 1)), (k*p) + j);
1280 
1281     t = gmul(t, gaddgs(x, k-1));
1282     s = gadd(s, gmul(gmul(gel(u, 1), gpowgs(gel(x,2), k)), t));
1283     if ((k&0xFL)==0) gerepileall(ltop, 3, &u,&s,&t);
1284   }
1285 
1286   return gneg(s);
1287 }
1288 
1289 /*Use Dwork expansion*/
1290 /*This is a O(p*e*log(pe)) algorithm, should be used when p small
1291  * If p==2 this is a O(pe) algorithm...
1292  */
1293 static GEN
gammap_Dwork(GEN x,long p)1294 gammap_Dwork(GEN x, long p)
1295 {
1296   pari_sp ltop = avma;
1297   long k = itos(gmodgs(x, p));
1298   GEN p1;
1299   long j;
1300   if (k)
1301   {
1302     x = gdivgs(gsubgs(x, k), p);
1303     k--;
1304     p1 = gadw(x, p);
1305     if (k%2==1) p1 = gneg(p1);
1306     for (j = 1; j <= k; ++j)
1307       p1 = gmul(p1, gaddgs(gmulsg(p, x), j));
1308   }
1309   else
1310     p1 = gneg(gadw(gdivgs(x, p), p));
1311   return gerepileupto(ltop, p1);
1312 }
1313 
1314 /*
1315  * Compute gammap using the definition. This is a O(x*M(log(pe))) algorithm.
1316  * This should be used if x is very small.
1317  */
1318 static GEN
gammap_Morita(long n,GEN p,long e)1319 gammap_Morita(long n, GEN p, long e)
1320 {
1321   pari_sp ltop=avma;
1322   GEN p2 = gaddsg((n&1)?-1:1, zeropadic(p, e+1));
1323   long i;
1324   long pp=is_bigint(p)? 0: itos(p);
1325   for (i = 2; i < n; i++)
1326     if (!pp || i%pp)
1327     {
1328       p2 = gmulgs(p2, i);
1329       if ((i&0xFL) == 0xFL)
1330         p2 = gerepileupto(ltop, p2);
1331     }
1332   return gerepileupto(ltop, p2);
1333 }
1334 
1335 /*
1336  * x\in\N: Gamma(-x)=(-1)^(1+x+x\p)*Gamma(1+x)
1337  */
1338 
1339 static GEN
gammap_neg_Morita(long n,GEN p,long e)1340 gammap_neg_Morita(long n, GEN p, long e)
1341 {
1342   GEN g = ginv(gammap_Morita(n+1,p,e));
1343   return ((n^sdivsi(n,p)) & 1)? g :gneg(g);
1344 }
1345 
1346 /* p-adic Gamma function for x a p-adic integer */
1347 /*
1348   There are three cases:
1349   n is small            : we use Morita definition.
1350   n is large, p is small: we use Dwork expansion.
1351   n is large, p is large: we don't know how to proceed.
1352   TODO: handle p=2 better (gammap_Dwork is very slow for p=2).
1353 */
1354 #define GAMMAP_DWORK_LIMIT 50000UL
1355 static GEN
gammap(GEN x)1356 gammap(GEN x)
1357 {
1358   GEN p = gel(x,2);
1359   long e= precp(x);
1360   GEN n,m,nm;
1361   if (valp(x)<0)
1362     pari_err(talker,"Gamma not defined for non-integral p-adic number");
1363   n = gtrunc(x);
1364   m = gtrunc(gneg(x));
1365   nm= cmpii(n,m)<=0?n:m;
1366   if (lgefint(nm)==3 && (is_bigint(p) || (ulong)nm[2]<GAMMAP_DWORK_LIMIT))
1367   {
1368     if(nm==n)
1369       return gammap_Morita(itos(n),p,e);
1370     else
1371       return gammap_neg_Morita(itos(m),p,e);
1372   }
1373   return gammap_Dwork(x, itos(p));
1374 }
1375 
1376 GEN
ggamma(GEN x,long prec)1377 ggamma(GEN x, long prec)
1378 {
1379   pari_sp av;
1380   long m;
1381   GEN y, z;
1382 
1383   switch(typ(x))
1384   {
1385     case t_INT:
1386       if (signe(x) <= 0) pari_err(talker,"non-positive integer argument in ggamma");
1387       if (cmpiu(x,481177) > 0) pari_err(talker,"argument too large in ggamma");
1388       return mpfactr(itos(x) - 1, prec);
1389 
1390     case t_REAL: case t_COMPLEX:
1391       return cxgamma(x, 0, prec);
1392 
1393     case t_FRAC:
1394       if (!equaliu(gel(x,2),2)) break;
1395       z = gel(x,1); /* true argument is z/2 */
1396       if (is_bigint(z) || labs(m = itos(z)) > 962354)
1397       {
1398         pari_err(talker, "argument too large in ggamma");
1399         return NULL; /* not reached */
1400       }
1401       return gammahs(m-1, prec);
1402 
1403     case t_PADIC: return gammap(x);
1404     case t_INTMOD: pari_err(typeer,"ggamma");
1405     default:
1406       av = avma; if (!(y = toser_i(x))) break;
1407       return gerepileupto(av, gexp(glngamma(y,prec),prec));
1408   }
1409   return transc(ggamma,x,prec);
1410 }
1411 
1412 GEN
mpfactr(long n,long prec)1413 mpfactr(long n, long prec)
1414 {
1415   GEN f = cgetr(prec);
1416   pari_sp av = avma;
1417 
1418   if (n+1 > 350 + 70*(prec-2)) /* heuristic */
1419     affrr(cxgamma(stor(n+1, prec), 0, prec), f);
1420   else
1421     affir(mpfact(n), f);
1422   avma = av; return f;
1423 }
1424 
1425 GEN
glngamma(GEN x,long prec)1426 glngamma(GEN x, long prec)
1427 {
1428   long i, n;
1429   pari_sp av;
1430   GEN a, y, p1;
1431 
1432   switch(typ(x))
1433   {
1434     case t_INT:
1435       if (signe(x) <= 0) pari_err(talker,"non-positive integer in glngamma");
1436       if (cmpiu(x,200 + 50*(prec-2)) > 0) /* heuristic */
1437 	return cxgamma(x, 1, prec);
1438       av = avma;
1439       return gerepileuptoleaf(av, logr_abs( itor(mpfact(itos(x) - 1), prec) ));
1440 
1441     case t_REAL: case t_COMPLEX:
1442       return cxgamma(x, 1, prec);
1443 
1444     default:
1445       av = avma; if (!(y = toser_i(x))) break;
1446       if (valp(y)) pari_err(negexper,"glngamma");
1447       p1 = gsubsg(1,y);
1448       if (!valp(p1)) pari_err(impl,"lngamma around a!=1");
1449       n = (lg(y)-3) / valp(p1);
1450       a = zeroser(varn(y), lg(y)-2);
1451       for (i=n; i>=2; i--) a = gmul(p1, gadd(a, gdivgs(szeta(i, prec),i)));
1452       a = gadd(a, mpeuler(prec));
1453       return gerepileupto(av, gmul(a, p1));
1454 
1455     case t_PADIC:  pari_err(impl,"p-adic lngamma function");
1456     case t_INTMOD: pari_err(typeer,"glngamma");
1457   }
1458   return transc(glngamma,x,prec);
1459 }
1460 /********************************************************************/
1461 /**                                                                **/
1462 /**                  PSI(x) = GAMMA'(x)/GAMMA(x)                   **/
1463 /**                                                                **/
1464 /********************************************************************/
1465 
1466 GEN
cxpsi(GEN s0,long prec)1467 cxpsi(GEN s0, long prec)
1468 {
1469   pari_sp av, av2;
1470   GEN sum,z,a,res,tes,in2,sig,s,unr;
1471   long lim,nn,k;
1472   const long la = 3;
1473   int funeq = 0;
1474 
1475   if (DEBUGLEVEL>2) (void)timer2();
1476   s = trans_fix_arg(&prec,&s0,&sig,&av,&res);
1477   if (signe(sig) <= 0) { funeq = 1; s = gsub(gen_1, s); sig = real_i(s); }
1478   if (typ(s0) == t_INT && signe(s0) <= 0)
1479     pari_err(talker,"non-positive integer argument in cxpsi");
1480 
1481   {
1482     double ssig = rtodbl(sig);
1483     double st = rtodbl(imag_i(s));
1484     double l;
1485     {
1486       double rlog, ilog; /* log (s - Euler) */
1487       dcxlog(ssig - 0.57721566, st, &rlog,&ilog);
1488       l = dnorm(rlog,ilog);
1489     }
1490     if (l < 0.000001) l = 0.000001;
1491     l = log(l) / 2.;
1492     lim = 2 + (long)ceil((bit_accuracy_mul(prec, LOG2) - l) / (2*(1+log((double)la))));
1493     if (lim < 2) lim = 2;
1494 
1495     l = (2*lim-1)*la / (2.*PI);
1496     l = l*l - st*st;
1497     if (l < 0.) l = 0.;
1498     nn = (long)ceil( sqrt(l) - ssig );
1499     if (nn < 1) nn = 1;
1500     if (DEBUGLEVEL>2) fprintferr("lim, nn: [%ld, %ld]\n",lim,nn);
1501   }
1502   prec++; unr = real_1(prec); /* one extra word of precision */
1503 
1504   a = gdiv(unr, gaddgs(s, nn)); /* 1 / (s+n) */
1505   av2 = avma; sum = gmul2n(a,-1);
1506   for (k = 0; k < nn; k++)
1507   {
1508     sum = gadd(sum, gdiv(unr, gaddgs(s, k)));
1509     if ((k & 127) == 0) sum = gerepileupto(av2, sum);
1510   }
1511   z = gsub(glog(gaddgs(s, nn), prec), sum);
1512   if (DEBUGLEVEL>2) msgtimer("sum from 0 to N-1");
1513 
1514   in2 = gsqr(a);
1515   av2 = avma; tes = divrs(bernreal(2*lim, prec), 2*lim);
1516   for (k=2*lim-2; k>=2; k-=2)
1517   {
1518     tes = gadd(gmul(in2,tes), divrs(bernreal(k, prec), k));
1519     if ((k & 255) == 0) tes = gerepileupto(av2, tes);
1520   }
1521   if (DEBUGLEVEL>2) msgtimer("Bernoulli sum");
1522   z = gsub(z, gmul(in2,tes));
1523   if (funeq)
1524   {
1525     GEN pi = mppi(prec);
1526     z = gadd(z, gmul(pi, gcotan(gmul(pi,s), prec)));
1527   }
1528   if (typ(z) == t_REAL) affr_fixlg(z, res);
1529   else
1530   {
1531     affr_fixlg(gel(z,1), gel(res,1));
1532     affr_fixlg(gel(z,2), gel(res,2));
1533   }
1534   avma = av; return res;
1535 }
1536 
1537 GEN
gpsi(GEN x,long prec)1538 gpsi(GEN x, long prec)
1539 {
1540   switch(typ(x))
1541   {
1542     case t_REAL: case t_COMPLEX: return cxpsi(x,prec);
1543     case t_INTMOD: case t_PADIC: pari_err(typeer,"gpsi");
1544     case t_SER: pari_err(impl,"psi of power series");
1545   }
1546   return transc(gpsi,x,prec);
1547 }
1548