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 to exceed 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 - trans1.c (source / functions) Hit Total Coverage
Test: PARI/GP v2.12.0 lcov report (development 23017-8c5e72c46) Lines: 2034 2096 97.0 %
Date: 2018-09-23 05:39:13 Functions: 151 153 98.7 %
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             : /**                                                                **/
      18             : /********************************************************************/
      19             : #include "pari.h"
      20             : #include "paripriv.h"
      21             : 
      22             : #ifdef LONG_IS_64BIT
      23             : static const long SQRTVERYBIGINT = 3037000500L; /* ceil(sqrt(LONG_MAX)) */
      24             : #else
      25             : static const long SQRTVERYBIGINT = 46341L;
      26             : #endif
      27             : 
      28             : static THREAD GEN gcatalan, geuler, glog2, gpi;
      29             : void
      30      113316 : pari_init_floats(void)
      31             : {
      32      113316 :   gcatalan = geuler = gpi = bernzone = glog2 = NULL;
      33      113316 : }
      34             : 
      35             : void
      36      112062 : pari_close_floats(void)
      37             : {
      38      112062 :   if (gcatalan) gunclone(gcatalan);
      39      112108 :   if (geuler) gunclone(geuler);
      40      112108 :   if (gpi) gunclone(gpi);
      41      112108 :   if (bernzone) gunclone(bernzone);
      42      112108 :   if (glog2) gunclone(glog2);
      43      112108 : }
      44             : 
      45             : /********************************************************************/
      46             : /**                   GENERIC BINARY SPLITTING                     **/
      47             : /**                    (Haible, Papanikolaou)                      **/
      48             : /********************************************************************/
      49             : void
      50        7102 : abpq_init(struct abpq *A, long n)
      51             : {
      52        7102 :   A->a = (GEN*)new_chunk(n+1);
      53        7102 :   A->b = (GEN*)new_chunk(n+1);
      54        7102 :   A->p = (GEN*)new_chunk(n+1);
      55        7102 :   A->q = (GEN*)new_chunk(n+1);
      56        7102 : }
      57             : static GEN
      58      653279 : mulii3(GEN a, GEN b, GEN c) { return mulii(mulii(a,b),c); }
      59             : static GEN
      60      193730 : mulii4(GEN a, GEN b, GEN c, GEN d) { return mulii(mulii(a,b),mulii(c,d)); }
      61             : 
      62             : /* T_{n1,n1+1}, given P = p[n1]p[n1+1] */
      63             : static GEN
      64      193730 : T2(struct abpq *A, long n1, GEN P)
      65             : {
      66      193730 :   GEN u1 = mulii4(A->a[n1], A->p[n1], A->b[n1+1], A->q[n1+1]);
      67      193730 :   GEN u2 = mulii3(A->b[n1],A->a[n1+1], P);
      68      193730 :   return addii(u1, u2);
      69             : }
      70             : 
      71             : /* assume n2 > n1. Compute sum_{n1 <= n < n2} a/b(n) p/q(n1)... p/q(n) */
      72             : void
      73      380444 : abpq_sum(struct abpq_res *r, long n1, long n2, struct abpq *A)
      74             : {
      75             :   struct abpq_res L, R;
      76             :   GEN u1, u2;
      77             :   pari_sp av;
      78             :   long n;
      79      380444 :   switch(n2 - n1)
      80             :   {
      81             :     GEN b, p, q;
      82             :     case 1:
      83          43 :       r->P = A->p[n1];
      84          43 :       r->Q = A->q[n1];
      85          43 :       r->B = A->b[n1];
      86          43 :       r->T = mulii(A->a[n1], A->p[n1]);
      87      193816 :       return;
      88             :     case 2:
      89      107523 :       r->P = mulii(A->p[n1], A->p[n1+1]);
      90      107523 :       r->Q = mulii(A->q[n1], A->q[n1+1]);
      91      107523 :       r->B = mulii(A->b[n1], A->b[n1+1]);
      92      107523 :       av = avma;
      93      107523 :       r->T = gerepileuptoint(av, T2(A, n1, r->P));
      94      107523 :       return;
      95             : 
      96             :     case 3:
      97       86207 :       p = mulii(A->p[n1+1], A->p[n1+2]);
      98       86207 :       q = mulii(A->q[n1+1], A->q[n1+2]);
      99       86207 :       b = mulii(A->b[n1+1], A->b[n1+2]);
     100       86207 :       r->P = mulii(A->p[n1], p);
     101       86207 :       r->Q = mulii(A->q[n1], q);
     102       86207 :       r->B = mulii(A->b[n1], b);
     103       86207 :       av = avma;
     104       86207 :       u1 = mulii3(b, q, A->a[n1]);
     105       86207 :       u2 = mulii(A->b[n1], T2(A, n1+1, p));
     106       86207 :       r->T = gerepileuptoint(av, mulii(A->p[n1], addii(u1, u2)));
     107       86207 :       return;
     108             :   }
     109             : 
     110      186671 :   av = avma;
     111      186671 :   n = (n1 + n2) >> 1;
     112      186671 :   abpq_sum(&L, n1, n, A);
     113      186671 :   abpq_sum(&R, n, n2, A);
     114             : 
     115      186671 :   r->P = mulii(L.P, R.P);
     116      186671 :   r->Q = mulii(L.Q, R.Q);
     117      186671 :   r->B = mulii(L.B, R.B);
     118      186671 :   u1 = mulii3(R.B,R.Q,L.T);
     119      186671 :   u2 = mulii3(L.B,L.P,R.T);
     120      186671 :   r->T = addii(u1,u2);
     121      186671 :   set_avma(av);
     122      186671 :   r->P = icopy(r->P);
     123      186671 :   r->Q = icopy(r->Q);
     124      186671 :   r->B = icopy(r->B);
     125      186671 :   r->T = icopy(r->T);
     126             : }
     127             : 
     128             : /********************************************************************/
     129             : /**                                                                **/
     130             : /**                               PI                               **/
     131             : /**                                                                **/
     132             : /********************************************************************/
     133             : /* replace *old clone by c. Protect against SIGINT */
     134             : static void
     135        5010 : swap_clone(GEN *old, GEN c)
     136             : {
     137        5010 :   GEN tmp = *old;
     138        5010 :   *old = c;
     139        5010 :   if (tmp) gunclone(tmp);
     140        5010 : }
     141             : 
     142             : /*                         ----
     143             :  *  53360 (640320)^(1/2)   \    (6n)! (545140134 n + 13591409)
     144             :  *  -------------------- = /    ------------------------------
     145             :  *        Pi               ----   (n!)^3 (3n)! (-640320)^(3n)
     146             :  *                         n>=0
     147             :  *
     148             :  * Ramanujan's formula + binary splitting */
     149             : static GEN
     150        2096 : pi_ramanujan(long prec)
     151             : {
     152        2096 :   const ulong B = 545140134, A = 13591409, C = 640320;
     153        2096 :   const double alpha2 = 47.11041314; /* 3log(C/12) / log(2) */
     154             :   long n, nmax, prec2;
     155             :   struct abpq_res R;
     156             :   struct abpq S;
     157             :   GEN D, u;
     158             : 
     159        2096 :   nmax = (long)(1 + prec2nbits(prec)/alpha2);
     160             : #ifdef LONG_IS_64BIT
     161        1792 :   D = utoipos(10939058860032000UL); /* C^3/24 */
     162             : #else
     163         304 :   D = uutoi(2546948UL,495419392UL);
     164             : #endif
     165        2096 :   abpq_init(&S, nmax);
     166        2096 :   S.a[0] = utoipos(A);
     167        2096 :   S.b[0] = S.p[0] = S.q[0] = gen_1;
     168       42052 :   for (n = 1; n <= nmax; n++)
     169             :   {
     170       39956 :     S.a[n] = addiu(muluu(B, n), A);
     171       39956 :     S.b[n] = gen_1;
     172       39956 :     S.p[n] = mulis(muluu(6*n-5, 2*n-1), 1-6*n);
     173       39956 :     S.q[n] = mulii(sqru(n), muliu(D,n));
     174             :   }
     175        2096 :   abpq_sum(&R, 0, nmax, &S); prec2 = prec+EXTRAPRECWORD;
     176        2096 :   u = itor(muliu(R.Q,C/12), prec2);
     177        2096 :   return rtor(mulrr(divri(u, R.T), sqrtr_abs(utor(C,prec2))), prec);
     178             : }
     179             : 
     180             : #if 0 /* Much slower than binary splitting at least up to prec = 10^8 */
     181             : /* Gauss - Brent-Salamin AGM iteration */
     182             : static GEN
     183             : pi_brent_salamin(long prec)
     184             : {
     185             :   GEN A, B, C;
     186             :   pari_sp av2;
     187             :   long i, G;
     188             : 
     189             :   G = - prec2nbits(prec);
     190             :   incrprec(prec);
     191             : 
     192             :   A = real2n(-1, prec);
     193             :   B = sqrtr_abs(A); /* = 1/sqrt(2) */
     194             :   setexpo(A, 0);
     195             :   C = real2n(-2, prec); av2 = avma;
     196             :   for (i = 0;; i++)
     197             :   {
     198             :     GEN y, a, b, B_A = subrr(B, A);
     199             :     pari_sp av3 = avma;
     200             :     if (expo(B_A) < G) break;
     201             :     a = addrr(A,B); shiftr_inplace(a, -1);
     202             :     b = mulrr(A,B);
     203             :     affrr(a, A);
     204             :     affrr(sqrtr_abs(b), B); set_avma(av3);
     205             :     y = sqrr(B_A); shiftr_inplace(y, i - 2);
     206             :     affrr(subrr(C, y), C); set_avma(av2);
     207             :   }
     208             :   shiftr_inplace(C, 2);
     209             :   return divrr(sqrr(addrr(A,B)), C);
     210             : }
     211             : #endif
     212             : 
     213             : GEN
     214    14410280 : constpi(long prec)
     215             : {
     216             :   pari_sp av;
     217             :   GEN tmp;
     218    14410280 :   if (gpi && realprec(gpi) >= prec) return gpi;
     219             : 
     220        2096 :   av = avma;
     221        2096 :   tmp = gclone(pi_ramanujan(prec));
     222        2096 :   swap_clone(&gpi,tmp);
     223        2096 :   set_avma(av); return gpi;
     224             : }
     225             : 
     226             : GEN
     227    14409589 : mppi(long prec) { return rtor(constpi(prec), prec); }
     228             : 
     229             : /* Pi * 2^n */
     230             : GEN
     231     5034459 : Pi2n(long n, long prec)
     232             : {
     233     5034459 :   GEN x = mppi(prec); shiftr_inplace(x, n);
     234     5034459 :   return x;
     235             : }
     236             : 
     237             : /* I * Pi * 2^n */
     238             : GEN
     239        8232 : PiI2n(long n, long prec) { retmkcomplex(gen_0, Pi2n(n, prec)); }
     240             : 
     241             : /* 2I * Pi */
     242             : GEN
     243        4361 : PiI2(long prec) { return PiI2n(1, prec); }
     244             : 
     245             : /********************************************************************/
     246             : /**                                                                **/
     247             : /**                       EULER CONSTANT                           **/
     248             : /**                                                                **/
     249             : /********************************************************************/
     250             : 
     251             : GEN
     252       31294 : consteuler(long prec)
     253             : {
     254             :   GEN u,v,a,b,tmpeuler;
     255             :   long l, n1, n, k, x;
     256             :   pari_sp av1, av2;
     257             : 
     258       31294 :   if (geuler && realprec(geuler) >= prec) return geuler;
     259             : 
     260         404 :   av1 = avma; tmpeuler = cgetr_block(prec);
     261             : 
     262         404 :   incrprec(prec);
     263             : 
     264         404 :   l = prec+EXTRAPRECWORD; x = (long) (1 + prec2nbits_mul(l, M_LN2/4));
     265         404 :   a = utor(x,l); u=logr_abs(a); setsigne(u,-1); affrr(u,a);
     266         404 :   b = real_1(l);
     267         404 :   v = real_1(l);
     268         404 :   n = (long)(1+3.591*x); /* z=3.591: z*[ ln(z)-1 ]=1 */
     269         404 :   n1 = minss(n, SQRTVERYBIGINT);
     270         404 :   if (x < SQRTVERYBIGINT)
     271             :   {
     272         404 :     ulong xx = x*x;
     273         404 :     av2 = avma;
     274      138045 :     for (k=1; k<n1; k++)
     275             :     {
     276      137641 :       affrr(divru(mulur(xx,b),k*k), b);
     277      137641 :       affrr(divru(addrr(divru(mulur(xx,a),k),b),k), a);
     278      137641 :       affrr(addrr(u,a), u);
     279      137641 :       affrr(addrr(v,b), v); set_avma(av2);
     280             :     }
     281         808 :     for (   ; k<=n; k++)
     282             :     {
     283         404 :       affrr(divru(divru(mulur(xx,b),k),k), b);
     284         404 :       affrr(divru(addrr(divru(mulur(xx,a),k),b),k), a);
     285         404 :       affrr(addrr(u,a), u);
     286         404 :       affrr(addrr(v,b), v); set_avma(av2);
     287             :     }
     288             :   }
     289             :   else
     290             :   {
     291           0 :     GEN xx = sqru(x);
     292           0 :     av2 = avma;
     293           0 :     for (k=1; k<n1; k++)
     294             :     {
     295           0 :       affrr(divru(mulir(xx,b),k*k), b);
     296           0 :       affrr(divru(addrr(divru(mulir(xx,a),k),b),k), a);
     297           0 :       affrr(addrr(u,a), u);
     298           0 :       affrr(addrr(v,b), v); set_avma(av2);
     299             :     }
     300           0 :     for (   ; k<=n; k++)
     301             :     {
     302           0 :       affrr(divru(divru(mulir(xx,b),k),k), b);
     303           0 :       affrr(divru(addrr(divru(mulir(xx,a),k),b),k), a);
     304           0 :       affrr(addrr(u,a), u);
     305           0 :       affrr(addrr(v,b), v); set_avma(av2);
     306             :     }
     307             :   }
     308         404 :   divrrz(u,v,tmpeuler);
     309         404 :   swap_clone(&geuler,tmpeuler);
     310         404 :   set_avma(av1); return geuler;
     311             : }
     312             : 
     313             : GEN
     314       31294 : mpeuler(long prec) { return rtor(consteuler(prec), prec); }
     315             : 
     316             : /********************************************************************/
     317             : /**                                                                **/
     318             : /**                       CATALAN CONSTANT                         **/
     319             : /**                                                                **/
     320             : /********************************************************************/
     321             : /* 8G = 3\sum_{n>=0} 1/(binomial(2n,n)(2n+1)^2) + Pi log(2+sqrt(3)) */
     322             : static GEN
     323          14 : catalan(long prec)
     324             : {
     325          14 :   long i, nmax = prec2nbits(prec) >> 1;
     326             :   struct abpq_res R;
     327             :   struct abpq A;
     328             :   GEN u, v;
     329          14 :   abpq_init(&A, nmax);
     330          14 :   A.a[0] = A.b[0] = A.p[0] = A.q[0] = gen_1;
     331        6510 :   for (i = 1; i <= nmax; i++)
     332             :   {
     333        6496 :     A.a[i] = gen_1;
     334        6496 :     A.b[i] = utoipos((i<<1)+1);
     335        6496 :     A.p[i] = utoipos(i);
     336        6496 :     A.q[i] = utoipos((i<<2)+2);
     337             :   }
     338          14 :   abpq_sum(&R, 0, nmax, &A);
     339          14 :   u = mulur(3, rdivii(R.T, mulii(R.B,R.Q),prec));
     340          14 :   v = mulrr(mppi(prec), logr_abs(addrs(sqrtr_abs(utor(3,prec)), 2)));
     341          14 :   u = addrr(u, v); shiftr_inplace(u, -3);
     342          14 :   return u;
     343             : }
     344             : 
     345             : GEN
     346          14 : constcatalan(long prec)
     347             : {
     348          14 :   pari_sp av = avma;
     349             :   GEN tmp;
     350          14 :   if (gcatalan && realprec(gcatalan) >= prec) return gcatalan;
     351          14 :   tmp = gclone(catalan(prec));
     352          14 :   swap_clone(&gcatalan,tmp);
     353          14 :   set_avma(av); return gcatalan;
     354             : }
     355             : 
     356             : GEN
     357          14 : mpcatalan(long prec) { return rtor(constcatalan(prec), prec); }
     358             : 
     359             : /********************************************************************/
     360             : /**                                                                **/
     361             : /**          TYPE CONVERSION FOR TRANSCENDENTAL FUNCTIONS          **/
     362             : /**                                                                **/
     363             : /********************************************************************/
     364             : static GEN
     365       73718 : transvec(GEN (*f)(GEN,long), GEN x, long prec)
     366             : {
     367             :   long i, l;
     368       73718 :   GEN y = cgetg_copy(x, &l);
     369       73718 :   for (i=1; i<l; i++) gel(y,i) = f(gel(x,i),prec);
     370       73711 :   return y;
     371             : }
     372             : 
     373             : GEN
     374      629045 : trans_eval(const char *fun, GEN (*f)(GEN,long), GEN x, long prec)
     375             : {
     376      629045 :   pari_sp av = avma;
     377      629045 :   if (prec < 3) pari_err_BUG("trans_eval [prec < 3]");
     378      629045 :   switch(typ(x))
     379             :   {
     380      470408 :     case t_INT:    x = f(itor(x,prec),prec); break;
     381       84863 :     case t_FRAC:   x = f(fractor(x, prec),prec); break;
     382           7 :     case t_QUAD:   x = f(quadtofp(x,prec),prec); break;
     383          14 :     case t_POLMOD: x = transvec(f, polmod_to_embed(x,prec), prec); break;
     384             :     case t_VEC:
     385             :     case t_COL:
     386       73704 :     case t_MAT: return transvec(f, x, prec);
     387          49 :     default: pari_err_TYPE(fun,x); return NULL;
     388             :   }
     389      555271 :   return gerepileupto(av, x);
     390             : }
     391             : 
     392             : /*******************************************************************/
     393             : /*                                                                 */
     394             : /*                            POWERING                             */
     395             : /*                                                                 */
     396             : /*******************************************************************/
     397             : /* x a t_REAL 0, return exp(x) */
     398             : static GEN
     399      207843 : mpexp0(GEN x)
     400             : {
     401      207843 :   long e = expo(x);
     402      207843 :   return e >= 0? real_0_bit(e): real_1_bit(-e);
     403             : }
     404             : static GEN
     405        2975 : powr0(GEN x)
     406        2975 : { return signe(x)? real_1(realprec(x)): mpexp0(x); }
     407             : 
     408             : /* x t_POL or t_SER, return scalarpol(Rg_get_1(x)) */
     409             : static GEN
     410      176073 : scalarpol_get_1(GEN x)
     411             : {
     412      176073 :   GEN y = cgetg(3,t_POL);
     413      176073 :   y[1] = evalvarn(varn(x)) | evalsigne(1);
     414      176073 :   gel(y,2) = Rg_get_1(x); return y;
     415             : }
     416             : /* to be called by the generic function gpowgs(x,s) when s = 0 */
     417             : static GEN
     418     3698920 : gpowg0(GEN x)
     419             : {
     420             :   long lx, i;
     421             :   GEN y;
     422             : 
     423     3698920 :   switch(typ(x))
     424             :   {
     425             :     case t_INT: case t_REAL: case t_FRAC: case t_PADIC:
     426      812009 :       return gen_1;
     427             : 
     428           7 :     case t_QUAD: x++; /*fall through*/
     429             :     case t_COMPLEX: {
     430       21270 :       pari_sp av = avma;
     431       21270 :       GEN a = gpowg0(gel(x,1));
     432       21270 :       GEN b = gpowg0(gel(x,2));
     433       21270 :       if (a == gen_1) return b;
     434          14 :       if (b == gen_1) return a;
     435           7 :       return gerepileupto(av, gmul(a,b));
     436             :     }
     437             :     case t_INTMOD:
     438         119 :       y = cgetg(3,t_INTMOD);
     439         119 :       gel(y,1) = icopy(gel(x,1));
     440         119 :       gel(y,2) = is_pm1(gel(x,1))? gen_0: gen_1;
     441         119 :       return y;
     442             : 
     443        5656 :     case t_FFELT: return FF_1(x);
     444             : 
     445             :     case t_POLMOD:
     446         133 :       retmkpolmod(scalarpol_get_1(gel(x,1)), gcopy(gel(x,1)));
     447             : 
     448             :     case t_RFRAC:
     449           7 :       return scalarpol_get_1(gel(x,2));
     450             :     case t_POL: case t_SER:
     451      175933 :       return scalarpol_get_1(x);
     452             : 
     453             :     case t_MAT:
     454          84 :       lx=lg(x); if (lx==1) return cgetg(1,t_MAT);
     455          77 :       if (lx != lgcols(x)) pari_err_DIM("gpow");
     456          77 :       y = matid(lx-1);
     457          77 :       for (i=1; i<lx; i++) gcoeff(y,i,i) = gpowg0(gcoeff(x,i,i));
     458          77 :       return y;
     459          14 :     case t_QFR: return qfr_1(x);
     460     2683681 :     case t_QFI: return qfi_1(x);
     461           7 :     case t_VECSMALL: return identity_perm(lg(x) - 1);
     462             :   }
     463           7 :   pari_err_TYPE("gpow",x);
     464             :   return NULL; /* LCOV_EXCL_LINE */
     465             : }
     466             : 
     467             : static GEN
     468    42563139 : _sqr(void *data /* ignored */, GEN x) { (void)data; return gsqr(x); }
     469             : static GEN
     470    17335421 : _mul(void *data /* ignored */, GEN x, GEN y) { (void)data; return gmul(x,y); }
     471             : static GEN
     472      110498 : _one(void *x) { return gpowg0((GEN) x); }
     473             : static GEN
     474     5318008 : _sqri(void *data /* ignored */, GEN x) { (void)data; return sqri(x); }
     475             : static GEN
     476     3648996 : _muli(void *data /* ignored */, GEN x, GEN y) { (void)data; return mulii(x,y); }
     477             : static GEN
     478     9480676 : _sqrr(void *data /* ignored */, GEN x) { (void)data; return sqrr(x); }
     479             : static GEN
     480     1543543 : _mulr(void *data /* ignored */, GEN x, GEN y) { (void)data; return mulrr(x,y); }
     481             : static GEN
     482       14505 : _oner(void *data /* prec */) { return real_1( *(long*) data); }
     483             : 
     484             : /* INTEGER POWERING (a^n for integer a != 0 and integer n > 0)
     485             :  *
     486             :  * Use left shift binary algorithm (RS is wasteful: multiplies big numbers,
     487             :  * with LS one of them is the base, hence small). Sign of result is set
     488             :  * to s (= 1,-1). Makes life easier for caller, which otherwise might do a
     489             :  * setsigne(gen_1 / gen_m1) */
     490             : static GEN
     491    37063754 : powiu_sign(GEN a, ulong N, long s)
     492             : {
     493             :   pari_sp av;
     494             :   GEN y;
     495             : 
     496    37063754 :   if (lgefint(a) == 3)
     497             :   { /* easy if |a| < 3 */
     498    33224578 :     ulong q = a[2];
     499    33224578 :     if (q == 1) return (s>0)? gen_1: gen_m1;
     500    24717676 :     if (q == 2) { a = int2u(N); setsigne(a,s); return a; }
     501    14376690 :     q = upowuu(q, N);
     502    14376690 :     if (q) return s>0? utoipos(q): utoineg(q);
     503             :   }
     504     4515214 :   if (N <= 2) {
     505     2896607 :     if (N == 2) return sqri(a);
     506       23839 :     a = icopy(a); setsigne(a,s); return a;
     507             :   }
     508     1618607 :   av = avma;
     509     1618607 :   y = gen_powu_i(a, N, NULL, &_sqri, &_muli);
     510     1618607 :   setsigne(y,s); return gerepileuptoint(av, y);
     511             : }
     512             : /* a^n */
     513             : GEN
     514    36936925 : powiu(GEN a, ulong n)
     515             : {
     516             :   long s;
     517    36936925 :   if (!n) return gen_1;
     518    36522271 :   s = signe(a);
     519    36522271 :   if (!s) return gen_0;
     520    36508844 :   return powiu_sign(a, n, (s < 0 && odd(n))? -1: 1);
     521             : }
     522             : GEN
     523    22576439 : powis(GEN a, long n)
     524             : {
     525             :   long s;
     526             :   GEN t, y;
     527    22576439 :   if (n >= 0) return powiu(a, n);
     528      106396 :   s = signe(a);
     529      106396 :   if (!s) pari_err_INV("powis",gen_0);
     530      106396 :   t = (s < 0 && odd(n))? gen_m1: gen_1;
     531      106396 :   if (is_pm1(a)) return t;
     532             :   /* n < 0, |a| > 1 */
     533      106102 :   y = cgetg(3,t_FRAC);
     534      106102 :   gel(y,1) = t;
     535      106102 :   gel(y,2) = powiu_sign(a, -n, 1); /* force denominator > 0 */
     536      106102 :   return y;
     537             : }
     538             : GEN
     539    25303010 : powuu(ulong p, ulong N)
     540             : {
     541    25303010 :   pari_sp av = avma;
     542    25303010 :   long P[] = {evaltyp(t_INT)|_evallg(3), evalsigne(1)|evallgefint(3),0};
     543             :   ulong pN;
     544             :   GEN y;
     545    25303010 :   if (N <= 2)
     546             :   {
     547    20068708 :     if (N == 2) return sqru(p);
     548    16424832 :     if (N == 1) return utoi(p);
     549     4260003 :     return gen_1;
     550             :   }
     551     5234302 :   if (!p) return gen_0;
     552     5234302 :   pN = upowuu(p, N);
     553     5234302 :   if (pN) return utoipos(pN);
     554      811889 :   if (p == 2) return int2u(N);
     555      807291 :   P[2] = p; av = avma;
     556      807291 :   y = gen_powu_i(P, N, NULL, &_sqri, &_muli);
     557      807291 :   return gerepileuptoint(av, y);
     558             : }
     559             : 
     560             : /* return 0 if overflow */
     561             : static ulong
     562     5319841 : usqru(ulong p) { return p & HIGHMASK? 0: p*p; }
     563             : ulong
     564    25353929 : upowuu(ulong p, ulong k)
     565             : {
     566             : #ifdef LONG_IS_64BIT
     567    21728165 :   const ulong CUTOFF3 = 2642245;
     568    21728165 :   const ulong CUTOFF4 = 65535;
     569    21728165 :   const ulong CUTOFF5 = 7131;
     570    21728165 :   const ulong CUTOFF6 = 1625;
     571    21728165 :   const ulong CUTOFF7 = 565;
     572    21728165 :   const ulong CUTOFF8 = 255;
     573    21728165 :   const ulong CUTOFF9 = 138;
     574    21728165 :   const ulong CUTOFF10 = 84;
     575    21728165 :   const ulong CUTOFF11 = 56;
     576    21728165 :   const ulong CUTOFF12 = 40;
     577    21728165 :   const ulong CUTOFF13 = 30;
     578    21728165 :   const ulong CUTOFF14 = 23;
     579    21728165 :   const ulong CUTOFF15 = 19;
     580    21728165 :   const ulong CUTOFF16 = 15;
     581    21728165 :   const ulong CUTOFF17 = 13;
     582    21728165 :   const ulong CUTOFF18 = 11;
     583    21728165 :   const ulong CUTOFF19 = 10;
     584    21728165 :   const ulong CUTOFF20 =  9;
     585             : #else
     586     3625764 :   const ulong CUTOFF3 = 1625;
     587     3625764 :   const ulong CUTOFF4 =  255;
     588     3625764 :   const ulong CUTOFF5 =   84;
     589     3625764 :   const ulong CUTOFF6 =   40;
     590     3625764 :   const ulong CUTOFF7 =   23;
     591     3625764 :   const ulong CUTOFF8 =   15;
     592     3625764 :   const ulong CUTOFF9 =   11;
     593     3625764 :   const ulong CUTOFF10 =   9;
     594     3625764 :   const ulong CUTOFF11 =   7;
     595     3625764 :   const ulong CUTOFF12 =   6;
     596     3625764 :   const ulong CUTOFF13 =   5;
     597     3625764 :   const ulong CUTOFF14 =   4;
     598     3625764 :   const ulong CUTOFF15 =   4;
     599     3625764 :   const ulong CUTOFF16 =   3;
     600     3625764 :   const ulong CUTOFF17 =   3;
     601     3625764 :   const ulong CUTOFF18 =   3;
     602     3625764 :   const ulong CUTOFF19 =   3;
     603     3625764 :   const ulong CUTOFF20 =   3;
     604             : #endif
     605             : 
     606    25353929 :   if (p <= 2)
     607             :   {
     608     1819324 :     if (p < 2) return p;
     609     1320324 :     return k < BITS_IN_LONG? 1UL<<k: 0;
     610             :   }
     611    23534605 :   switch(k)
     612             :   {
     613             :     ulong p2, p3, p4, p5, p8;
     614     1828283 :     case 0:  return 1;
     615     8445774 :     case 1:  return p;
     616     5319841 :     case 2:  return usqru(p);
     617     2559917 :     case 3:  if (p > CUTOFF3) return 0; return p*p*p;
     618     1052069 :     case 4:  if (p > CUTOFF4) return 0; p2=p*p; return p2*p2;
     619     1222863 :     case 5:  if (p > CUTOFF5) return 0; p2=p*p; return p2*p2*p;
     620      710468 :     case 6:  if (p > CUTOFF6) return 0; p2=p*p; return p2*p2*p2;
     621      138485 :     case 7:  if (p > CUTOFF7) return 0; p2=p*p; return p2*p2*p2*p;
     622      130746 :     case 8:  if (p > CUTOFF8) return 0; p2=p*p; p4=p2*p2; return p4*p4;
     623      232681 :     case 9:  if (p > CUTOFF9) return 0; p2=p*p; p4=p2*p2; return p4*p4*p;
     624       73673 :     case 10: if (p > CUTOFF10)return 0; p2=p*p; p4=p2*p2; return p4*p4*p2;
     625       38241 :     case 11: if (p > CUTOFF11)return 0; p2=p*p; p4=p2*p2; return p4*p4*p2*p;
     626       56293 :     case 12: if (p > CUTOFF12)return 0; p2=p*p; p4=p2*p2; return p4*p4*p4;
     627       64464 :     case 13: if (p > CUTOFF13)return 0; p2=p*p; p4=p2*p2; return p4*p4*p4*p;
     628       69600 :     case 14: if (p > CUTOFF14)return 0; p2=p*p; p4=p2*p2; return p4*p4*p4*p2;
     629       88021 :     case 15: if (p > CUTOFF15)return 0;
     630       20218 :       p2=p*p; p3=p2*p; p5=p3*p2; return p5*p5*p5;
     631       61516 :     case 16: if (p > CUTOFF16)return 0;
     632       15979 :       p2=p*p; p4=p2*p2; p8=p4*p4; return p8*p8;
     633       35440 :     case 17: if (p > CUTOFF17)return 0;
     634       10754 :       p2=p*p; p4=p2*p2; p8=p4*p4; return p*p8*p8;
     635       24755 :     case 18: if (p > CUTOFF18)return 0;
     636        7571 :       p2=p*p; p4=p2*p2; p8=p4*p4; return p2*p8*p8;
     637      797053 :     case 19: if (p > CUTOFF19)return 0;
     638      783774 :       p2=p*p; p4=p2*p2; p8=p4*p4; return p*p2*p8*p8;
     639       13138 :     case 20: if (p > CUTOFF20)return 0;
     640        4193 :       p2=p*p; p4=p2*p2; p8=p4*p4; return p4*p8*p8;
     641             :   }
     642             : #ifdef LONG_IS_64BIT
     643      489088 :   switch(p)
     644             :   {
     645       54942 :     case 3: if (k > 40) return 0;
     646       23628 :       break;
     647       16206 :     case 4: if (k > 31) return 0;
     648         222 :       return 1UL<<(2*k);
     649       26812 :     case 5: if (k > 27) return 0;
     650        3880 :       break;
     651       13206 :     case 6: if (k > 24) return 0;
     652          54 :       break;
     653       17478 :     case 7: if (k > 22) return 0;
     654        1074 :       break;
     655      360444 :     default: return 0;
     656             :   }
     657             :   /* no overflow */
     658             :   {
     659       28636 :     ulong q = upowuu(p, k >> 1);
     660       28636 :     q *= q ;
     661       28636 :     return odd(k)? q*p: q;
     662             :   }
     663             : #else
     664       82196 :   return 0;
     665             : #endif
     666             : }
     667             : 
     668             : typedef struct {
     669             :   long prec, a;
     670             :   GEN (*sqr)(GEN);
     671             :   GEN (*mulug)(ulong,GEN);
     672             : } sr_muldata;
     673             : 
     674             : static GEN
     675      443422 : _rpowuu_sqr(void *data, GEN x)
     676             : {
     677      443422 :   sr_muldata *D = (sr_muldata *)data;
     678      443422 :   if (typ(x) == t_INT && lgefint(x) >= D->prec)
     679             :   { /* switch to t_REAL */
     680       15500 :     D->sqr   = &sqrr;
     681       15500 :     D->mulug = &mulur; x = itor(x, D->prec);
     682             :   }
     683      443422 :   return D->sqr(x);
     684             : }
     685             : 
     686             : static GEN
     687      164754 : _rpowuu_msqr(void *data, GEN x)
     688             : {
     689      164754 :   GEN x2 = _rpowuu_sqr(data, x);
     690      164754 :   sr_muldata *D = (sr_muldata *)data;
     691      164754 :   return D->mulug(D->a, x2);
     692             : }
     693             : 
     694             : /* return a^n as a t_REAL of precision prec. Assume a > 0, n > 0 */
     695             : GEN
     696      155582 : rpowuu(ulong a, ulong n, long prec)
     697             : {
     698             :   pari_sp av;
     699             :   GEN y, z;
     700             :   sr_muldata D;
     701             : 
     702      155582 :   if (a == 1) return real_1(prec);
     703      155582 :   if (a == 2) return real2n(n, prec);
     704      155582 :   if (n == 1) return utor(a, prec);
     705      153930 :   z = cgetr(prec);
     706      153930 :   av = avma;
     707      153930 :   D.sqr   = &sqri;
     708      153930 :   D.mulug = &mului;
     709      153930 :   D.prec = prec;
     710      153930 :   D.a = (long)a;
     711      153930 :   y = gen_powu_fold_i(utoipos(a), n, (void*)&D, &_rpowuu_sqr, &_rpowuu_msqr);
     712      153930 :   mpaff(y, z); set_avma(av); return z;
     713             : }
     714             : 
     715             : GEN
     716     6934259 : powrs(GEN x, long n)
     717             : {
     718     6934259 :   pari_sp av = avma;
     719             :   GEN y;
     720     6934259 :   if (!n) return powr0(x);
     721     6934259 :   y = gen_powu_i(x, (ulong)labs(n), NULL, &_sqrr, &_mulr);
     722     6934259 :   if (n < 0) y = invr(y);
     723     6934259 :   return gerepileuptoleaf(av,y);
     724             : }
     725             : GEN
     726      751159 : powru(GEN x, ulong n)
     727             : {
     728      751159 :   pari_sp av = avma;
     729             :   GEN y;
     730      751159 :   if (!n) return powr0(x);
     731      748695 :   y = gen_powu_i(x, n, NULL, &_sqrr, &_mulr);
     732      748695 :   return gerepileuptoleaf(av,y);
     733             : }
     734             : 
     735             : GEN
     736       14505 : powersr(GEN x, long n)
     737             : {
     738       14505 :   long prec = realprec(x);
     739       14505 :   return gen_powers(x, n, 1, &prec, &_sqrr, &_mulr, &_oner);
     740             : }
     741             : 
     742             : /* x^(s/2), assume x t_REAL */
     743             : GEN
     744           0 : powrshalf(GEN x, long s)
     745             : {
     746           0 :   if (s & 1) return sqrtr(powrs(x, s));
     747           0 :   return powrs(x, s>>1);
     748             : }
     749             : /* x^(s/2), assume x t_REAL */
     750             : GEN
     751       23640 : powruhalf(GEN x, ulong s)
     752             : {
     753       23640 :   if (s & 1) return sqrtr(powru(x, s));
     754        6321 :   return powru(x, s>>1);
     755             : }
     756             : /* x^(n/d), assume x t_REAL, return t_REAL */
     757             : GEN
     758         511 : powrfrac(GEN x, long n, long d)
     759             : {
     760             :   long z;
     761         511 :   if (!n) return powr0(x);
     762           0 :   z = cgcd(n, d); if (z > 1) { n /= z; d /= z; }
     763           0 :   if (d == 1) return powrs(x, n);
     764           0 :   x = powrs(x, n);
     765           0 :   if (d == 2) return sqrtr(x);
     766           0 :   return sqrtnr(x, d);
     767             : }
     768             : 
     769             : /* assume x != 0 */
     770             : static GEN
     771      565722 : pow_monome(GEN x, long n)
     772             : {
     773      565722 :   long i, d, dx = degpol(x);
     774             :   GEN A, b, y;
     775             : 
     776      565722 :   if (n < 0) { n = -n; y = cgetg(3, t_RFRAC); } else y = NULL;
     777             : 
     778      565722 :   if (HIGHWORD(dx) || HIGHWORD(n))
     779           8 :   {
     780             :     LOCAL_HIREMAINDER;
     781           9 :     d = (long)mulll((ulong)dx, (ulong)n);
     782           9 :     if (hiremainder || (d &~ LGBITS)) d = LGBITS; /* overflow */
     783           9 :     d += 2;
     784             :   }
     785             :   else
     786      565713 :     d = dx*n + 2;
     787      565722 :   if ((d + 1) & ~LGBITS) pari_err(e_OVERFLOW,"pow_monome [degree]");
     788      565715 :   A = cgetg(d+1, t_POL); A[1] = x[1];
     789      565715 :   for (i=2; i < d; i++) gel(A,i) = gen_0;
     790      565715 :   b = gpowgs(gel(x,dx+2), n); /* not memory clean if (n < 0) */
     791      565715 :   if (!y) y = A;
     792             :   else {
     793       20440 :     GEN c = denom_i(b);
     794       20440 :     gel(y,1) = c; if (c != gen_1) b = gmul(b,c);
     795       20440 :     gel(y,2) = A;
     796             :   }
     797      565715 :   gel(A,d) = b; return y;
     798             : }
     799             : 
     800             : /* x t_PADIC */
     801             : static GEN
     802       41307 : powps(GEN x, long n)
     803             : {
     804       41307 :   long e = n*valp(x), v;
     805       41307 :   GEN t, y, mod, p = gel(x,2);
     806             :   pari_sp av;
     807             : 
     808       41307 :   if (!signe(gel(x,4))) {
     809          84 :     if (n < 0) pari_err_INV("powps",x);
     810          77 :     return zeropadic(p, e);
     811             :   }
     812       41223 :   v = z_pval(n, p);
     813             : 
     814       41223 :   y = cgetg(5,t_PADIC);
     815       41223 :   mod = gel(x,3);
     816       41223 :   if (v == 0) mod = icopy(mod);
     817             :   else
     818             :   {
     819       40341 :     if (precp(x) == 1 && absequaliu(p, 2)) v++;
     820       40341 :     mod = mulii(mod, powiu(p,v));
     821       40341 :     mod = gerepileuptoint((pari_sp)y, mod);
     822             :   }
     823       41223 :   y[1] = evalprecp(precp(x) + v) | evalvalp(e);
     824       41223 :   gel(y,2) = icopy(p);
     825       41223 :   gel(y,3) = mod;
     826             : 
     827       41223 :   av = avma; t = gel(x,4);
     828       41223 :   if (n < 0) { t = Fp_inv(t, mod); n = -n; }
     829       41223 :   t = Fp_powu(t, n, mod);
     830       41223 :   gel(y,4) = gerepileuptoint(av, t);
     831       41223 :   return y;
     832             : }
     833             : /* x t_PADIC */
     834             : static GEN
     835         161 : powp(GEN x, GEN n)
     836             : {
     837             :   long v;
     838         161 :   GEN y, mod, p = gel(x,2);
     839             : 
     840         161 :   if (valp(x)) pari_err_OVERFLOW("valp()");
     841             : 
     842         161 :   if (!signe(gel(x,4))) {
     843          14 :     if (signe(n) < 0) pari_err_INV("powp",x);
     844           7 :     return zeropadic(p, 0);
     845             :   }
     846         147 :   v = Z_pval(n, p);
     847             : 
     848         147 :   y = cgetg(5,t_PADIC);
     849         147 :   mod = gel(x,3);
     850         147 :   if (v == 0) mod = icopy(mod);
     851             :   else
     852             :   {
     853          70 :     mod = mulii(mod, powiu(p,v));
     854          70 :     mod = gerepileuptoint((pari_sp)y, mod);
     855             :   }
     856         147 :   y[1] = evalprecp(precp(x) + v) | _evalvalp(0);
     857         147 :   gel(y,2) = icopy(p);
     858         147 :   gel(y,3) = mod;
     859         147 :   gel(y,4) = Fp_pow(gel(x,4), n, mod);
     860         147 :   return y;
     861             : }
     862             : static GEN
     863       29707 : pow_polmod(GEN x, GEN n)
     864             : {
     865       29707 :   GEN z = cgetg(3, t_POLMOD), a = gel(x,2), T = gel(x,1);
     866       29707 :   gel(z,1) = gcopy(T);
     867       29707 :   if (typ(a) != t_POL || varn(a) != varn(T) || lg(a) <= 3)
     868        2751 :     a = powgi(a, n);
     869             :   else {
     870       26956 :     pari_sp av = avma;
     871       26956 :     GEN p = NULL;
     872       26956 :     if (RgX_is_FpX(T, &p) && RgX_is_FpX(a, &p) && p)
     873             :     {
     874        7602 :       T = RgX_to_FpX(T, p); a = RgX_to_FpX(a, p);
     875        7602 :       if (lgefint(p) == 3)
     876             :       {
     877        7595 :         ulong pp = p[2];
     878        7595 :         a = Flxq_pow(ZX_to_Flx(a, pp), n, ZX_to_Flx(T, pp), pp);
     879        7595 :         a = Flx_to_ZX(a);
     880             :       }
     881             :       else
     882           7 :         a = FpXQ_pow(a, n, T, p);
     883        7602 :       a = FpX_to_mod(a, p);
     884        7602 :       a = gerepileupto(av, a);
     885             :     }
     886             :     else
     887             :     {
     888       19354 :       set_avma(av);
     889       19354 :       a = RgXQ_pow(a, n, gel(z,1));
     890             :     }
     891             :   }
     892       29707 :   gel(z,2) = a; return z;
     893             : }
     894             : 
     895             : GEN
     896   127949872 : gpowgs(GEN x, long n)
     897             : {
     898             :   long m;
     899             :   pari_sp av;
     900             :   GEN y;
     901             : 
     902   127949872 :   if (n == 0) return gpowg0(x);
     903   124404165 :   if (n == 1)
     904    72463544 :     switch (typ(x)) {
     905      802271 :       case t_QFI: return redimag(x);
     906          14 :       case t_QFR: return redreal(x);
     907    71661259 :       default: return gcopy(x);
     908             :     }
     909    51940621 :   if (n ==-1) return ginv(x);
     910    46773135 :   switch(typ(x))
     911             :   {
     912    22528854 :     case t_INT: return powis(x,n);
     913     6933349 :     case t_REAL: return powrs(x,n);
     914             :     case t_INTMOD:
     915       21560 :       y = cgetg(3,t_INTMOD); gel(y,1) = icopy(gel(x,1));
     916       21560 :       gel(y,2) = Fp_pows(gel(x,2), n, gel(x,1));
     917       21560 :       return y;
     918             :     case t_FRAC:
     919             :     {
     920      224400 :       GEN a = gel(x,1), b = gel(x,2);
     921      224400 :       long s = (signe(a) < 0 && odd(n))? -1: 1;
     922      224400 :       if (n < 0) {
     923          49 :         n = -n;
     924          49 :         if (is_pm1(a)) return powiu_sign(b, n, s); /* +-1/x[2] inverts to t_INT */
     925          42 :         swap(a, b);
     926             :       }
     927      224393 :       y = cgetg(3, t_FRAC);
     928      224393 :       gel(y,1) = powiu_sign(a, n, s);
     929      224393 :       gel(y,2) = powiu_sign(b, n, 1);
     930      224393 :       return y;
     931             :     }
     932       41307 :     case t_PADIC: return powps(x, n);
     933             :     case t_RFRAC:
     934             :     {
     935      249011 :       av = avma; y = cgetg(3, t_RFRAC); m = labs(n);
     936      249011 :       gel(y,1) = gpowgs(gel(x,1),m);
     937      249011 :       gel(y,2) = gpowgs(gel(x,2),m);
     938      249011 :       if (n < 0) y = ginv(y);
     939      249011 :       return gerepileupto(av,y);
     940             :     }
     941             :     case t_POLMOD: {
     942       29700 :       long N[] = {evaltyp(t_INT) | _evallg(3),0,0};
     943       29700 :       affsi(n,N); return pow_polmod(x, N);
     944             :     }
     945             :     case t_POL:
     946      610823 :       if (RgX_is_monomial(x)) return pow_monome(x, n);
     947             :     default: {
     948    16179232 :       pari_sp av = avma;
     949    16179232 :       y = gen_powu_i(x, (ulong)labs(n), NULL, &_sqr, &_mul);
     950    16179226 :       if (n < 0) y = ginv(y);
     951    16179226 :       return gerepileupto(av,y);
     952             :     }
     953             :   }
     954             : }
     955             : 
     956             : /* n a t_INT */
     957             : GEN
     958   114741232 : powgi(GEN x, GEN n)
     959             : {
     960             :   GEN y;
     961             : 
     962   114741232 :   if (!is_bigint(n)) return gpowgs(x, itos(n));
     963             :   /* probable overflow for non-modular types (typical exception: (X^0)^N) */
     964         341 :   switch(typ(x))
     965             :   {
     966             :     case t_INTMOD:
     967          49 :       y = cgetg(3,t_INTMOD); gel(y,1) = icopy(gel(x,1));
     968          49 :       gel(y,2) = Fp_pow(gel(x,2), n, gel(x,1));
     969          49 :       return y;
     970          59 :     case t_FFELT: return FF_pow(x,n);
     971         161 :     case t_PADIC: return powp(x, n);
     972             : 
     973             :     case t_INT:
     974          35 :       if (is_pm1(x)) return (signe(x) < 0 && mpodd(n))? gen_m1: gen_1;
     975          14 :       if (signe(x)) pari_err_OVERFLOW("lg()");
     976           7 :       if (signe(n) < 0) pari_err_INV("powgi",gen_0);
     977           7 :       return gen_0;
     978             :     case t_FRAC:
     979           7 :       pari_err_OVERFLOW("lg()");
     980             : 
     981          14 :     case t_QFR: return qfrpow(x, n);
     982           7 :     case t_POLMOD: return pow_polmod(x, n);
     983             :     default: {
     984           9 :       pari_sp av = avma;
     985           9 :       y = gen_pow(x, n, NULL, &_sqr, &_mul);
     986           9 :       if (signe(n) < 0) y = ginv(y);
     987           9 :       return gerepileupto(av,y);
     988             :     }
     989             :   }
     990             : }
     991             : 
     992             : /* Assume x = 1 + O(t), n a scalar. Return x^n */
     993             : static GEN
     994        7756 : ser_pow_1(GEN x, GEN n)
     995             : {
     996             :   long lx, mi, i, j, d;
     997        7756 :   GEN y = cgetg_copy(x, &lx), X = x+2, Y = y + 2;
     998        7756 :   y[1] = evalsigne(1) | _evalvalp(0) | evalvarn(varn(x));
     999        7756 :   d = mi = lx-3; while (mi>=1 && isrationalzero(gel(X,mi))) mi--;
    1000        7756 :   gel(Y,0) = gen_1;
    1001      109305 :   for (i=1; i<=d; i++)
    1002             :   {
    1003      101549 :     pari_sp av = avma;
    1004      101549 :     GEN s = gen_0;
    1005      479178 :     for (j=1; j<=minss(i,mi); j++)
    1006             :     {
    1007      377629 :       GEN t = gsubgs(gmulgs(n,j),i-j);
    1008      377629 :       s = gadd(s, gmul(gmul(t, gel(X,j)), gel(Y,i-j)));
    1009             :     }
    1010      101549 :     gel(Y,i) = gerepileupto(av, gdivgs(s,i));
    1011             :   }
    1012        7756 :   return y;
    1013             : }
    1014             : 
    1015             : /* we suppose n != 0, valp(x) = 0 and leading-term(x) != 0. Not stack clean */
    1016             : static GEN
    1017        7875 : ser_pow(GEN x, GEN n, long prec)
    1018             : {
    1019             :   GEN y, c, lead;
    1020        7875 :   if (varncmp(gvar(n), varn(x)) <= 0) return gexp(gmul(n, glog(x,prec)), prec);
    1021        7756 :   lead = gel(x,2);
    1022        7756 :   if (gequal1(lead)) return ser_pow_1(x, n);
    1023        7434 :   x = ser_normalize(x);
    1024        7434 :   if (typ(n) == t_FRAC && !isinexact(lead) && ispower(lead, gel(n,2), &c))
    1025          77 :     c = powgi(c, gel(n,1));
    1026             :   else
    1027        7357 :     c = gpow(lead,n, prec);
    1028        7434 :   y = gmul(c, ser_pow_1(x, n));
    1029             :   /* gpow(t_POLMOD,n) can be a t_COL [conjvec] */
    1030        7434 :   if (typ(y) != t_SER) pari_err_TYPE("gpow", y);
    1031        7434 :   return y;
    1032             : }
    1033             : 
    1034             : static long
    1035        7770 : val_from_i(GEN E)
    1036             : {
    1037        7770 :   if (is_bigint(E)) pari_err_OVERFLOW("sqrtn [valuation]");
    1038        7763 :   return itos(E);
    1039             : }
    1040             : 
    1041             : /* return x^q, assume typ(x) = t_SER, typ(q) = t_INT/t_FRAC and q != 0 */
    1042             : static GEN
    1043        7777 : ser_powfrac(GEN x, GEN q, long prec)
    1044             : {
    1045        7777 :   GEN y, E = gmulsg(valp(x), q);
    1046             :   long e;
    1047             : 
    1048        7777 :   if (!signe(x))
    1049             :   {
    1050          21 :     if (gsigne(q) < 0) pari_err_INV("gpow", x);
    1051          21 :     return zeroser(varn(x), val_from_i(gfloor(E)));
    1052             :   }
    1053        7756 :   if (typ(E) != t_INT)
    1054           7 :     pari_err_DOMAIN("sqrtn", "valuation", "!=", mkintmod(gen_0, gel(q,2)), x);
    1055        7749 :   e = val_from_i(E);
    1056        7749 :   y = leafcopy(x); setvalp(y, 0);
    1057        7749 :   y = ser_pow(y, q, prec);
    1058        7749 :   setvalp(y, e); return y;
    1059             : }
    1060             : 
    1061             : static GEN
    1062         147 : gpow0(GEN x, GEN n, long prec)
    1063             : {
    1064         147 :   pari_sp av = avma;
    1065             :   long i, lx;
    1066             :   GEN y;
    1067         147 :   switch(typ(n))
    1068             :   {
    1069             :     case t_INT: case t_REAL: case t_FRAC: case t_COMPLEX: case t_QUAD:
    1070         105 :       break;
    1071             :     case t_VEC: case t_COL: case t_MAT:
    1072          35 :       y = cgetg_copy(n, &lx);
    1073          35 :       for (i=1; i<lx; i++) gel(y,i) = gpow0(x,gel(n,i),prec);
    1074          35 :       return y;
    1075           7 :     default: pari_err_TYPE("gpow(0,n)", n);
    1076             :   }
    1077         105 :   n = real_i(n);
    1078         105 :   if (gsigne(n) <= 0) pari_err_DOMAIN("gpow(0,n)", "n", "<=", gen_0, n);
    1079          98 :   if (!precision(x)) return gcopy(x);
    1080             : 
    1081          35 :   x = ground(gmulsg(gexpo(x),n));
    1082          35 :   if (is_bigint(x) || uel(x,2) >= HIGHEXPOBIT)
    1083           7 :     pari_err_OVERFLOW("gpow");
    1084          28 :   set_avma(av); return real_0_bit(itos(x));
    1085             : }
    1086             : 
    1087             : GEN
    1088    16621983 : gpow(GEN x, GEN n, long prec)
    1089             : {
    1090    16621983 :   long prec0, i, lx, tx, tn = typ(n);
    1091             :   pari_sp av;
    1092             :   GEN y;
    1093             : 
    1094    16621983 :   if (tn == t_INT) return powgi(x,n);
    1095     4852354 :   tx = typ(x);
    1096     4852354 :   if (is_matvec_t(tx))
    1097             :   {
    1098          49 :     y = cgetg_copy(x, &lx);
    1099          49 :     for (i=1; i<lx; i++) gel(y,i) = gpow(gel(x,i),n,prec);
    1100          49 :     return y;
    1101             :   }
    1102     4852305 :   av = avma;
    1103     4852305 :   switch (tx)
    1104             :   {
    1105          28 :     case t_POL: case t_RFRAC: x = toser_i(x); /* fall through */
    1106             :     case t_SER:
    1107        7574 :       if (tn == t_FRAC) return gerepileupto(av, ser_powfrac(x, n, prec));
    1108         154 :       if (valp(x))
    1109          21 :         pari_err_DOMAIN("gpow [irrational exponent]",
    1110             :                         "valuation", "!=", gen_0, x);
    1111         133 :       if (lg(x) == 2) return gerepilecopy(av, x); /* O(1) */
    1112         126 :       return gerepileupto(av, ser_pow(x, n, prec));
    1113             :   }
    1114     4844731 :   if (gequal0(x)) return gpow0(x, n, prec);
    1115     4844654 :   if (tn == t_FRAC)
    1116             :   {
    1117     4425723 :     GEN z, a = gel(n,1), d = gel(n,2);
    1118             :     long D;
    1119     4425723 :     switch (tx)
    1120             :     {
    1121             :     case t_INT:
    1122             :     case t_FRAC:
    1123     4427242 :       if (ispower(x, d, &z)) return powgi(z, a);
    1124       90484 :       break;
    1125             :     case t_INTMOD:
    1126             :       {
    1127          21 :         GEN p = gel(x,1);
    1128          21 :         if (!BPSW_psp(p)) pari_err_PRIME("gpow",p);
    1129          14 :         y = cgetg(3,t_INTMOD); gel(y,1) = icopy(p);
    1130          14 :         av = avma;
    1131          14 :         z = Fp_sqrtn(gel(x,2), d, p, NULL);
    1132          14 :         if (!z) pari_err_SQRTN("gpow",x);
    1133           7 :         gel(y,2) = gerepileuptoint(av, Fp_pow(z, a, p));
    1134           7 :         return y;
    1135             :       }
    1136             : 
    1137             :     case t_PADIC:
    1138          14 :       z = Qp_sqrtn(x, d, NULL); if (!z) pari_err_SQRTN("gpow",x);
    1139           7 :       return gerepileupto(av, powgi(z, a));
    1140             : 
    1141             :     case t_FFELT:
    1142          21 :       return gerepileupto(av,FF_pow(FF_sqrtn(x,d,NULL),a));
    1143             :     default:
    1144     4332548 :       D = itos_or_0(d);
    1145     4332548 :       if (!D) break;
    1146     4332548 :       if (D == 2)
    1147             :       {
    1148     3980515 :         GEN y = gsqrt(x,prec), t = shifti(subiu(a,1), -1);
    1149     3980515 :         if (signe(t)) y = gmul(y, powgi(x,t));
    1150     3980515 :         return gerepileupto(av, y);
    1151             :       }
    1152      352033 :       if (is_real_t(tx) && gsigne(x) > 0)
    1153             :       {
    1154      350938 :         prec += nbits2extraprec(expi(a));
    1155      350938 :         if (tx != t_REAL) x = gtofp(x, prec);
    1156      350938 :         z = sqrtnr(x, D);
    1157      350938 :         if (!equali1(a)) z = powgi(z, a);
    1158      350938 :         return gerepileuptoleaf(av, z);
    1159             :       }
    1160             :     }
    1161             :   }
    1162      510510 :   i = precision(n);
    1163      510510 :   if (i) prec = i;
    1164      510510 :   prec0 = prec;
    1165      510510 :   if (!gprecision(x))
    1166             :   {
    1167      132729 :     long e = gexpo_safe(n); /* avoided if n = 0 or gexpo not defined */
    1168      132729 :     if (e > 2) prec += nbits2extraprec(e);
    1169             :   }
    1170      510510 :   y = gmul(n, glog(x,prec));
    1171      510454 :   y = gexp(y,prec);
    1172      510454 :   if (prec0 == prec) return gerepileupto(av, y);
    1173       33320 :   return gerepilecopy(av, gprec_wtrunc(y,prec0));
    1174             : }
    1175             : 
    1176             : GEN
    1177      163612 : gpowers0(GEN x, long n, GEN x0)
    1178             : {
    1179             :   long i, l;
    1180             :   GEN V;
    1181      163612 :   if (!x0) return gpowers(x,n);
    1182      150914 :   if (n < 0) return cgetg(1,t_VEC);
    1183      150914 :   l = n+2; V = cgetg(l, t_VEC); gel(V,1) = gcopy(x0);
    1184      150914 :   for (i = 2; i < l; i++) gel(V,i) = gmul(gel(V,i-1),x);
    1185      150914 :   return V;
    1186             : }
    1187             : 
    1188             : GEN
    1189      110505 : gpowers(GEN x, long n)
    1190             : {
    1191      110505 :   if (n < 0) return cgetg(1,t_VEC);
    1192      110498 :   return gen_powers(x, n, 1, (void*)x, &_sqr, &_mul, &_one);
    1193             : }
    1194             : 
    1195             : /* return [q^1,q^4,...,q^{n^2}] */
    1196             : GEN
    1197       25348 : gsqrpowers(GEN q, long n)
    1198             : {
    1199       25348 :   pari_sp av = avma;
    1200       25348 :   GEN L = gpowers0(gsqr(q), n, q); /* L[i] = q^(2i - 1), i <= n+1 */
    1201       25348 :   GEN v = cgetg(n+1, t_VEC);
    1202             :   long i;
    1203       25348 :   gel(v, 1) = gcopy(q);
    1204       25348 :   for (i = 2; i <= n ; ++i) gel(v, i) = q = gmul(q, gel(L,i)); /* q^(i^2) */
    1205       25348 :   return gerepileupto(av, v);
    1206             : }
    1207             : 
    1208             : /* 4 | N. returns a vector RU which contains exp(2*i*k*Pi/N), k=0..N-1 */
    1209             : static GEN
    1210       93259 : grootsof1_4(long N, long prec)
    1211             : {
    1212       93259 :   GEN z, RU = cgetg(N+1,t_VEC), *v  = ((GEN*)RU) + 1;
    1213       93259 :   long i, N2 = (N>>1), N4 = (N>>2), N8 = (N>>3);
    1214             :   /* z^N2 = -1, z^N4 = I; if z^k = a+I*b, then z^(N4-k) = I*conj(z) = b+a*I */
    1215             : 
    1216       93259 :   v[0] = gen_1; v[1] = z = rootsof1u_cx(N, prec);
    1217       93259 :   if (odd(N4)) N8++;
    1218      124842 :   for (i=1; i<N8; i++)
    1219             :   {
    1220       31583 :     GEN t = v[i];
    1221       31583 :     v[i+1] = gmul(z, t);
    1222       31583 :     v[N4-i] = mkcomplex(gel(t,2), gel(t,1));
    1223             :   }
    1224       93259 :   for (i=0; i<N4; i++) v[i+N4] = mulcxI(v[i]);
    1225       93259 :   for (i=0; i<N2; i++) v[i+N2] = gneg(v[i]);
    1226       93259 :   return RU;
    1227             : }
    1228             : 
    1229             : /* as above, N arbitrary */
    1230             : GEN
    1231      128250 : grootsof1(long N, long prec)
    1232             : {
    1233             :   GEN z, RU, *v;
    1234             :   long i, k;
    1235             : 
    1236      128250 :   if ((N & 3) == 0) return grootsof1_4(N, prec);
    1237       34991 :   if (N <= 2) return N == 1? mkvec(gen_1): mkvec2(gen_1, gen_m1);
    1238        9192 :   k = (N+3)>>1;
    1239        9192 :   RU = cgetg(N+1,t_VEC);
    1240        9192 :   v  = ((GEN*)RU) + 1;
    1241        9192 :   v[0] = gen_1; v[1] = z = rootsof1u_cx(N, prec);
    1242        9192 :   if (odd(N))
    1243        7766 :     for (i=2; i<k; i++) v[i] = gmul(z, v[i-1]);
    1244             :   else
    1245             :   {
    1246        1426 :     for (i=2; i<k-1; i++) v[i] = gmul(z, v[i-1]);
    1247        1426 :     v[i++] = gen_m1; /*avoid loss of accuracy*/
    1248             :   }
    1249        9192 :   for (   ; i<N; i++) v[i] = gconj(v[N-i]);
    1250        9192 :   return RU;
    1251             : }
    1252             : 
    1253             : /********************************************************************/
    1254             : /**                                                                **/
    1255             : /**                        RACINE CARREE                           **/
    1256             : /**                                                                **/
    1257             : /********************************************************************/
    1258             : /* assume x unit, e = precp(x) */
    1259             : GEN
    1260        1477 : Z2_sqrt(GEN x, long e)
    1261             : {
    1262        1477 :   ulong r = signe(x)>=0?mod16(x):16-mod16(x);
    1263             :   GEN z;
    1264             :   long ez;
    1265             :   pari_sp av;
    1266             : 
    1267        1477 :   switch(e)
    1268             :   {
    1269           7 :     case 1: return gen_1;
    1270         119 :     case 2: return (r & 3UL) == 1? gen_1: NULL;
    1271          14 :     case 3: return (r & 7UL) == 1? gen_1: NULL;
    1272          70 :     case 4: if (r == 1) return gen_1;
    1273          35 :             else return (r == 9)? utoipos(3): NULL;
    1274        1267 :     default: if ((r&7UL) != 1) return NULL;
    1275             :   }
    1276        1267 :   av = avma; z = (r==1)? gen_1: utoipos(3);
    1277        1267 :   ez = 3; /* number of correct bits in z (compared to sqrt(x)) */
    1278             :   for(;;)
    1279        3367 :   {
    1280             :     GEN mod;
    1281        4634 :     ez = (ez<<1) - 1;
    1282        4634 :     if (ez > e) ez = e;
    1283        4634 :     mod = int2n(ez);
    1284        4634 :     z = addii(z, remi2n(mulii(x, Fp_inv(z,mod)), ez));
    1285        4634 :     z = shifti(z, -1); /* (z + x/z) / 2 */
    1286        4634 :     if (e == ez) return gerepileuptoint(av, z);
    1287        3367 :     if (ez < e) ez--;
    1288        3367 :     if (gc_needed(av,2))
    1289             :     {
    1290           0 :       if (DEBUGMEM > 1) pari_warn(warnmem,"Qp_sqrt");
    1291           0 :       z = gerepileuptoint(av,z);
    1292             :     }
    1293             :   }
    1294             : }
    1295             : 
    1296             : /* x unit defined modulo p^e, e > 0 */
    1297             : GEN
    1298        1764 : Qp_sqrt(GEN x)
    1299             : {
    1300        1764 :   long pp, e = valp(x);
    1301        1764 :   GEN z,y,mod, p = gel(x,2);
    1302             : 
    1303        1764 :   if (gequal0(x)) return zeropadic(p, (e+1) >> 1);
    1304        1764 :   if (e & 1) return NULL;
    1305             : 
    1306        1750 :   y = cgetg(5,t_PADIC);
    1307        1750 :   pp = precp(x);
    1308        1750 :   mod = gel(x,3);
    1309        1750 :   z   = gel(x,4); /* lift to t_INT */
    1310        1750 :   e >>= 1;
    1311        1750 :   z = Zp_sqrt(z, p, pp);
    1312        1750 :   if (!z) return NULL;
    1313        1701 :   if (absequaliu(p,2))
    1314             :   {
    1315         805 :     pp  = (pp <= 3) ? 1 : pp-1;
    1316         805 :     mod = int2n(pp);
    1317             :   }
    1318         896 :   else mod = icopy(mod);
    1319        1701 :   y[1] = evalprecp(pp) | evalvalp(e);
    1320        1701 :   gel(y,2) = icopy(p);
    1321        1701 :   gel(y,3) = mod;
    1322        1701 :   gel(y,4) = z; return y;
    1323             : }
    1324             : 
    1325             : GEN
    1326         350 : Zn_sqrt(GEN d, GEN fn)
    1327             : {
    1328         350 :   pari_sp ltop = avma, btop;
    1329         350 :   GEN b = gen_0, m = gen_1;
    1330             :   long j, np;
    1331         350 :   if (typ(d) != t_INT) pari_err_TYPE("Zn_sqrt",d);
    1332         350 :   if (typ(fn) == t_INT)
    1333           0 :     fn = absZ_factor(fn);
    1334         350 :   else if (!is_Z_factorpos(fn))
    1335           0 :     pari_err_TYPE("Zn_sqrt",fn);
    1336         350 :   np = nbrows(fn);
    1337         350 :   btop = avma;
    1338        2800 :   for (j = 1; j <= np; ++j)
    1339             :   {
    1340             :     GEN  bp, mp, pr, r;
    1341        1050 :     GEN  p = gcoeff(fn, j, 1);
    1342        1050 :     long e = itos(gcoeff(fn, j, 2));
    1343        1050 :     long v = Z_pvalrem(d,p,&r);
    1344        1050 :     if (v >= e) bp =gen_0;
    1345             :     else
    1346             :     {
    1347         952 :       if (odd(v)) return NULL;
    1348         952 :       bp = Zp_sqrt(r, p, e-v);
    1349         952 :       if (!bp)    return NULL;
    1350         952 :       if (v) bp = mulii(bp, powiu(p, v>>1L));
    1351             :     }
    1352        1050 :     mp = powiu(p, e);
    1353        1050 :     pr = mulii(m, mp);
    1354        1050 :     b = Z_chinese_coprime(b, bp, m, mp, pr);
    1355        1050 :     m = pr;
    1356        1050 :     if (gc_needed(btop, 1))
    1357           0 :       gerepileall(btop, 2, &b, &m);
    1358             :   }
    1359         350 :   return gerepileupto(ltop, b);
    1360             : }
    1361             : 
    1362             : static GEN
    1363       17045 : sqrt_ser(GEN b, long prec)
    1364             : {
    1365       17045 :   long e = valp(b), vx = varn(b), lx, lold, j;
    1366             :   ulong mask;
    1367             :   GEN a, x, lta, ltx;
    1368             : 
    1369       17045 :   if (!signe(b)) return zeroser(vx, e>>1);
    1370       17045 :   a = leafcopy(b);
    1371       17045 :   x = cgetg_copy(b, &lx);
    1372       17045 :   if (e & 1)
    1373          14 :     pari_err_DOMAIN("sqrtn", "valuation", "!=", mkintmod(gen_0, gen_2), b);
    1374       17031 :   a[1] = x[1] = evalsigne(1) | evalvarn(0) | _evalvalp(0);
    1375       17031 :   lta = gel(a,2);
    1376       17031 :   if (gequal1(lta)) ltx = lta;
    1377       14791 :   else if (!issquareall(lta,&ltx)) ltx = gsqrt(lta,prec);
    1378       17024 :   gel(x,2) = ltx;
    1379       17024 :   for (j = 3; j < lx; j++) gel(x,j) = gen_0;
    1380       17024 :   setlg(x,3);
    1381       17024 :   mask = quadratic_prec_mask(lx - 2);
    1382       17024 :   lold = 1;
    1383       98801 :   while (mask > 1)
    1384             :   {
    1385       64753 :     GEN y, x2 = gmul2n(x,1);
    1386       64753 :     long l = lold << 1, lx;
    1387             : 
    1388       64753 :     if (mask & 1) l--;
    1389       64753 :     mask >>= 1;
    1390       64753 :     setlg(a, l + 2);
    1391       64753 :     setlg(x, l + 2);
    1392       64753 :     y = sqr_ser_part(x, lold, l-1) - lold;
    1393       64753 :     for (j = lold+2; j < l+2; j++) gel(y,j) = gsub(gel(y,j), gel(a,j));
    1394       64753 :     y += lold; setvalp(y, lold);
    1395       64753 :     y = normalize(y);
    1396       64753 :     y = gsub(x, gdiv(y, x2)); /* = gmul2n(gsub(x, gdiv(a,x)), -1); */
    1397       64753 :     lx = minss(l+2, lg(y));
    1398       64753 :     for (j = lold+2; j < lx; j++) gel(x,j) = gel(y,j);
    1399       64753 :     lold = l;
    1400             :   }
    1401       17024 :   x[1] = evalsigne(1) | evalvarn(vx) | _evalvalp(e >> 1);
    1402       17024 :   return x;
    1403             : }
    1404             : 
    1405             : GEN
    1406    11930550 : gsqrt(GEN x, long prec)
    1407             : {
    1408             :   pari_sp av;
    1409             :   GEN y;
    1410             : 
    1411    11930550 :   switch(typ(x))
    1412             :   {
    1413             :     case t_INT:
    1414     3605555 :       if (!signe(x)) return real_0(prec); /* no loss of accuracy */
    1415     3605548 :       x = itor(x,prec); /* fall through */
    1416    11304566 :     case t_REAL: return sqrtr(x);
    1417             : 
    1418             :     case t_INTMOD:
    1419             :     {
    1420          35 :       GEN p = gel(x,1), a;
    1421          35 :       y = cgetg(3,t_INTMOD); gel(y,1) = icopy(p);
    1422          35 :       a = Fp_sqrt(gel(x,2),p);
    1423          21 :       if (!a)
    1424             :       {
    1425           7 :         if (!BPSW_psp(p)) pari_err_PRIME("sqrt [modulus]",p);
    1426           7 :         pari_err_SQRTN("gsqrt",x);
    1427             :       }
    1428          14 :       gel(y,2) = a; return y;
    1429             :     }
    1430             : 
    1431             :     case t_COMPLEX:
    1432             :     { /* (u+iv)^2 = a+ib <=> u^2+v^2 = sqrt(a^2+b^2), u^2-v^2=a, 2uv=b */
    1433      524202 :       GEN a = gel(x,1), b = gel(x,2), r, u, v;
    1434      524202 :       if (isrationalzero(b)) return gsqrt(a, prec);
    1435      524202 :       y = cgetg(3,t_COMPLEX); av = avma;
    1436             : 
    1437      524202 :       r = cxnorm(x);
    1438      524202 :       if (typ(r) == t_INTMOD) pari_err_IMPL("sqrt(complex of t_INTMODs)");
    1439      524202 :       r = gsqrt(r, prec); /* t_REAL, |a+Ib| */
    1440      524202 :       if (!signe(r))
    1441         117 :         u = v = gerepileuptoleaf(av, sqrtr(r));
    1442      524085 :       else if (gsigne(a) < 0)
    1443             :       {
    1444             :         /* v > 0 since r > 0, a < 0, rounding errors can't make the sum of two
    1445             :          * positive numbers = 0 */
    1446       50347 :         v = sqrtr( gmul2n(gsub(r,a), -1) );
    1447       50347 :         if (gsigne(b) < 0) togglesign(v);
    1448       50347 :         v = gerepileuptoleaf(av, v); av = avma;
    1449             :         /* v = 0 is impossible */
    1450       50347 :         u = gerepileuptoleaf(av, gdiv(b, shiftr(v,1)));
    1451             :       } else {
    1452      473738 :         u = sqrtr( gmul2n(gadd(r,a), -1) );
    1453      473738 :         u = gerepileuptoleaf(av, u); av = avma;
    1454      473738 :         if (!signe(u)) /* possible if a = 0.0, e.g. sqrt(0.e-10+1e-10*I) */
    1455           7 :           v = u;
    1456             :         else
    1457      473731 :           v = gerepileuptoleaf(av, gdiv(b, shiftr(u,1)));
    1458             :       }
    1459      524202 :       gel(y,1) = u;
    1460      524202 :       gel(y,2) = v; return y;
    1461             :     }
    1462             : 
    1463             :     case t_PADIC:
    1464          63 :       y = Qp_sqrt(x);
    1465          63 :       if (!y) pari_err_SQRTN("Qp_sqrt",x);
    1466          42 :       return y;
    1467             : 
    1468          70 :     case t_FFELT: return FF_sqrt(x);
    1469             : 
    1470             :     default:
    1471      101607 :       av = avma; if (!(y = toser_i(x))) break;
    1472       17045 :       return gerepilecopy(av, sqrt_ser(y, prec));
    1473             :   }
    1474       84562 :   return trans_eval("sqrt",gsqrt,x,prec);
    1475             : }
    1476             : /********************************************************************/
    1477             : /**                                                                **/
    1478             : /**                          N-th ROOT                             **/
    1479             : /**                                                                **/
    1480             : /********************************************************************/
    1481             : static void
    1482           7 : bug_logp(GEN p)
    1483             : {
    1484           7 :   if (!BPSW_psp(p)) pari_err_PRIME("p-adic log",p);
    1485           0 :   pari_err_BUG("log_p");
    1486           0 : }
    1487             : /* Let x = 1 mod p and y := (x-1)/(x+1) = 0 (p). Then
    1488             :  * log(x) = log(1+y) - log(1-y) = 2 \sum_{k odd} y^k / k.
    1489             :  * palogaux(x) returns the last sum (not multiplied by 2) */
    1490             : static GEN
    1491       24885 : palogaux(GEN x)
    1492             : {
    1493             :   long i, k, e, pp, t;
    1494       24885 :   GEN y,s,y2, p = gel(x,2);
    1495       24885 :   int is2 = absequaliu(p,2);
    1496             : 
    1497       24885 :   y = subiu(gel(x,4), 1);
    1498       24885 :   if (!signe(y))
    1499             :   {
    1500         686 :     long v = valp(x)+precp(x);
    1501         686 :     if (is2) v--;
    1502         686 :     return zeropadic(p, v);
    1503             :   }
    1504             :   /* optimize t: log(x) = log(x^(p^t)) / p^t */
    1505       24199 :   e = Z_pval(y, p); /* valp(y) = e >= 1; precp(y) = precp(x)-e */
    1506       24199 :   if (!e) bug_logp(p);
    1507       24192 :   if (is2)
    1508        2499 :     t = sqrt( (double)(precp(x)-e) / e ); /* instead of (2*e) */
    1509             :   else
    1510       21693 :     t = sqrt( (double)(precp(x)-e) / (e * (expi(p) + hammingweight(p))) );
    1511       24192 :   for (i = 0; i < t; i++) x = gpow(x, p, 0);
    1512             : 
    1513       24192 :   y = gdiv(gaddgs(x,-1), gaddgs(x,1));
    1514       24192 :   e = valp(y); /* > 0 */
    1515       24192 :   if (e <= 0) bug_logp(p);
    1516       24192 :   pp = precp(y) + e;
    1517       24192 :   if (is2) pp--;
    1518             :   else
    1519             :   {
    1520             :     GEN p1;
    1521       21693 :     for (p1=utoipos(e); abscmpui(pp,p1) > 0; pp++) p1 = mulii(p1, p);
    1522       21693 :     pp -= 2;
    1523             :   }
    1524       24192 :   k = pp/e; if (!odd(k)) k--;
    1525       24192 :   if (DEBUGLEVEL>5)
    1526           0 :     err_printf("logp: [pp,k,e,t] = [%ld,%ld,%ld,%ld]\n",pp,k,e,t);
    1527       24192 :   if (k > 1)
    1528             :   {
    1529       22799 :     y2 = gsqr(y); s = gdivgs(gen_1,k);
    1530      123207 :     while (k > 2)
    1531             :     {
    1532       77609 :       k -= 2;
    1533       77609 :       s = gadd(gmul(y2,s), gdivgs(gen_1,k));
    1534             :     }
    1535       22799 :     y = gmul(s,y);
    1536             :   }
    1537       24192 :   if (t) setvalp(y, valp(y) - t);
    1538       24192 :   return y;
    1539             : }
    1540             : 
    1541             : GEN
    1542       24892 : Qp_log(GEN x)
    1543             : {
    1544       24892 :   pari_sp av = avma;
    1545       24892 :   GEN y, p = gel(x,2), a = gel(x,4);
    1546             : 
    1547       24892 :   if (!signe(a)) pari_err_DOMAIN("Qp_log", "argument", "=", gen_0, x);
    1548       24885 :   y = leafcopy(x); setvalp(y,0);
    1549       24885 :   if (absequaliu(p,2))
    1550        2821 :     y = palogaux(gsqr(y));
    1551       22064 :   else if (gequal1(modii(a, p)))
    1552       12988 :     y = gmul2n(palogaux(y), 1);
    1553             :   else
    1554             :   { /* compute log(x^(p-1)) / (p-1) */
    1555        9076 :     GEN mod = gel(y,3), p1 = subiu(p,1);
    1556        9076 :     gel(y,4) = Fp_pow(a, p1, mod);
    1557        9076 :     p1 = diviiexact(subsi(1,mod), p1); /* 1/(p-1) */
    1558        9076 :     y = gmul(palogaux(y), shifti(p1,1));
    1559             :   }
    1560       24878 :   return gerepileupto(av,y);
    1561             : }
    1562             : 
    1563             : static GEN Qp_exp_safe(GEN x);
    1564             : 
    1565             : /*compute the p^e th root of x p-adic, assume x != 0 */
    1566             : static GEN
    1567         854 : Qp_sqrtn_ram(GEN x, long e)
    1568             : {
    1569         854 :   pari_sp ltop=avma;
    1570         854 :   GEN a, p = gel(x,2), n = powiu(p,e);
    1571         854 :   long v = valp(x), va;
    1572         854 :   if (v)
    1573             :   {
    1574             :     long z;
    1575         161 :     v = sdivsi_rem(v, n, &z);
    1576         161 :     if (z) return NULL;
    1577          91 :     x = leafcopy(x);
    1578          91 :     setvalp(x,0);
    1579             :   }
    1580             :   /*If p = 2, -1 is a root of 1 in U1: need extra check*/
    1581         784 :   if (absequaliu(p, 2) && mod8(gel(x,4)) != 1) return NULL;
    1582         749 :   a = Qp_log(x);
    1583         749 :   va = valp(a) - e;
    1584         749 :   if (va <= 0)
    1585             :   {
    1586         287 :     if (signe(gel(a,4))) return NULL;
    1587             :     /* all accuracy lost */
    1588         119 :     a = cvtop(remii(gel(x,4),p), p, 1);
    1589             :   }
    1590             :   else
    1591             :   {
    1592         462 :     setvalp(a, va); /* divide by p^e */
    1593         462 :     a = Qp_exp_safe(a);
    1594         462 :     if (!a) return NULL;
    1595             :     /* n=p^e and a^n=z*x where z is a (p-1)th-root of 1.
    1596             :      * Since z^n=z, we have (a/z)^n = x. */
    1597         462 :     a = gdiv(x, powgi(a,subiu(n,1))); /* = a/z = x/a^(n-1)*/
    1598         462 :     if (v) setvalp(a,v);
    1599             :   }
    1600         581 :   return gerepileupto(ltop,a);
    1601             : }
    1602             : 
    1603             : /*compute the nth root of x p-adic p prime with n*/
    1604             : static GEN
    1605         616 : Qp_sqrtn_unram(GEN x, GEN n, GEN *zetan)
    1606             : {
    1607             :   pari_sp av;
    1608         616 :   GEN Z, a, r, p = gel(x,2);
    1609         616 :   long v = valp(x);
    1610         616 :   if (v)
    1611             :   {
    1612             :     long z;
    1613          84 :     v = sdivsi_rem(v,n,&z);
    1614          84 :     if (z) return NULL;
    1615             :   }
    1616         609 :   r = cgetp(x); setvalp(r,v);
    1617         609 :   Z = NULL; /* -Wall */
    1618         609 :   if (zetan) Z = cgetp(x);
    1619         609 :   av = avma; a = Fp_sqrtn(gel(x,4), n, p, zetan);
    1620         609 :   if (!a) return NULL;
    1621         595 :   affii(Zp_sqrtnlift(gel(x,4), n, a, p, precp(x)), gel(r,4));
    1622         595 :   if (zetan)
    1623             :   {
    1624          14 :     affii(Zp_sqrtnlift(gen_1, n, *zetan, p, precp(x)), gel(Z,4));
    1625          14 :     *zetan = Z;
    1626             :   }
    1627         595 :   set_avma(av); return r;
    1628             : }
    1629             : 
    1630             : GEN
    1631        1183 : Qp_sqrtn(GEN x, GEN n, GEN *zetan)
    1632             : {
    1633             :   pari_sp av, tetpil;
    1634             :   GEN q, p;
    1635             :   long e;
    1636        1183 :   if (absequaliu(n, 2))
    1637             :   {
    1638          70 :     if (zetan) *zetan = gen_m1;
    1639          70 :     if (signe(n) < 0) x = ginv(x);
    1640          63 :     return Qp_sqrt(x);
    1641             :   }
    1642        1113 :   av = avma; p = gel(x,2);
    1643        1113 :   if (!signe(gel(x,4)))
    1644             :   {
    1645         203 :     if (signe(n) < 0) pari_err_INV("Qp_sqrtn", x);
    1646         203 :     q = divii(addis(n, valp(x)-1), n);
    1647         203 :     if (zetan) *zetan = gen_1;
    1648         203 :     set_avma(av); return zeropadic(p, itos(q));
    1649             :   }
    1650             :   /* treat the ramified part using logarithms */
    1651         910 :   e = Z_pvalrem(n, p, &q);
    1652         910 :   if (e) { x = Qp_sqrtn_ram(x,e); if (!x) return NULL; }
    1653         637 :   if (is_pm1(q))
    1654             :   { /* finished */
    1655          21 :     if (signe(q) < 0) x = ginv(x);
    1656          21 :     x = gerepileupto(av, x);
    1657          21 :     if (zetan)
    1658          28 :       *zetan = (e && absequaliu(p, 2))? gen_m1 /*-1 in Q_2*/
    1659          21 :                                    : gen_1;
    1660          21 :     return x;
    1661             :   }
    1662         616 :   tetpil = avma;
    1663             :   /* use hensel lift for unramified case */
    1664         616 :   x = Qp_sqrtn_unram(x, q, zetan);
    1665         616 :   if (!x) return NULL;
    1666         595 :   if (zetan)
    1667             :   {
    1668             :     GEN *gptr[2];
    1669          14 :     if (e && absequaliu(p, 2))/*-1 in Q_2*/
    1670             :     {
    1671           7 :       tetpil = avma; x = gcopy(x); *zetan = gneg(*zetan);
    1672             :     }
    1673          14 :     gptr[0] = &x; gptr[1] = zetan;
    1674          14 :     gerepilemanysp(av,tetpil,gptr,2);
    1675          14 :     return x;
    1676             :   }
    1677         581 :   return gerepile(av,tetpil,x);
    1678             : }
    1679             : 
    1680             : GEN
    1681       23927 : sqrtnint(GEN a, long n)
    1682             : {
    1683       23927 :   pari_sp ltop = avma;
    1684             :   GEN x, b, q;
    1685             :   long s, k, e;
    1686       23927 :   const ulong nm1 = n - 1;
    1687       23927 :   if (typ(a) != t_INT) pari_err_TYPE("sqrtnint",a);
    1688       23927 :   if (n <= 0) pari_err_DOMAIN("sqrtnint", "n", "<=", gen_0, stoi(n));
    1689       23920 :   if (n == 1) return icopy(a);
    1690       23920 :   s = signe(a);
    1691       23920 :   if (s < 0) pari_err_DOMAIN("sqrtnint", "x", "<", gen_0, a);
    1692       23913 :   if (!s) return gen_0;
    1693       23913 :   if (lgefint(a) == 3) return utoi(usqrtn(itou(a), n));
    1694       23504 :   e = expi(a); k = e/(2*n);
    1695       23504 :   if (k == 0)
    1696             :   {
    1697             :     long flag;
    1698         291 :     if (n > e) {set_avma(ltop); return gen_1;}
    1699         291 :     flag = cmpii(a, powuu(3, n)); set_avma(ltop);
    1700         291 :     return (flag < 0) ? gen_2: stoi(3);
    1701             :   }
    1702       23213 :   if (e < n*BITS_IN_LONG - 1)
    1703             :   {
    1704             :     ulong xs, qs;
    1705        7084 :     b = itor(a, (2*e < n*BITS_IN_LONG)? DEFAULTPREC: MEDDEFAULTPREC);
    1706        7084 :     x = mpexp(divru(logr_abs(b), n));
    1707        7084 :     xs = itou(floorr(x)) + 1; /* >= a^(1/n) */
    1708             :     for(;;) {
    1709       21216 :       q = divii(a, powuu(xs, nm1));
    1710       14150 :       if (lgefint(q) > 3) break;
    1711       14143 :       qs = itou(q); if (qs >= xs) break;
    1712        7066 :       xs -= (xs - qs + nm1)/n;
    1713             :     }
    1714        7084 :     return utoi(xs);
    1715             :   }
    1716       16129 :   b = addui(1, shifti(a, -n*k));
    1717       16129 :   x = shifti(addui(1, sqrtnint(b, n)), k);
    1718       16129 :   q = divii(a, powiu(x, nm1));
    1719       50230 :   while (cmpii(q, x) < 0) /* a priori one iteration, no GC necessary */
    1720             :   {
    1721       17972 :     x = subii(x, divis(addui(nm1, subii(x, q)), n));
    1722       17972 :     q = divii(a, powiu(x, nm1));
    1723             :   }
    1724       16129 :   return gerepileuptoleaf(ltop, x);
    1725             : }
    1726             : 
    1727             : ulong
    1728         409 : usqrtn(ulong a, ulong n)
    1729             : {
    1730             :   ulong x, s, q;
    1731         409 :   const ulong nm1 = n - 1;
    1732         409 :   if (!n) pari_err_DOMAIN("sqrtnint", "n", "=", gen_0, utoi(n));
    1733         409 :   if (n == 1 || a == 0) return a;
    1734         409 :   s = 1 + expu(a)/n; x = 1UL << s;
    1735         409 :   q = (nm1*s >= BITS_IN_LONG)? 0: a >> (nm1*s);
    1736        1622 :   while (q < x) {
    1737             :     ulong X;
    1738         804 :     x -= (x - q + nm1)/n;
    1739         804 :     X = upowuu(x, nm1);
    1740         804 :     q = X? a/X: 0;
    1741             :   }
    1742         409 :   return x;
    1743             : }
    1744             : 
    1745             : static ulong
    1746      292535 : cubic_prec_mask(long n)
    1747             : {
    1748      292535 :   long a = n, i;
    1749      292535 :   ulong mask = 0;
    1750     1648344 :   for(i = 1;; i++, mask *= 3)
    1751     1355809 :   {
    1752     1648344 :     long c = a%3;
    1753     1648344 :     if (c) mask += 3 - c;
    1754     1648344 :     a = (a+2)/3;
    1755     1940879 :     if (a==1) return mask + upowuu(3, i);
    1756             :   }
    1757             : }
    1758             : 
    1759             : /* cubic Newton iteration, |a|^(1/n), assuming a != 0 */
    1760             : GEN
    1761      412540 : sqrtnr_abs(GEN a, long n)
    1762             : {
    1763             :   pari_sp av;
    1764             :   GEN x, b;
    1765             :   long eextra, eold, n1, n2, prec, B, v;
    1766             :   ulong mask;
    1767             : 
    1768      412540 :   if (n == 1) return mpabs(a);
    1769      411981 :   if (n == 2) return sqrtr_abs(a);
    1770             : 
    1771      391725 :   prec = realprec(a);
    1772      391725 :   B = prec2nbits(prec);
    1773      391725 :   eextra = expu(n)-1;
    1774      391725 :   n1 = n+1;
    1775      391725 :   n2 = 2*n; av = avma;
    1776      391725 :   v = expo(a) / n;
    1777      391725 :   if (v) a = shiftr(a, -n*v);
    1778             : 
    1779      391725 :   b = rtor(a, DEFAULTPREC);
    1780      391725 :   x = mpexp(divru(logr_abs(b), n));
    1781      391725 :   if (prec == DEFAULTPREC)
    1782             :   {
    1783      116784 :     if (v) shiftr_inplace(x, v);
    1784      116784 :     return gerepileuptoleaf(av, x);
    1785             :   }
    1786      274941 :   mask = cubic_prec_mask(B + 63);
    1787      274941 :   eold = 1;
    1788             :   for(;;)
    1789     1091390 :   { /* reach 64 */
    1790     1366331 :     long enew = eold * 3;
    1791     1366331 :     enew -= mask % 3;
    1792     1366331 :     if (enew > 64) break; /* back up one step */
    1793     1091390 :     mask /= 3;
    1794     1091390 :     eold = enew;
    1795             :   }
    1796             :   for(;;)
    1797      184831 :   {
    1798      459772 :     long pr, enew = eold * 3;
    1799             :     GEN y, z;
    1800      459772 :     enew -= mask % 3;
    1801      459772 :     mask /= 3;
    1802      459772 :     pr = nbits2prec(enew + eextra);
    1803      459772 :     b = rtor(a, pr); setsigne(b,1);
    1804      459772 :     x = rtor(x, pr);
    1805      459772 :     y = subrr(powru(x, n), b);
    1806      459772 :     z = divrr(y, addrr(mulur(n1, y), mulur(n2, b)));
    1807      459772 :     shiftr_inplace(z,1);
    1808      459772 :     x = mulrr(x, subsr(1,z));
    1809      459772 :     if (mask == 1)
    1810             :     {
    1811      274941 :       if (v) shiftr_inplace(x, v);
    1812      274941 :       return gerepileuptoleaf(av, gprec_wtrunc(x,prec));
    1813             :     }
    1814      184831 :     eold = enew;
    1815             :   }
    1816             : }
    1817             : 
    1818             : static void
    1819       28829 : shiftc_inplace(GEN z, long d)
    1820             : {
    1821       28829 :   shiftr_inplace(gel(z,1), d);
    1822       28829 :   shiftr_inplace(gel(z,2), d);
    1823       28829 : }
    1824             : 
    1825             : /* exp(2*Pi*I/n), same iteration as sqrtnr_abs, different initial point */
    1826             : static GEN
    1827       93332 : sqrtnof1(ulong n, long prec)
    1828             : {
    1829             :   pari_sp av;
    1830             :   GEN x;
    1831             :   long eold, n1, n2, B;
    1832             :   ulong mask;
    1833             : 
    1834       93332 :   B = prec2nbits(prec);
    1835       93332 :   n1 = n+1;
    1836       93332 :   n2 = 2*n; av = avma;
    1837             : 
    1838       93332 :   x = expIr(divru(Pi2n(1, LOWDEFAULTPREC), n));
    1839       93332 :   if (prec == LOWDEFAULTPREC) return gerepileupto(av, x);
    1840       17594 :   mask = cubic_prec_mask(B + BITS_IN_LONG-1);
    1841       17594 :   eold = 1;
    1842             :   for(;;)
    1843       68353 :   { /* reach BITS_IN_LONG */
    1844       85947 :     long enew = eold * 3;
    1845       85947 :     enew -= mask % 3;
    1846       85947 :     if (enew > BITS_IN_LONG) break; /* back up one step */
    1847       68353 :     mask /= 3;
    1848       68353 :     eold = enew;
    1849             :   }
    1850             :   for(;;)
    1851       11235 :   {
    1852       28829 :     long pr, enew = eold * 3;
    1853             :     GEN y, z;
    1854       28829 :     enew -= mask % 3;
    1855       28829 :     mask /= 3;
    1856       28829 :     pr = nbits2prec(enew);
    1857       28829 :     x = cxtofp(x, pr);
    1858       28829 :     y = gsub(gpowgs(x, n), gen_1);
    1859       28829 :     z = gdiv(y, gaddgs(gmulsg(n1, y), n2));
    1860       28829 :     shiftc_inplace(z,1);
    1861       28829 :     x = gmul(x, gsubsg(1, z));
    1862       28829 :     if (mask == 1) return gerepileupto(av, gprec_w(x,prec));
    1863       11235 :     eold = enew;
    1864             :   }
    1865             : }
    1866             : 
    1867             : /* exp(2iPi/d) */
    1868             : GEN
    1869      227853 : rootsof1u_cx(ulong n, long prec)
    1870             : {
    1871      227853 :   switch(n)
    1872             :   {
    1873       20559 :     case 1: return gen_1;
    1874        2534 :     case 2: return gen_m1;
    1875       39886 :     case 4: return gen_I();
    1876             :     case 3: case 6: case 12:
    1877             :     {
    1878        6071 :       pari_sp av = avma;
    1879        6071 :       GEN a = (n == 3)? mkfrac(gen_m1,gen_2): ghalf;
    1880        6071 :       GEN sq3 = sqrtr_abs(utor(3, prec));
    1881        6071 :       shiftr_inplace(sq3, -1);
    1882        6071 :       a = (n == 12)? mkcomplex(sq3, a): mkcomplex(a, sq3);
    1883        6071 :       return gerepilecopy(av, a);
    1884             :     }
    1885             :     case 8:
    1886             :     {
    1887       65471 :       pari_sp av = avma;
    1888       65471 :       GEN sq2 = sqrtr_abs(utor(2, prec));
    1889       65471 :       shiftr_inplace(sq2,-1);
    1890       65471 :       return gerepilecopy(av, mkcomplex(sq2, sq2));
    1891             :     }
    1892             :   }
    1893       93332 :   return sqrtnof1(n, prec);
    1894             : }
    1895             : /* e(a/b) */
    1896             : GEN
    1897       16380 : rootsof1q_cx(long a, long b, long prec)
    1898             : {
    1899       16380 :   long g = cgcd(a,b);
    1900             :   GEN z;
    1901       16380 :   if (g != 1) { a /= g; b /= g; }
    1902       16380 :   if (b < 0) { b = -b; a = -a; }
    1903       16380 :   z = rootsof1u_cx(b, prec);
    1904       16380 :   if (a < 0) { z = conj_i(z); a = -a; }
    1905       16380 :   return gpowgs(z, a);
    1906             : }
    1907             : 
    1908             : /* initializes powers of e(a/b) */
    1909             : GEN
    1910       12593 : rootsof1powinit(long a, long b, long prec)
    1911             : {
    1912       12593 :   long g = cgcd(a,b);
    1913       12593 :   if (g != 1) { a /= g; b /= g; }
    1914       12593 :   if (b < 0) { b = -b; a = -a; }
    1915       12593 :   a %= b; if (a < 0) a += b;
    1916       12593 :   return mkvec2(grootsof1(b,prec), mkvecsmall2(a,b));
    1917             : }
    1918             : /* T = rootsof1powinit(a,b); return  e(a/b)^c */
    1919             : GEN
    1920    10186659 : rootsof1pow(GEN T, long c)
    1921             : {
    1922    10186659 :   GEN vz = gel(T,1), ab = gel(T,2);
    1923    10186659 :   long a = ab[1], b = ab[2]; /* a >= 0, b > 0 */
    1924    10186659 :   c %= b; if (c < 0) c += b;
    1925    10186659 :   a = Fl_mul(a, c, b);
    1926    10186659 :   return gel(vz, a + 1);
    1927             : }
    1928             : 
    1929             : /* exp(2iPi/d), assume d a t_INT */
    1930             : GEN
    1931        3780 : rootsof1_cx(GEN d, long prec)
    1932             : {
    1933        3780 :   if (lgefint(d) == 3) return rootsof1u_cx((ulong)d[2], prec);
    1934           0 :   return expIr(divri(Pi2n(1,prec), d));
    1935             : }
    1936             : 
    1937             : GEN
    1938        3693 : gsqrtn(GEN x, GEN n, GEN *zetan, long prec)
    1939             : {
    1940             :   long i, lx, tx;
    1941             :   pari_sp av;
    1942             :   GEN y, z;
    1943        3693 :   if (typ(n)!=t_INT) pari_err_TYPE("sqrtn",n);
    1944        3693 :   if (!signe(n)) pari_err_DOMAIN("sqrtn", "n", "=", gen_0, n);
    1945        3693 :   if (is_pm1(n))
    1946             :   {
    1947          70 :     if (zetan) *zetan = gen_1;
    1948          70 :     return (signe(n) > 0)? gcopy(x): ginv(x);
    1949             :   }
    1950        3623 :   if (zetan) *zetan = gen_0;
    1951        3623 :   tx = typ(x);
    1952        3623 :   if (is_matvec_t(tx))
    1953             :   {
    1954           7 :     y = cgetg_copy(x, &lx);
    1955           7 :     for (i=1; i<lx; i++) gel(y,i) = gsqrtn(gel(x,i),n,NULL,prec);
    1956           7 :     return y;
    1957             :   }
    1958        3616 :   av = avma;
    1959        3616 :   switch(tx)
    1960             :   {
    1961             :   case t_INTMOD:
    1962             :     {
    1963         168 :       GEN p = gel(x,1), s;
    1964         168 :       z = gen_0;
    1965         168 :       y = cgetg(3,t_INTMOD);  gel(y,1) = icopy(p);
    1966         168 :       if (zetan) { z = cgetg(3,t_INTMOD); gel(z,1) = gel(y,1); }
    1967         168 :       s = Fp_sqrtn(gel(x,2),n,p,zetan);
    1968         147 :       if (!s) {
    1969          14 :         if (zetan) {set_avma(av); return gen_0;}
    1970           7 :         if (!BPSW_psp(p)) pari_err_PRIME("sqrtn [modulus]",p);
    1971           7 :         pari_err_SQRTN("gsqrtn",x);
    1972             :       }
    1973         133 :       gel(y,2) = s;
    1974         133 :       if (zetan) { gel(z,2) = *zetan; *zetan = z; }
    1975         133 :       return y;
    1976             :     }
    1977             : 
    1978             :   case t_PADIC:
    1979          56 :     y = Qp_sqrtn(x,n,zetan);
    1980          49 :     if (!y) {
    1981           7 :       if (zetan) return gen_0;
    1982           7 :       pari_err_SQRTN("gsqrtn",x);
    1983             :     }
    1984          42 :     return y;
    1985             : 
    1986         196 :   case t_FFELT: return FF_sqrtn(x,n,zetan);
    1987             : 
    1988             :   case t_INT: case t_FRAC: case t_REAL: case t_COMPLEX:
    1989        2832 :     i = precision(x); if (i) prec = i;
    1990        2832 :     if (isint1(x))
    1991           7 :       y = real_1(prec);
    1992        2825 :     else if (gequal0(x))
    1993             :     {
    1994             :       long b;
    1995          21 :       if (signe(n) < 0) pari_err_INV("gsqrtn",x);
    1996          21 :       if (isinexactreal(x))
    1997          14 :         b = sdivsi(gexpo(x), n);
    1998             :       else
    1999           7 :         b = -prec2nbits(prec);
    2000          21 :       if (typ(x) == t_COMPLEX)
    2001             :       {
    2002           7 :         y = cgetg(3,t_COMPLEX);
    2003           7 :         gel(y,1) = gel(y,2) = real_0_bit(b);
    2004             :       }
    2005             :       else
    2006          14 :         y = real_0_bit(b);
    2007             :     }
    2008             :     else
    2009             :     {
    2010        2804 :       long nn = itos_or_0(n);
    2011        2804 :       if (tx == t_INT) { x = itor(x,prec); tx = t_REAL; }
    2012        2804 :       if (nn > 0 && tx == t_REAL && signe(x) > 0)
    2013        1090 :         y = sqrtnr(x, nn);
    2014             :       else
    2015        1714 :         y = gexp(gdiv(glog(x,prec), n), prec);
    2016        2804 :       y = gerepileupto(av, y);
    2017             :     }
    2018        2832 :     if (zetan) *zetan = rootsof1_cx(n, prec);
    2019        2832 :     return y;
    2020             : 
    2021             :   case t_QUAD:
    2022           7 :     return gsqrtn(quadtofp(x, prec), n, zetan, prec);
    2023             : 
    2024             :   default:
    2025         357 :     av = avma; if (!(y = toser_i(x))) break;
    2026         357 :     return gerepileupto(av, ser_powfrac(y, ginv(n), prec));
    2027             :   }
    2028           0 :   pari_err_TYPE("sqrtn",x);
    2029             :   return NULL;/* LCOV_EXCL_LINE */
    2030             : }
    2031             : 
    2032             : /********************************************************************/
    2033             : /**                                                                **/
    2034             : /**                             EXP(X) - 1                         **/
    2035             : /**                                                                **/
    2036             : /********************************************************************/
    2037             : /* exp(|x|) - 1, assume x != 0.
    2038             :  * For efficiency, x should be reduced mod log(2): if so, we have a < 0 */
    2039             : GEN
    2040    10338178 : exp1r_abs(GEN x)
    2041             : {
    2042    10338178 :   long l = realprec(x), a = expo(x), b = prec2nbits(l), L, i, n, m, B;
    2043             :   GEN y, p2, X;
    2044             :   pari_sp av;
    2045             :   double d;
    2046             : 
    2047    10338178 :   if (b + a <= 0) return mpabs(x);
    2048             : 
    2049    10324017 :   y = cgetr(l); av = avma;
    2050    10324017 :   B = b/3 + BITS_IN_LONG + (BITS_IN_LONG*BITS_IN_LONG)/ b;
    2051    10324017 :   d = a/2.; m = (long)(d + sqrt(d*d + B)); /* >= 0 */
    2052    10324017 :   if (m < (-a) * 0.1) m = 0; /* not worth it */
    2053    10324017 :   L = l + nbits2extraprec(m);
    2054             :  /* Multiplication is quadratic in this range (l is small, otherwise we
    2055             :   * use logAGM + Newton). Set Y = 2^(-e-a) x, compute truncated series
    2056             :   * sum x^k/k!: this costs roughly
    2057             :   *    m b^2 + sum_{k <= n} (k e + BITS_IN_LONG)^2
    2058             :   * bit operations with |x| <  2^(1+a), |Y| < 2^(1-e) , m = e+a and b bits of
    2059             :   * accuracy needed, so
    2060             :   *    B := (b / 3 + BITS_IN_LONG + BITS_IN_LONG^2 / b) ~ m(m-a)
    2061             :   * we want b ~ 3 m (m-a) or m~b+a hence
    2062             :   *     m = min( a/2 + sqrt(a^2/4 + B),  b + a )
    2063             :   * NB: e ~ (b/3)^(1/2) as b -> oo
    2064             :   *
    2065             :   * Truncate the sum at k = n (>= 1), the remainder is
    2066             :   *   sum_{k >= n+1} Y^k / k! < Y^(n+1) / (n+1)! (1-Y) < Y^(n+1) / n!
    2067             :   * We want Y^(n+1) / n! <= Y 2^-b, hence -n log_2 |Y| + log_2 n! >= b
    2068             :   *   log n! ~ (n + 1/2) log(n+1) - (n+1) + log(2Pi)/2,
    2069             :   * error bounded by 1/6(n+1) <= 1/12. Finally, we want
    2070             :   * n (-1/log(2) -log_2 |Y| + log_2(n+1)) >= b  */
    2071    10324017 :   b += m;
    2072    10324017 :   d = m-dbllog2(x)-1/M_LN2; /* ~ -log_2 Y - 1/log(2) */
    2073    10324017 :   n = (long)(b / d);
    2074    10324017 :   if (n > 1)
    2075    10308904 :     n = (long)(b / (d + log2((double)n+1))); /* log~constant in small ranges */
    2076    10324017 :   while (n*(d+log2((double)n+1)) < b) n++; /* expect few corrections */
    2077             : 
    2078    10324017 :   X = rtor(x,L); shiftr_inplace(X, -m); setsigne(X, 1);
    2079    10324017 :   if (n == 1) p2 = X;
    2080             :   else
    2081             :   {
    2082    10324017 :     long s = 0, l1 = nbits2prec((long)(d + n + 16));
    2083    10324017 :     GEN unr = real_1(L);
    2084             :     pari_sp av2;
    2085             : 
    2086    10324017 :     p2 = cgetr(L); av2 = avma;
    2087   145904068 :     for (i=n; i>=2; i--, set_avma(av2))
    2088             :     { /* compute X^(n-1)/n! + ... + X/2 + 1 */
    2089             :       GEN p1, p3;
    2090   135580051 :       setprec(X,l1); p3 = divru(X,i);
    2091   135580051 :       l1 += dvmdsBIL(s - expo(p3), &s); if (l1>L) l1=L;
    2092   135580051 :       setprec(unr,l1); p1 = addrr_sign(unr,1, i == n? p3: mulrr(p3,p2),1);
    2093   135580051 :       setprec(p2,l1); affrr(p1,p2); /* p2 <- 1 + (X/i)*p2 */
    2094             :     }
    2095    10324017 :     setprec(X,L); p2 = mulrr(X,p2);
    2096             :   }
    2097             : 
    2098   113393342 :   for (i=1; i<=m; i++)
    2099             :   {
    2100   103069325 :     if (realprec(p2) > L) setprec(p2,L);
    2101   103069325 :     p2 = mulrr(p2, addsr(2,p2));
    2102             :   }
    2103    10324017 :   affrr_fixlg(p2,y); set_avma(av); return y;
    2104             : }
    2105             : 
    2106             : GEN
    2107        8991 : mpexpm1(GEN x)
    2108             : {
    2109        8991 :   const long s = 6;
    2110        8991 :   long l, sx = signe(x);
    2111             :   GEN y, z;
    2112             :   pari_sp av;
    2113        8991 :   if (!sx) return real_0_bit(expo(x));
    2114        8984 :   l = realprec(x);
    2115        8984 :   if (l > maxss(EXPNEWTON_LIMIT, (1L<<s) + 2))
    2116             :   {
    2117           0 :     long e = expo(x);
    2118           0 :     if (e < 0) x = rtor(x, l + nbits2extraprec(-e));
    2119           0 :     return subrs(mpexp(x), 1);
    2120             :   }
    2121        8984 :   if (sx > 0) return exp1r_abs(x);
    2122             :   /* compute exp(x) * (1 - exp(-x)) */
    2123        4149 :   av = avma; y = exp1r_abs(x);
    2124        4149 :   z = addsr(1, y); setsigne(z, -1);
    2125        4149 :   return gerepileupto(av, divrr(y, z));
    2126             : }
    2127             : 
    2128             : static GEN serexp(GEN x, long prec);
    2129             : GEN
    2130       10573 : gexpm1(GEN x, long prec)
    2131             : {
    2132       10573 :   switch(typ(x))
    2133             :   {
    2134        4060 :     case t_REAL: return mpexpm1(x);
    2135        4994 :     case t_COMPLEX: return cxexpm1(x,prec);
    2136          14 :     case t_PADIC: return gsubgs(Qp_exp(x), 1);
    2137             :     default:
    2138             :     {
    2139        1505 :       pari_sp av = avma;
    2140             :       long ey;
    2141             :       GEN y;
    2142        1505 :       if (!(y = toser_i(x))) break;
    2143        1484 :       ey = valp(y);
    2144        1484 :       if (ey < 0) pari_err_DOMAIN("expm1","valuation", "<", gen_0, x);
    2145        1484 :       if (gequal0(y)) return gcopy(y);
    2146        1477 :       if (ey)
    2147         343 :         return gerepileupto(av, gsubgs(serexp(y,prec), 1));
    2148             :       else
    2149             :       {
    2150        1134 :         GEN e1 = gexpm1(gel(y,2), prec), e = gaddgs(e1,1);
    2151        1134 :         y = gmul(e, serexp(serchop0(y),prec));
    2152        1134 :         gel(y,2) = e1;
    2153        1134 :         return gerepilecopy(av, y);
    2154             :       }
    2155             :     }
    2156             :   }
    2157          21 :   return trans_eval("expm1",gexpm1,x,prec);
    2158             : }
    2159             : /********************************************************************/
    2160             : /**                                                                **/
    2161             : /**                             EXP(X)                             **/
    2162             : /**                                                                **/
    2163             : /********************************************************************/
    2164             : 
    2165             : /* centermod(x, log(2)), set *sh to the quotient */
    2166             : static GEN
    2167    10277676 : modlog2(GEN x, long *sh)
    2168             : {
    2169    10277676 :   double d = rtodbl(x), da = fabs(d);
    2170    10277676 :   long q = (long) ((da + (M_LN2/2))/M_LN2);
    2171    10277676 :   if (da > M_LN2 * LONG_MAX)
    2172          14 :     pari_err_OVERFLOW("expo()"); /* avoid overflow in  q */
    2173    10277662 :   if (d < 0) q = -q;
    2174    10277662 :   *sh = q;
    2175    10277662 :   if (q) {
    2176     9818305 :     long l = realprec(x) + 1;
    2177     9818305 :     x = subrr(rtor(x,l), mulsr(q, mplog2(l)));
    2178     9818305 :     if (!signe(x)) return NULL;
    2179             :   }
    2180    10277662 :   return x;
    2181             : }
    2182             : 
    2183             : static GEN
    2184    10276985 : mpexp_basecase(GEN x)
    2185             : {
    2186    10276985 :   pari_sp av = avma;
    2187    10276985 :   long sh, l = realprec(x);
    2188             :   GEN y, z;
    2189             : 
    2190    10276985 :   y = modlog2(x, &sh);
    2191    10276971 :   if (!y) { set_avma(av); return real2n(sh, l); }
    2192    10276971 :   z = addsr(1, exp1r_abs(y));
    2193    10276971 :   if (signe(y) < 0) z = invr(z);
    2194    10276971 :   if (sh) {
    2195     9817614 :     shiftr_inplace(z, sh);
    2196     9817614 :     if (realprec(z) > l) z = rtor(z, l); /* spurious precision increase */
    2197             :   }
    2198             : #ifdef DEBUG
    2199             : {
    2200             :   GEN t = mplog(z), u = divrr(subrr(x, t),x);
    2201             :   if (signe(u) && expo(u) > 5-prec2nbits(minss(l,realprec(t))))
    2202             :     pari_err_BUG("exp");
    2203             : }
    2204             : #endif
    2205    10276971 :   return gerepileuptoleaf(av, z); /* NOT affrr, precision often increases */
    2206             : }
    2207             : 
    2208             : GEN
    2209    10484828 : mpexp(GEN x)
    2210             : {
    2211    10484828 :   const long s = 6; /*Initial steps using basecase*/
    2212    10484828 :   long i, p, l = realprec(x), sh;
    2213             :   GEN a, t, z;
    2214             :   ulong mask;
    2215             : 
    2216    10484828 :   if (l <= maxss(EXPNEWTON_LIMIT, (1L<<s) + 2))
    2217             :   {
    2218    10484137 :     if (!signe(x)) return mpexp0(x);
    2219    10276294 :     return mpexp_basecase(x);
    2220             :   }
    2221         691 :   z = cgetr(l); /* room for result */
    2222         691 :   x = modlog2(x, &sh);
    2223         691 :   if (!x) { avma = (pari_sp)(z+lg(z)); return real2n(sh, l); }
    2224         691 :   constpi(l); /* precompute for later logr_abs() */
    2225         691 :   mask = quadratic_prec_mask(prec2nbits(l)+BITS_IN_LONG);
    2226         691 :   for(i=0, p=1; i<s+TWOPOTBITS_IN_LONG; i++) { p <<= 1; if (mask & 1) p-=1; mask >>= 1; }
    2227         691 :   a = mpexp_basecase(rtor(x, nbits2prec(p)));
    2228         691 :   x = addrs(x,1);
    2229         691 :   if (realprec(x) < l+EXTRAPRECWORD) x = rtor(x, l+EXTRAPRECWORD);
    2230         691 :   a = rtor(a, l+EXTRAPRECWORD); /*append 0s */
    2231         691 :   t = NULL;
    2232             :   for(;;)
    2233             :   {
    2234         765 :     p <<= 1; if (mask & 1) p--;
    2235         728 :     mask >>= 1;
    2236         728 :     setprec(x, nbits2prec(p));
    2237         728 :     setprec(a, nbits2prec(p));
    2238         728 :     t = mulrr(a, subrr(x, logr_abs(a))); /* a (x - log(a)) */
    2239         728 :     if (mask == 1) break;
    2240          37 :     affrr(t, a); avma = (pari_sp)a;
    2241             :   }
    2242         691 :   affrr(t,z);
    2243         691 :   if (sh) shiftr_inplace(z, sh);
    2244         691 :   avma = (pari_sp)z; return z;
    2245             : }
    2246             : 
    2247             : static long
    2248       20148 : Qp_exp_prec(GEN x)
    2249             : {
    2250       20148 :   long k, e = valp(x), n = e + precp(x);
    2251       20148 :   GEN p = gel(x,2);
    2252       20148 :   int is2 = absequaliu(p,2);
    2253       20148 :   if (e < 1 || (e == 1 && is2)) return -1;
    2254       20120 :   if (is2)
    2255             :   {
    2256        6323 :     n--; e--; k = n/e;
    2257        6323 :     if (n%e == 0) k--;
    2258             :   }
    2259             :   else
    2260             :   { /* e > 0, n > 0 */
    2261       13797 :     GEN r, t = subiu(p, 1);
    2262       13797 :     k = itos(dvmdii(subiu(muliu(t,n), 1), subiu(muliu(t,e), 1), &r));
    2263       13797 :     if (!signe(r)) k--;
    2264             :   }
    2265       20120 :   return k;
    2266             : }
    2267             : 
    2268             : static GEN
    2269       21576 : Qp_exp_safe(GEN x)
    2270             : {
    2271             :   long k;
    2272             :   pari_sp av;
    2273             :   GEN y;
    2274             : 
    2275       21576 :   if (gequal0(x)) return gaddgs(x,1);
    2276       20078 :   k = Qp_exp_prec(x);
    2277       20077 :   if (k < 0) return NULL;
    2278       20070 :   av = avma;
    2279       20070 :   for (y=gen_1; k; k--) y = gaddsg(1, gdivgs(gmul(y,x), k));
    2280       20072 :   return gerepileupto(av, y);
    2281             : }
    2282             : 
    2283             : GEN
    2284       21114 : Qp_exp(GEN x)
    2285             : {
    2286       21114 :   GEN y = Qp_exp_safe(x);
    2287       21115 :   if (!y) pari_err_DOMAIN("gexp(t_PADIC)","argument","",gen_0,x);
    2288       21108 :   return y;
    2289             : }
    2290             : 
    2291             : static GEN
    2292          49 : cos_p(GEN x)
    2293             : {
    2294             :   long k;
    2295             :   pari_sp av;
    2296             :   GEN x2, y;
    2297             : 
    2298          49 :   if (gequal0(x)) return gaddgs(x,1);
    2299          28 :   k = Qp_exp_prec(x);
    2300          28 :   if (k < 0) return NULL;
    2301          21 :   av = avma; x2 = gsqr(x);
    2302          21 :   if (k & 1) k--;
    2303         105 :   for (y=gen_1; k; k-=2)
    2304             :   {
    2305          84 :     GEN t = gdiv(gmul(y,x2), muluu(k, k-1));
    2306          84 :     y = gsubsg(1, t);
    2307             :   }
    2308          21 :   return gerepileupto(av, y);
    2309             : }
    2310             : static GEN
    2311          63 : sin_p(GEN x)
    2312             : {
    2313             :   long k;
    2314             :   pari_sp av;
    2315             :   GEN x2, y;
    2316             : 
    2317          63 :   if (gequal0(x)) return gcopy(x);
    2318          42 :   k = Qp_exp_prec(x);
    2319          42 :   if (k < 0) return NULL;
    2320          28 :   av = avma; x2 = gsqr(x);
    2321          28 :   if (k & 1) k--;
    2322         133 :   for (y=gen_1; k; k-=2)
    2323             :   {
    2324         105 :     GEN t = gdiv(gmul(y,x2), muluu(k, k+1));
    2325         105 :     y = gsubsg(1, t);
    2326             :   }
    2327          28 :   return gerepileupto(av, gmul(y, x));
    2328             : }
    2329             : 
    2330             : static GEN
    2331     1351890 : cxexp(GEN x, long prec)
    2332             : {
    2333     1351890 :   GEN r, p1, p2, y = cgetg(3,t_COMPLEX);
    2334     1351890 :   pari_sp av = avma, tetpil;
    2335             :   long l;
    2336     1351890 :   l = precision(x); if (l > prec) prec = l;
    2337     1351890 :   r = gexp(gel(x,1),prec);
    2338     1351890 :   if (gequal0(r)) { gel(y,1) = r; gel(y,2) = r; return y; }
    2339     1351890 :   gsincos(gel(x,2),&p2,&p1,prec);
    2340     1351890 :   tetpil = avma;
    2341     1351890 :   gel(y,1) = gmul(r,p1);
    2342     1351890 :   gel(y,2) = gmul(r,p2);
    2343     1351890 :   gerepilecoeffssp(av,tetpil,y+1,2);
    2344     1351890 :   return y;
    2345             : }
    2346             : 
    2347             : /* given a t_SER x^v s(x), with s(0) != 0, return x^v(s - s(0)), shallow */
    2348             : GEN
    2349       28553 : serchop0(GEN s)
    2350             : {
    2351       28553 :   long i, l = lg(s);
    2352             :   GEN y;
    2353       28553 :   if (l == 2) return s;
    2354       28553 :   if (l == 3 && isexactzero(gel(s,2))) return s;
    2355       28553 :   y = cgetg(l, t_SER); y[1] = s[1];
    2356       28553 :   gel(y,2) = gen_0; for (i=3; i <l; i++) gel(y,i) = gel(s,i);
    2357       28553 :   return normalize(y);
    2358             : }
    2359             : 
    2360             : GEN
    2361          42 : serchop_i(GEN s, long n)
    2362             : {
    2363          42 :   long i, m, l = lg(s);
    2364             :   GEN y;
    2365          42 :   if (l == 2 || (l == 3 && isexactzero(gel(s,2))))
    2366             :   {
    2367          14 :     if (valp(s) < n) { s = shallowcopy(s); setvalp(s,n); }
    2368          14 :     return s;
    2369             :   }
    2370          28 :   m = n - valp(s); if (m < 0) return s;
    2371          21 :   if (l-m <= 2) return zeroser(varn(s), n);
    2372          14 :   y = cgetg(l-m, t_SER); y[1] = s[1]; setvalp(y, valp(y)+m);
    2373          14 :   for (i=m+2; i < l; i++) gel(y,i-m) = gel(s,i);
    2374          14 :   return normalize(y);
    2375             : }
    2376             : GEN
    2377          42 : serchop(GEN s, long n)
    2378             : {
    2379          42 :   pari_sp av = avma;
    2380          42 :   if (typ(s) != t_SER) pari_err_TYPE("serchop",s);
    2381          42 :   return gerepilecopy(av, serchop_i(s,n));
    2382             : }
    2383             : 
    2384             : static GEN
    2385       75376 : serexp(GEN x, long prec)
    2386             : {
    2387             :   pari_sp av;
    2388             :   long i,j,lx,ly,ex,mi;
    2389             :   GEN p1,y,xd,yd;
    2390             : 
    2391       75376 :   ex = valp(x);
    2392       75376 :   if (ex < 0) pari_err_DOMAIN("exp","valuation", "<", gen_0, x);
    2393       75369 :   if (gequal0(x)) return gaddsg(1,x);
    2394       55440 :   lx = lg(x);
    2395       55440 :   if (ex)
    2396             :   {
    2397       39732 :     ly = lx+ex; y = cgetg(ly,t_SER);
    2398       39732 :     mi = lx-1; while (mi>=3 && isrationalzero(gel(x,mi))) mi--;
    2399       39732 :     mi += ex-2;
    2400       39732 :     y[1] = evalsigne(1) | _evalvalp(0) | evalvarn(varn(x));
    2401             :     /* zd[i] = coefficient of X^i in z */
    2402       39732 :     xd = x+2-ex; yd = y+2; ly -= 2;
    2403       39732 :     gel(yd,0) = gen_1;
    2404       39732 :     for (i=1; i<ex; i++) gel(yd,i) = gen_0;
    2405      419160 :     for (   ; i<ly; i++)
    2406             :     {
    2407      379428 :       av = avma; p1 = gen_0;
    2408     2706200 :       for (j=ex; j<=minss(i,mi); j++)
    2409     2326772 :         p1 = gadd(p1, gmulgs(gmul(gel(xd,j),gel(yd,i-j)),j));
    2410      379428 :       gel(yd,i) = gerepileupto(av, gdivgs(p1,i));
    2411             :     }
    2412       39732 :     return y;
    2413             :   }
    2414       15708 :   av = avma;
    2415       15708 :   return gerepileupto(av, gmul(gexp(gel(x,2),prec), serexp(serchop0(x),prec)));
    2416             : }
    2417             : 
    2418             : GEN
    2419     8532676 : gexp(GEN x, long prec)
    2420             : {
    2421     8532676 :   switch(typ(x))
    2422             :   {
    2423     6939939 :     case t_REAL: return mpexp(x);
    2424     1351890 :     case t_COMPLEX: return cxexp(x,prec);
    2425          56 :     case t_PADIC: return Qp_exp(x);
    2426             :     default:
    2427             :     {
    2428      240791 :       pari_sp av = avma;
    2429             :       GEN y;
    2430      240791 :       if (!(y = toser_i(x))) break;
    2431       58191 :       return gerepileupto(av, serexp(y,prec));
    2432             :     }
    2433             :   }
    2434      182600 :   return trans_eval("exp",gexp,x,prec);
    2435             : }
    2436             : 
    2437             : /********************************************************************/
    2438             : /**                                                                **/
    2439             : /**                           AGM(X, Y)                            **/
    2440             : /**                                                                **/
    2441             : /********************************************************************/
    2442             : static int
    2443     4855897 : agmr_gap(GEN a, GEN b, long L)
    2444             : {
    2445     4855897 :   GEN d = subrr(b, a);
    2446     4855897 :   return (signe(d) && expo(d) - expo(b) >= L);
    2447             : }
    2448             : /* assume x > 0 */
    2449             : static GEN
    2450      342710 : agm1r_abs(GEN x)
    2451             : {
    2452      342710 :   long l = realprec(x), L = 5-prec2nbits(l);
    2453      342710 :   GEN a1, b1, y = cgetr(l);
    2454      342710 :   pari_sp av = avma;
    2455             : 
    2456      342710 :   a1 = addrr(real_1(l), x); shiftr_inplace(a1, -1);
    2457      342710 :   b1 = sqrtr_abs(x);
    2458     5198607 :   while (agmr_gap(a1,b1,L))
    2459             :   {
    2460     4513187 :     GEN a = a1;
    2461     4513187 :     a1 = addrr(a,b1); shiftr_inplace(a1, -1);
    2462     4513187 :     b1 = sqrtr_abs(mulrr(a,b1));
    2463             :   }
    2464      342710 :   affrr_fixlg(a1,y); set_avma(av); return y;
    2465             : }
    2466             : 
    2467             : struct agmcx_gap_t { long L, ex, cnt; };
    2468             : 
    2469             : static void
    2470       14014 : agmcx_init(GEN x, long *prec, struct agmcx_gap_t *S)
    2471             : {
    2472       14014 :   long l = precision(x);
    2473       14014 :   if (l) *prec = l;
    2474       14014 :   S->L = 1-prec2nbits(*prec);
    2475       14014 :   S->cnt = 0;
    2476       14014 :   S->ex = LONG_MAX;
    2477       14014 : }
    2478             : 
    2479             : static long
    2480       14014 : agmcx_a_b(GEN x, GEN *a1, GEN *b1, long prec)
    2481             : {
    2482       14014 :   long rotate = 0;
    2483       14014 :   if (gsigne(real_i(x))<0)
    2484             :   { /* Rotate by +/-Pi/2, so that the choice of the principal square
    2485             :      * root gives the optimal AGM. So a1 = +/-I*a1, b1=sqrt(-x). */
    2486        1477 :     if (gsigne(imag_i(x))<0) { *a1=mulcxI(*a1);  rotate=-1; }
    2487         987 :     else                     { *a1=mulcxmI(*a1); rotate=1; }
    2488        1477 :     x = gneg(x);
    2489             :   }
    2490       14014 :   *b1 = gsqrt(x, prec);
    2491       14014 :   return rotate;
    2492             : }
    2493             : /* return 0 if we must stop the AGM loop (a=b or a ~ b), 1 otherwise */
    2494             : static int
    2495      197913 : agmcx_gap(GEN a, GEN b, struct agmcx_gap_t *S)
    2496             : {
    2497      197913 :   GEN d = gsub(b, a);
    2498      197913 :   long ex = S->ex;
    2499      197913 :   S->ex = gexpo(d);
    2500      197913 :   if (gequal0(d) || S->ex - gexpo(b) < S->L) return 0;
    2501             :   /* if (S->ex >= ex) we're no longer making progress; twice in a row */
    2502      187210 :   if (S->ex < ex) S->cnt = 0;
    2503             :   else
    2504        6671 :     if (S->cnt++) return 0;
    2505      183899 :   return 1;
    2506             : }
    2507             : static GEN
    2508       13972 : agm1cx(GEN x, long prec)
    2509             : {
    2510             :   struct agmcx_gap_t S;
    2511             :   GEN a1, b1;
    2512       13972 :   pari_sp av = avma;
    2513             :   long rotate;
    2514       13972 :   agmcx_init(x, &prec, &S);
    2515       13972 :   a1 = gtofp(gmul2n(gadd(real_1(prec), x), -1), prec);
    2516       13972 :   rotate = agmcx_a_b(x, &a1, &b1, prec);
    2517      211585 :   while (agmcx_gap(a1,b1,&S))
    2518             :   {
    2519      183641 :     GEN a = a1;
    2520      183641 :     a1 = gmul2n(gadd(a,b1),-1);
    2521      183641 :     b1 = gsqrt(gmul(a,b1), prec);
    2522             :   }
    2523       13972 :   if (rotate) a1 = rotate>0 ? mulcxI(a1):mulcxmI(a1);
    2524       13972 :   return gerepilecopy(av,a1);
    2525             : }
    2526             : 
    2527             : GEN
    2528          42 : zellagmcx(GEN a0, GEN b0, GEN r, GEN t, long prec)
    2529             : {
    2530             :   struct agmcx_gap_t S;
    2531          42 :   pari_sp av = avma;
    2532          42 :   GEN x = gdiv(a0, b0), a1, b1;
    2533             :   long rotate;
    2534          42 :   agmcx_init(x, &prec, &S);
    2535          42 :   a1 = gtofp(gmul2n(gadd(real_1(prec), x), -1), prec);
    2536          42 :   r = gsqrt(gdiv(gmul(a1,gaddgs(r, 1)),gadd(r, x)), prec);
    2537          42 :   t = gmul(r, t);
    2538          42 :   rotate = agmcx_a_b(x, &a1, &b1, prec);
    2539         342 :   while (agmcx_gap(a1,b1,&S))
    2540             :   {
    2541         258 :     GEN a = a1, b = b1;
    2542         258 :     a1 = gmul2n(gadd(a,b),-1);
    2543         258 :     b1 = gsqrt(gmul(a,b), prec);
    2544         258 :     r = gsqrt(gdiv(gmul(a1,gaddgs(r, 1)),gadd(gmul(b, r), a )), prec);
    2545         258 :     t = gmul(r, t);
    2546             :   }
    2547          42 :   if (rotate) a1 = rotate>0 ? mulcxI(a1):mulcxmI(a1);
    2548          42 :   a1 = gmul(a1, b0);
    2549          42 :   t = gatan(gdiv(a1,t), prec);
    2550             :   /* send t to the fundamental domain if necessary */
    2551          42 :   if (gsigne(real_i(t))<0) t = gadd(t, mppi(prec));
    2552          42 :   return gerepileupto(av,gdiv(t,a1));
    2553             : }
    2554             : 
    2555             : static long
    2556          49 : ser_cmp_expo(GEN A, GEN B)
    2557             : {
    2558          49 :   long e = -(long)HIGHEXPOBIT, d = valp(B) - valp(A);
    2559          49 :   long i, la = lg(A), v = varn(B);
    2560        9849 :   for (i = 2; i < la; i++)
    2561             :   {
    2562        9800 :     GEN a = gel(A,i), b;
    2563             :     long ei;
    2564        9800 :     if (isexactzero(a)) continue;
    2565        9800 :     b = polcoef_i(B, i-2 + d, v);
    2566        9800 :     ei = gexpo(a);
    2567        9800 :     if (!isexactzero(b)) ei -= gexpo(b);
    2568        9800 :     e = maxss(e, ei);
    2569             :   }
    2570          49 :   return e;
    2571             : }
    2572             : 
    2573             : static GEN
    2574          14 : ser_agm1(GEN y, long prec)
    2575             : {
    2576          14 :   GEN a1 = y, b1 = gen_1;
    2577          14 :   long l = lg(y)-2, l2 = 6-prec2nbits(prec), eold = LONG_MAX;
    2578             :   for(;;)
    2579          63 :   {
    2580          77 :     GEN a = a1, p1;
    2581          77 :     a1 = gmul2n(gadd(a,b1),-1);
    2582          77 :     b1 = gsqrt(gmul(a,b1), prec);
    2583          77 :     p1 = gsub(b1,a1);
    2584          77 :     if (isinexactreal(p1))
    2585             :     {
    2586          49 :       long e = ser_cmp_expo(p1, b1);
    2587          49 :       if (e < l2 || e >= eold) break;
    2588          42 :       eold = e;
    2589             :     }
    2590             :     else
    2591             :     {
    2592          28 :       long ep = valp(p1)-valp(b1);
    2593          28 :       if (ep >= l && gequal0(p1)) break;
    2594             :     }
    2595             :   }
    2596          14 :   return a1;
    2597             : }
    2598             : 
    2599             : /* agm(1,x) */
    2600             : static GEN
    2601       10663 : agm1(GEN x, long prec)
    2602             : {
    2603             :   GEN y;
    2604             :   pari_sp av;
    2605             : 
    2606       10663 :   if (gequal0(x)) return gcopy(x);
    2607       10663 :   switch(typ(x))
    2608             :   {
    2609             :     case t_INT:
    2610          28 :       if (!is_pm1(x)) break;
    2611          21 :       return (signe(x) > 0)? real_1(prec): real_0(prec);
    2612             : 
    2613        5140 :     case t_REAL: return signe(x) > 0? agm1r_abs(x): agm1cx(x, prec);
    2614             : 
    2615             :     case t_COMPLEX:
    2616        5362 :       if (gequal0(gel(x,2))) return agm1(gel(x,1), prec);
    2617        5173 :       return agm1cx(x, prec);
    2618             : 
    2619             :     case t_PADIC:
    2620             :     {
    2621          14 :       GEN a1 = x, b1 = gen_1;
    2622          14 :       long l = precp(x);
    2623          14 :       av = avma;
    2624             :       for(;;)
    2625          14 :       {
    2626          28 :         GEN a = a1, p1;
    2627             :         long ep;
    2628          28 :         a1 = gmul2n(gadd(a,b1),-1);
    2629          28 :         a = gmul(a,b1);
    2630          28 :         b1 = Qp_sqrt(a); if (!b1) pari_err_SQRTN("Qp_sqrt",a);
    2631          21 :         p1 = gsub(b1,a1); ep = valp(p1)-valp(b1);
    2632          21 :         if (ep<=0) { b1 = gneg_i(b1); p1 = gsub(b1,a1); ep=valp(p1)-valp(b1); }
    2633          21 :         if (ep >= l || gequal0(p1)) return gerepilecopy(av,a1);
    2634             :       }
    2635             :     }
    2636             : 
    2637             :     default:
    2638         119 :       av = avma; if (!(y = toser_i(x))) break;
    2639          14 :       return gerepilecopy(av, ser_agm1(y, prec));
    2640             :   }
    2641         112 :   return trans_eval("agm",agm1,x,prec);
    2642             : }
    2643             : 
    2644             : GEN
    2645       10362 : agm(GEN x, GEN y, long prec)
    2646             : {
    2647             :   pari_sp av;
    2648       10362 :   if (is_matvec_t(typ(y)))
    2649             :   {
    2650          14 :     if (is_matvec_t(typ(x))) pari_err_TYPE2("agm",x,y);
    2651           7 :     swap(x, y);
    2652             :   }
    2653       10355 :   if (gequal0(y)) return gcopy(y);
    2654       10355 :   av = avma;
    2655       10355 :   return gerepileupto(av, gmul(y, agm1(gdiv(x,y), prec)));
    2656             : }
    2657             : 
    2658             : /********************************************************************/
    2659             : /**                                                                **/
    2660             : /**                             LOG(X)                             **/
    2661             : /**                                                                **/
    2662             : /********************************************************************/
    2663             : /* atanh(u/v) using binary splitting */
    2664             : static GEN
    2665        4992 : atanhQ_split(ulong u, ulong v, long prec)
    2666             : {
    2667             :   long i, nmax;
    2668        4992 :   GEN u2 = sqru(u), v2 = sqru(v);
    2669        4992 :   double d = ((double)v) / u;
    2670             :   struct abpq_res R;
    2671             :   struct abpq A;
    2672             :   /* satisfies (2n+1) (v/u)^2n > 2^bitprec */
    2673        4992 :   nmax = (long)(prec2nbits(prec) / (2*log2(d)));
    2674        4992 :   abpq_init(&A, nmax);
    2675        4992 :   A.a[0] = A.b[0] = gen_1;
    2676        4992 :   A.p[0] = utoipos(u);
    2677        4992 :   A.q[0] = utoipos(v);
    2678      432250 :   for (i = 1; i <= nmax; i++)
    2679             :   {
    2680      427258 :     A.a[i] = gen_1;
    2681      427258 :     A.b[i] = utoipos((i<<1)+1);
    2682      427258 :     A.p[i] = u2;
    2683      427258 :     A.q[i] = v2;
    2684             :   }
    2685        4992 :   abpq_sum(&R, 0, nmax, &A);
    2686        4992 :   return rdivii(R.T, mulii(R.B,R.Q),prec);
    2687             : }
    2688             : /* log(2) = 10*atanh(1/17)+4*atanh(13/499); faster than logagmr_abs()
    2689             :  * and Pi/2M(1,4/2^n) ~ n log(2) */
    2690             : static GEN
    2691        2496 : log2_split(long prec)
    2692             : {
    2693        2496 :   GEN u = atanhQ_split(1, 17, prec);
    2694        2496 :   GEN v = atanhQ_split(13, 499, prec);
    2695        2496 :   shiftr_inplace(v, 2);
    2696        2496 :   return addrr(mulur(10, u), v);
    2697             : }
    2698             : GEN
    2699    13310123 : constlog2(long prec)
    2700             : {
    2701             :   pari_sp av;
    2702             :   GEN tmp;
    2703    13310123 :   if (glog2 && realprec(glog2) >= prec) return glog2;
    2704             : 
    2705        2496 :   tmp = cgetr_block(prec);
    2706        2496 :   av = avma;
    2707        2496 :   affrr(log2_split(prec+EXTRAPRECWORD), tmp);
    2708        2496 :   swap_clone(&glog2,tmp);
    2709        2496 :   set_avma(av); return glog2;
    2710             : }
    2711             : 
    2712             : GEN
    2713    13310123 : mplog2(long prec) { return rtor(constlog2(prec), prec); }
    2714             : 
    2715             : /* dont check that q != 2^expo(q), done in logr_abs */
    2716             : static GEN
    2717      337612 : logagmr_abs(GEN q)
    2718             : {
    2719      337612 :   long prec = realprec(q), e = expo(q), lim;
    2720      337612 :   GEN z = cgetr(prec), y, Q, _4ovQ;
    2721      337612 :   pari_sp av = avma;
    2722             : 
    2723      337612 :   incrprec(prec);
    2724      337612 :   lim = prec2nbits(prec) >> 1;
    2725      337612 :   Q = rtor(q,prec);
    2726      337612 :   shiftr_inplace(Q,lim-e); setsigne(Q,1);
    2727             : 
    2728      337612 :   _4ovQ = invr(Q); shiftr_inplace(_4ovQ, 2); /* 4/Q */
    2729             :   /* Pi / 2agm(1, 4/Q) ~ log(Q), q = Q * 2^(e-lim) */
    2730      337612 :   y = divrr(Pi2n(-1, prec), agm1r_abs(_4ovQ));
    2731      337612 :   y = addrr(y, mulsr(e - lim, mplog2(prec)));
    2732      337612 :   affrr_fixlg(y, z); set_avma(av); return z;
    2733             : }
    2734             : 
    2735             : /* sum_{k >= 0} y^(2k+1) / (2k+1), y close to 0 */
    2736             : static GEN
    2737     3146110 : logr_aux(GEN y)
    2738             : {
    2739     3146110 :   long k, L = realprec(y); /* should be ~ l+1 - (k-2) */
    2740             :   /* log(x) = log(1+y) - log(1-y) = 2 sum_{k odd} y^k / k
    2741             :    * Truncate the sum at k = 2n+1, the remainder is
    2742             :    *   2 sum_{k >= 2n+3} y^k / k < 2y^(2n+3) / (2n+3)(1-y) < y^(2n+3)
    2743             :    * We want y^(2n+3) < y 2^(-prec2nbits(L)), hence
    2744             :    *   n+1 > -prec2nbits(L) /-log_2(y^2) */
    2745     3146110 :   double d = -2*dbllog2r(y); /* ~ -log_2(y^2) */
    2746     3146110 :   k = (long)(2*(prec2nbits(L) / d));
    2747     3146110 :   k |= 1;
    2748     3146110 :   if (k >= 3)
    2749             :   {
    2750     3126671 :     GEN T, S = cgetr(L), y2 = sqrr(y), unr = real_1(L);
    2751     3126671 :     pari_sp av = avma;
    2752     3126671 :     long s = 0, incs = (long)d, l1 = nbits2prec((long)d);
    2753     3126671 :     setprec(S,  l1);
    2754     3126671 :     setprec(unr,l1); affrr(divru(unr,k), S);
    2755    51094505 :     for (k -= 2;; k -= 2) /* k = 2n+1, ..., 1 */
    2756             :     { /* S = y^(2n+1-k)/(2n+1) + ... + 1 / k */
    2757    99062339 :       setprec(y2, l1); T = mulrr(S,y2);
    2758    51094505 :       if (k == 1) break;
    2759             : 
    2760    47967834 :       l1 += dvmdsBIL(s + incs, &s); if (l1>L) l1=L;
    2761    47967834 :       setprec(S, l1);
    2762    47967834 :       setprec(unr,l1);
    2763    47967834 :       affrr(addrr(divru(unr, k), T), S); set_avma(av);
    2764             :     }
    2765             :     /* k = 1 special-cased for eficiency */
    2766     3126671 :     y = mulrr(y, addsr(1,T)); /* = log(X)/2 */
    2767             :   }
    2768     3146110 :   return y;
    2769             : }
    2770             : /*return log(|x|), assuming x != 0 */
    2771             : GEN
    2772     3819586 : logr_abs(GEN X)
    2773             : {
    2774     3819586 :   long EX, L, m, k, a, b, l = realprec(X);
    2775             :   GEN z, x, y;
    2776             :   ulong u;
    2777             :   double d;
    2778             : 
    2779             :  /* Assuming 1 < x < 2, we want delta = x-1, 1-x/2, 1-1/x, or 2/x-1 small.
    2780             :   * We have 2/x-1 > 1-x/2, 1-1/x < x-1. So one should be choosing between
    2781             :   * 1-1/x and 1-x/2 ( crossover sqrt(2), worse ~ 0.29 ). To avoid an inverse,
    2782             :   * we choose between x-1 and 1-x/2 ( crossover 4/3, worse ~ 0.33 ) */
    2783     3819586 :   EX = expo(X);
    2784     3819586 :   u = uel(X,2);
    2785     3819586 :   k = 2;
    2786     3819586 :   if (u > (~0UL / 3) * 2) { /* choose 1-x/2 */
    2787     2031403 :     EX++; u = ~u;
    2788     2031403 :     while (!u && ++k < l) { u = uel(X,k); u = ~u; }
    2789             :   } else { /* choose x - 1 */
    2790     1788183 :     u &= ~HIGHBIT; /* u - HIGHBIT, assuming HIGHBIT set */
    2791     1788183 :     while (!u && ++k < l) u = uel(X,k);
    2792             :   }
    2793     3819586 :   if (k == l) return EX? mulsr(EX, mplog2(l)): real_0(l);
    2794     3483680 :   a = prec2nbits(k) + bfffo(u); /* ~ -log2 |1-x| */
    2795     3483680 :   L = l+1;
    2796     3483680 :   b = prec2nbits(L - (k-2)); /* take loss of accuracy into account */
    2797     5970011 :   if (b > 24*a*log2(L) &&
    2798     2823943 :       prec2nbits(l) > prec2nbits(LOGAGM_LIMIT)) return logagmr_abs(X);
    2799             : 
    2800     3146068 :   z = cgetr(EX? l: l - (k-2));
    2801             : 
    2802             :  /* Multiplication is quadratic in this range (l is small, otherwise we
    2803             :   * use AGM). Set Y = x^(1/2^m), y = (Y - 1) / (Y + 1) and compute truncated
    2804             :   * series sum y^(2k+1)/(2k+1): the costs is less than
    2805             :   *    m b^2 + sum_{k <= n} ((2k+1) e + BITS_IN_LONG)^2
    2806             :   * bit operations with |x-1| <  2^(1-a), |Y| < 2^(1-e) , m = e-a and b bits of
    2807             :   * accuracy needed (+ BITS_IN_LONG since bit accuracies increase by
    2808             :   * increments of BITS_IN_LONG), so
    2809             :   * 4n^3/3 e^2 + n^2 2e BITS_IN_LONG+ n BITS_IN_LONG ~ m b^2, with n ~ b/2e
    2810             :   * or b/6e + BITS_IN_LONG/2e + BITS_IN_LONG/2be ~ m
    2811             :   *    B := (b / 6 + BITS_IN_LONG/2 + BITS_IN_LONG^2 / 2b) ~ m(m+a)
    2812             :   *     m = min( -a/2 + sqrt(a^2/4 + B),  b - a )
    2813             :   * NB: e ~ (b/6)^(1/2) as b -> oo
    2814             :   * Instead of the above pessimistic estimate for the cost of the sum, use
    2815             :   * optimistic estimate (BITS_IN_LONG -> 0) */
    2816     3146068 :   d = -a/2.; m = (long)(d + sqrt(d*d + b/6)); /* >= 0 */
    2817             : 
    2818     3146068 :   if (m > b-a) m = b-a;
    2819     3146068 :   if (m < 0.2*a) m = 0; else L += nbits2extraprec(m);
    2820     3146068 :   x = rtor(X,L);
    2821     3146068 :   setsigne(x,1); shiftr_inplace(x,-EX);
    2822             :   /* 2/3 < x < 4/3 */
    2823     3146068 :   for (k=1; k<=m; k++) x = sqrtr_abs(x);
    2824             : 
    2825     3146068 :   y = divrr(subrs(x,1), addrs(x,1)); /* = (x-1) / (x+1), close to 0 */
    2826     3146068 :   y = logr_aux(y); /* log(1+y) - log(1-y) = log(x) */
    2827     3146068 :   shiftr_inplace(y, m + 1);
    2828     3146068 :   if (EX) y = addrr(y, mulsr(EX, mplog2(l+1)));
    2829     3146068 :   affrr_fixlg(y, z); avma = (pari_sp)z; return z;
    2830             : }
    2831             : 
    2832             : /* assume Im(q) != 0 and precision(q) >= prec. Compute log(q) with accuracy
    2833             :  * prec [disregard input accuracy] */
    2834             : GEN
    2835        8757 : logagmcx(GEN q, long prec)
    2836             : {
    2837        8757 :   GEN z = cgetc(prec), y, Q, a, b;
    2838             :   long lim, e, ea, eb;
    2839        8757 :   pari_sp av = avma;
    2840        8757 :   int neg = 0;
    2841             : 
    2842        8757 :   incrprec(prec);
    2843        8757 :   if (gsigne(gel(q,1)) < 0) { q = gneg(q); neg = 1; }
    2844        8757 :   lim = prec2nbits(prec) >> 1;
    2845        8757 :   Q = gtofp(q, prec);
    2846        8757 :   a = gel(Q,1);
    2847        8757 :   b = gel(Q,2);
    2848        8757 :   if (gequal0(a)) {
    2849           0 :     affrr_fixlg(logr_abs(b), gel(z,1));
    2850           0 :     y = Pi2n(-1, prec);
    2851           0 :     if (signe(b) < 0) setsigne(y, -1);
    2852           0 :     affrr_fixlg(y, gel(z,2)); set_avma(av); return z;
    2853             :   }
    2854        8757 :   ea = expo(a);
    2855        8757 :   eb = expo(b);
    2856        8757 :   e = ea <= eb ? lim - eb : lim - ea;
    2857        8757 :   shiftr_inplace(a, e);
    2858        8757 :   shiftr_inplace(b, e);
    2859             : 
    2860             :   /* Pi / 2agm(1, 4/Q) ~ log(Q), q = Q * 2^e */
    2861        8757 :   y = gdiv(Pi2n(-1, prec), agm1cx( gdivsg(4, Q), prec ));
    2862        8757 :   a = gel(y,1);
    2863        8757 :   b = gel(y,2);
    2864        8757 :   a = addrr(a, mulsr(-e, mplog2(prec)));
    2865        8757 :   if (realprec(a) <= LOWDEFAULTPREC) a = real_0_bit(expo(a));
    2866       15189 :   if (neg) b = gsigne(b) <= 0? gadd(b, mppi(prec))
    2867        6432 :                              : gsub(b, mppi(prec));
    2868        8757 :   affrr_fixlg(a, gel(z,1));
    2869        8757 :   affrr_fixlg(b, gel(z,2)); set_avma(av); return z;
    2870             : }
    2871             : 
    2872             : GEN
    2873      107159 : mplog(GEN x)
    2874             : {
    2875      107159 :   if (signe(x)<=0) pari_err_DOMAIN("mplog", "argument", "<=", gen_0, x);
    2876      107159 :   return logr_abs(x);
    2877             : }
    2878             : 
    2879             : /* pe = p^e, p prime, 0 < x < pe a t_INT coprime to p. Return the (p-1)-th
    2880             :  * root of 1 in (Z/pe)^* congruent to x mod p, resp x mod 4 if p = 2.
    2881             :  * Simplified form of Zp_sqrtnlift: 1/(p-1) is trivial to compute */
    2882             : GEN
    2883        1344 : Zp_teichmuller(GEN x, GEN p, long e, GEN pe)
    2884             : {
    2885             :   GEN q, z, p1;
    2886             :   pari_sp av;
    2887             :   ulong mask;
    2888        1344 :   if (absequaliu(p,2)) return (mod4(x) & 2)? subiu(pe,1): gen_1;
    2889        1344 :   if (e == 1) return icopy(x);
    2890        1344 :   av = avma;
    2891        1344 :   p1 = subiu(p, 1);
    2892        1344 :   mask = quadratic_prec_mask(e);
    2893        1344 :   q = p; z = remii(x, p);
    2894        7868 :   while (mask > 1)
    2895             :   { /* Newton iteration solving z^{1 - p} = 1, z = x (mod p) */
    2896        5180 :     GEN w, t, qold = q;
    2897        5180 :     if (mask <= 3) /* last iteration */
    2898        1344 :       q = pe;
    2899             :     else
    2900             :     {
    2901        3836 :       q = sqri(q);
    2902        3836 :       if (mask & 1) q = diviiexact(q, p);
    2903             :     }
    2904        5180 :     mask >>= 1;
    2905             :     /* q <= qold^2 */
    2906        5180 :     if (lgefint(q) == 3)
    2907             :     {
    2908        5032 :       ulong Z = uel(z,2), Q = uel(q,2), P1 = uel(p1,2);
    2909        5032 :       ulong W = (Q-1) / P1; /* -1/(p-1) + O(qold) */
    2910        5032 :       ulong T = Fl_mul(W, Fl_powu(Z,P1,Q) - 1, Q);
    2911        5032 :       Z = Fl_mul(Z, 1 + T, Q);
    2912        5032 :       z = utoi(Z);
    2913             :     }
    2914             :     else
    2915             :     {
    2916         148 :       w = diviiexact(subiu(qold,1),p1); /* -1/(p-1) + O(qold) */
    2917         148 :       t = Fp_mul(w, subiu(Fp_pow(z,p1,q), 1), q);
    2918         148 :       z = Fp_mul(z, addui(1,t), q);
    2919             :     }
    2920             :   }
    2921        1344 :   return gerepileuptoint(av, z);
    2922             : }
    2923             : 
    2924             : GEN
    2925        1225 : teichmullerinit(long p, long n)
    2926             : {
    2927             :   GEN t, pn, g, v;
    2928             :   ulong gp, tp;
    2929             :   long a, m;
    2930             : 
    2931        1225 :   if (p == 2) return mkvec(gen_1);
    2932        1225 :   if (!uisprime(p)) pari_err_PRIME("teichmullerinit",utoipos(p));
    2933             : 
    2934        1225 :   m = p >> 1; /* (p-1)/2 */
    2935        1225 :   tp= gp= pgener_Fl(p); /* order (p-1), gp^m = -1 */
    2936        1225 :   pn = powuu(p, n);
    2937        1225 :   v = cgetg(p, t_VEC);
    2938        1225 :   t = g = Zp_teichmuller(utoipos(gp), utoipos(p), n, pn);
    2939        1225 :   gel(v, 1) = gen_1;
    2940        1225 :   gel(v, p-1) = subiu(pn,1);
    2941        3031 :   for (a = 1; a < m; a++)
    2942             :   {
    2943        1806 :     gel(v, tp) = t;
    2944        1806 :     gel(v, p - tp) = Fp_neg(t, pn); /* g^(m+a) = -g^a */
    2945        1806 :     if (a < m-1)
    2946             :     {
    2947        1029 :       t = Fp_mul(t, g, pn); /* g^(a+1) */
    2948        1029 :       tp = Fl_mul(tp, gp, p); /* t mod p  */
    2949             :     }
    2950             :   }
    2951        1225 :   return v;
    2952             : }
    2953             : 
    2954             : /* tab from teichmullerinit or NULL */
    2955             : GEN
    2956         238 : teichmuller(GEN x, GEN tab)
    2957             : {
    2958             :   GEN p, q, y, z;
    2959         238 :   long n, tx = typ(x);
    2960             : 
    2961         238 :   if (!tab)
    2962             :   {
    2963         126 :     if (tx == t_VEC && lg(x) == 3)
    2964             :     {
    2965           7 :       p = gel(x,1);
    2966           7 :       q = gel(x,2);
    2967           7 :       if (typ(p) == t_INT && typ(q) == t_INT)
    2968           7 :         return teichmullerinit(itos(p), itos(q));
    2969             :     }
    2970             :   }
    2971         112 :   else if (typ(tab) != t_VEC) pari_err_TYPE("teichmuller",tab);
    2972         231 :   if (tx!=t_PADIC) pari_err_TYPE("teichmuller",x);
    2973         231 :   z = gel(x,4);
    2974         231 :   if (!signe(z)) return gcopy(x);
    2975         231 :   p = gel(x,2);
    2976         231 :   q = gel(x,3);
    2977         231 :   n = precp(x);
    2978         231 :   y = cgetg(5,t_PADIC);
    2979         231 :   y[1] = evalprecp(n) | _evalvalp(0);
    2980         231 :   gel(y,2) = icopy(p);
    2981         231 :   gel(y,3) = icopy(q);
    2982         231 :   if (tab)
    2983             :   {
    2984         112 :     ulong pp = itou_or_0(p);
    2985         112 :     if (lg(tab) != (long)pp) pari_err_TYPE("teichmuller",tab);
    2986         112 :     z = gel(tab, umodiu(z, pp));
    2987         112 :     if (typ(z) != t_INT) pari_err_TYPE("teichmuller",tab);
    2988         112 :     z = remii(z, q);
    2989             :   }
    2990             :   else
    2991         119 :     z = Zp_teichmuller(z, p, n, q);
    2992         231 :   gel(y,4) = z;
    2993         231 :   return y;
    2994             : }
    2995             : GEN
    2996           0 : teich(GEN x) { return teichmuller(x, NULL); }
    2997             : 
    2998             : GEN
    2999     3288406 : glog(GEN x, long prec)
    3000             : {
    3001             :   pari_sp av, tetpil;
    3002             :   GEN y, p1;
    3003             :   long l;
    3004             : 
    3005     3288406 :   switch(typ(x))
    3006             :   {
    3007             :     case t_REAL:
    3008     2211150 :       if (signe(x) >= 0)
    3009             :       {
    3010     1937283 :         if (!signe(x)) pari_err_DOMAIN("log", "argument", "=", gen_0, x);
    3011     1937269 :         return logr_abs(x);
    3012             :       }
    3013      273867 :       retmkcomplex(logr_abs(x), mppi(realprec(x)));
    3014             : 
    3015             :     case t_FRAC:
    3016             :     {
    3017             :       GEN a, b;
    3018             :       long e1, e2;
    3019      204068 :       av = avma;
    3020      204068 :       a = gel(x,1);
    3021      204068 :       b = gel(x,2);
    3022      204068 :       e1 = expi(subii(a,b)); e2 = expi(b);
    3023      204068 :       if (e2 > e1) prec += nbits2nlong(e2 - e1);
    3024      204068 :       x = fractor(x, prec);
    3025      204068 :       return gerepileupto(av, glog(x, prec));
    3026             :     }
    3027             :     case t_COMPLEX:
    3028      533257 :       if (ismpzero(gel(x,2))) return glog(gel(x,1), prec);
    3029      528752 :       l = precision(x); if (l > prec) prec = l;
    3030      528752 :       if (ismpzero(gel(x,1)))
    3031             :       {
    3032        5393 :         GEN a = gel(x,2), b;
    3033        5393 :         av = avma; b = Pi2n(-1,prec);
    3034        5393 :         if (gsigne(a) < 0) { setsigne(b, -1); a = gabs(a,prec); }
    3035        5393 :         a = isint1(a) ? gen_0: glog(a,prec);
    3036        5393 :         return gerepilecopy(av, mkcomplex(a, b));
    3037             :       }
    3038      523359 :       if (prec >= LOGAGMCX_LIMIT) return logagmcx(x, prec);
    3039      514798 :       y = cgetg(3,t_COMPLEX);
    3040      514798 :       gel(y,2) = garg(x,prec);
    3041      514798 :       av = avma; p1 = glog(cxnorm(x),prec); tetpil = avma;
    3042      514798 :       gel(y,1) = gerepile(av,tetpil,gmul2n(p1,-1)); return y;
    3043             : 
    3044         287 :     case t_PADIC: return Qp_log(x);
    3045             :     default:
    3046      339644 :       av = avma; if (!(y = toser_i(x))) break;
    3047         182 :       if (!signe(y)) pari_err_DOMAIN("log", "argument", "=", gen_0, x);
    3048         182 :       if (valp(y)) pari_err_DOMAIN("log", "series valuation", "!=", gen_0, x);
    3049         175 :       p1 = integser(gdiv(derivser(y), y)); /* log(y)' = y'/y */
    3050         175 :       if (!gequal1(gel(y,2))) p1 = gadd(p1, glog(gel(y,2),prec));
    3051         175 :       return gerepileupto(av, p1);
    3052             :   }
    3053      339462 :   return trans_eval("log",glog,x,prec);
    3054             : }
    3055             : 
    3056             : static GEN
    3057          42 : mplog1p(GEN x)
    3058             : {
    3059             :   long ex, a, b, l, L;
    3060          42 :   if (!signe(x)) return rcopy(x);
    3061          42 :   ex = expo(x); if (ex >= 0) return glog(addrs(x,1), 0);
    3062          42 :   a = -ex;
    3063          42 :   b = realprec(x); L = b+1;
    3064          42 :   if (b > a*log2(L) && prec2nbits(b) > prec2nbits(LOGAGM_LIMIT))
    3065             :   {
    3066           0 :     x = addrs(x,1); l = b + nbits2extraprec(a);
    3067           0 :     if (realprec(x) < l) x = rtor(x,l);
    3068           0 :     return logagmr_abs(x);
    3069             :   }
    3070          42 :   x = rtor(x, L);
    3071          42 :   x = logr_aux(divrr(x, addrs(x,2)));
    3072          42 :   if (realprec(x) > b) fixlg(x, b);
    3073          42 :   shiftr_inplace(x,1); return x;
    3074             : }
    3075             : 
    3076             : static GEN log1p_i(GEN x, long prec);
    3077             : static GEN
    3078          14 : cxlog1p(GEN x, long prec)
    3079             : {
    3080             :   pari_sp av;
    3081          14 :   GEN z, a, b = gel(x,2);
    3082             :   long l;
    3083          14 :   if (ismpzero(b)) return log1p_i(gel(x,1), prec);
    3084          14 :   l = precision(x); if (l > prec) prec = l;
    3085          14 :   if (prec >= LOGAGMCX_LIMIT) return logagmcx(gaddgs(x,1), prec);
    3086          14 :   a = gel(x,1);
    3087          14 :   z = cgetg(3,t_COMPLEX); av = avma;
    3088          14 :   a = gadd(gadd(gmul2n(a,1), gsqr(a)), gsqr(b));
    3089          14 :   a = log1p_i(a, prec); shiftr_inplace(a,-1);
    3090          14 :   gel(z,1) = gerepileupto(av, a);
    3091          14 :   gel(z,2) = garg(gaddgs(x,1),prec); return z;
    3092             : }
    3093             : static GEN
    3094          91 : log1p_i(GEN x, long prec)
    3095             : {
    3096          91 :   switch(typ(x))
    3097             :   {
    3098          42 :     case t_REAL: return mplog1p(x);
    3099          14 :     case t_COMPLEX: return cxlog1p(x, prec);
    3100           7 :     case t_PADIC: return Qp_log(gaddgs(x,1));
    3101             :     default:
    3102             :     {
    3103             :       long ey;
    3104             :       GEN y;
    3105          28 :       if (!(y = toser_i(x))) break;
    3106          21 :       ey = valp(y);
    3107          21 :       if (ey < 0) pari_err_DOMAIN("log1p","valuation", "<", gen_0, x);
    3108          21 :       if (gequal0(y)) return gcopy(y);
    3109          14 :       if (ey)
    3110           7 :         return glog(gaddgs(y,1),prec);
    3111             :       else
    3112             :       {
    3113           7 :         GEN a = gel(y,2), a1 = gaddgs(a,1);
    3114           7 :         y = gdiv(y, a1); gel(y,2) = gen_1;
    3115           7 :         return gadd(glog1p(a,prec), glog(y, prec));
    3116             :       }
    3117             :     }
    3118             :   }
    3119           7 :   return trans_eval("log1p",glog1p,x,prec);
    3120             : }
    3121             : GEN
    3122          77 : glog1p(GEN x, long prec)
    3123             : {
    3124          77 :   pari_sp av = avma;
    3125          77 :   return gerepileupto(av, log1p_i(x, prec));
    3126             : }
    3127             : /********************************************************************/
    3128             : /**                                                                **/
    3129             : /**                        SINE, COSINE                            **/
    3130             : /**                                                                **/
    3131             : /********************************************************************/
    3132             : 
    3133             : /* Reduce x0 mod Pi/2 to x in [-Pi/4, Pi/4]. Return cos(x)-1 */
    3134             : static GEN
    3135     4448570 : mpcosm1(GEN x, long *ptmod8)
    3136             : {
    3137     4448570 :   long a = expo(x), l = realprec(x), b, L, i, n, m, B;
    3138             :   GEN y, p2, x2;
    3139             :   double d;
    3140             : 
    3141     4448570 :   n = 0;
    3142     4448570 :   if (a >= 0)
    3143             :   {
    3144             :     long p;
    3145             :     GEN q;
    3146     3533442 :     if (a > 30)
    3147             :     {
    3148          92 :       GEN z, pitemp = Pi2n(-2, nbits2prec(a + 32));
    3149          92 :       z = addrr(x,pitemp); /* = x + Pi/4 */
    3150          92 :       if (expo(z) >= bit_prec(z) + 3) pari_err_PREC("mpcosm1");
    3151          92 :       shiftr_inplace(pitemp, 1);
    3152          92 :       q = floorr( divrr(z,pitemp) ); /* round ( x / (Pi/2) ) */
    3153          92 :       p = l+EXTRAPRECWORD; x = rtor(x,p);
    3154             :     } else {
    3155     3533350 :       q = stoi((long)floor(rtodbl(x) / (M_PI/2) + 0.5));
    3156     3533350 :       p = l;
    3157             :     }
    3158     3533442 :     if (signe(q))
    3159             :     {
    3160     3533442 :       x = subrr(x, mulir(q, Pi2n(-1,p))); /* x mod Pi/2  */
    3161     3533442 :       a = expo(x);
    3162     3533442 :       if (!signe(x) && a >= 0) pari_err_PREC("mpcosm1");
    3163     3533442 :       n = mod4(q); if (n && signe(q) < 0) n = 4 - n;
    3164             :     }
    3165             :   }
    3166             :   /* a < 0 */
    3167     4448570 :   b = signe(x); *ptmod8 = (b < 0)? 4 + n: n;
    3168     4448570 :   if (!b) return real_0_bit(expo(x)*2 - 1);
    3169             : 
    3170     4270830 :   b = prec2nbits(l);
    3171     4270830 :   if (b + 2*a <= 0) {
    3172      192644 :     y = sqrr(x); shiftr_inplace(y, -1); setsigne(y, -1);
    3173      192644 :     return y;
    3174             :   }
    3175             : 
    3176     4078186 :   y = cgetr(l);
    3177     4078186 :   B = b/6 + BITS_IN_LONG + (BITS_IN_LONG*BITS_IN_LONG/2)/ b;
    3178     4078186 :   d = a/2.; m = (long)(d + sqrt(d*d + B)); /* >= 0 ,*/
    3179     4078186 :   if (m < (-a) * 0.1) m = 0; /* not worth it */
    3180     4078186 :   L = l + nbits2extraprec(m);
    3181             : 
    3182     4078186 :   b += m;
    3183     4078186 :   d = 2.0 * (m-dbllog2r(x)-1/M_LN2); /* ~ 2( - log_2 Y - 1/log(2) ) */
    3184     4078186 :   n = (long)(b / d);
    3185     4078186 :   if (n > 1)
    3186     4047591 :     n = (long)(b / (d + log2((double)n+1))); /* log~constant in small ranges */
    3187     4078186 :   while (n*(d+log2((double)n+1)) < b) n++; /* expect few corrections */
    3188             : 
    3189             :  /* Multiplication is quadratic in this range (l is small, otherwise we
    3190             :   * use logAGM + Newton). Set Y = 2^(-e-a) x, compute truncated series
    3191             :   * sum Y^2k/(2k)!: this costs roughly
    3192             :   *   m b^2 + sum_{k <= n} (2k e + BITS_IN_LONG)^2
    3193             :   *   ~ floor(b/2e) b^2 / 3  + m b^2
    3194             :   * bit operations with |x| <  2^(1+a), |Y| < 2^(1-e) , m = e+a and b bits of
    3195             :   * accuracy needed, so
    3196             :   *    B := ( b / 6 + BITS_IN_LONG + BITS_IN_LONG^2 / 2b) ~ m(m-a)
    3197             :   * we want b ~ 6 m (m-a) or m~b+a hence
    3198             :   *     m = min( a/2 + sqrt(a^2/4 + b/6),  b/2 + a )
    3199             :   * NB1: e ~ (b/6)^(1/2) or b/2.
    3200             :   * NB2: We use b/4 instead of b/6 in the formula above: hand-optimized...
    3201             :   *
    3202             :   * Truncate the sum at k = n (>= 1), the remainder is
    3203             :   * < sum_{k >= n+1} Y^2k / 2k! < Y^(2n+2) / (2n+2)!(1-Y^2) < Y^(2n+2)/(2n+1)!
    3204             :   * We want ... <= Y^2 2^-b, hence -2n log_2 |Y| + log_2 (2n+1)! >= b
    3205             :   *   log n! ~ (n + 1/2) log(n+1) - (n+1) + log(2Pi)/2,
    3206             :   * error bounded by 1/6(n+1) <= 1/12. Finally, we want
    3207             :   * 2n (-1/log(2) - log_2 |Y| + log_2(2n+2)) >= b  */
    3208     4078186 :   x = rtor(x, L); shiftr_inplace(x, -m); setsigne(x, 1);
    3209     4078186 :   x2 = sqrr(x);
    3210     4078186 :   if (n == 1) { p2 = x2; shiftr_inplace(p2, -1); setsigne(p2, -1); } /*-Y^2/2*/
    3211             :   else
    3212             :   {
    3213     4078186 :     GEN unr = real_1(L);
    3214             :     pari_sp av;
    3215     4078186 :     long s = 0, l1 = nbits2prec((long)(d + n + 16));
    3216             : 
    3217     4078186 :     p2 = cgetr(L); av = avma;
    3218    26195505 :     for (i=n; i>=2; i--)
    3219             :     {
    3220             :       GEN p1;
    3221    22117319 :       setprec(x2,l1); p1 = divrunu(x2, 2*i-1);
    3222    22117319 :       l1 += dvmdsBIL(s - expo(p1), &s); if (l1>L) l1=L;
    3223    22117319 :       if (i != n) p1 = mulrr(p1,p2);
    3224    22117319 :       setprec(unr,l1); p1 = addrr_sign(unr,1, p1,-signe(p1));
    3225    22117319 :       setprec(p2,l1); affrr(p1,p2); set_avma(av);
    3226             :     }
    3227     4078186 :     shiftr_inplace(p2, -1); togglesign(p2); /* p2 := -p2/2 */
    3228     4078186 :     setprec(x2,L); p2 = mulrr(x2,p2);
    3229             :   }
    3230             :   /* Now p2 = sum {1<= i <=n} (-1)^i x^(2i) / (2i)! ~ cos(x) - 1 */
    3231    37428510 :   for (i=1; i<=m; i++)
    3232             :   { /* p2 = cos(x)-1 --> cos(2x)-1 */
    3233    33350324 :     p2 = mulrr(p2, addsr(2,p2));
    3234    33350324 :     shiftr_inplace(p2, 1);
    3235    33350324 :     if ((i & 31) == 0) p2 = gerepileuptoleaf((pari_sp)y, p2);
    3236             :   }
    3237     4078186 :   affrr_fixlg(p2,y); return y;
    3238             : }
    3239             : 
    3240             : /* sqrt (|1 - (1+x)^2|) = sqrt(|x*(x+2)|). Sends cos(x)-1 to |sin(x)| */
    3241             : static GEN
    3242     2988922 : mpaut(GEN x)
    3243             : {
    3244     2988922 :   pari_sp av = avma;
    3245     2988922 :   GEN t = mulrr(x, addsr(2,x)); /* != 0 */
    3246     2988922 :   if (!signe(t)) return real_0_bit(expo(t) >> 1);
    3247     2815844 :   return gerepileuptoleaf(av, sqrtr_abs(t));
    3248             : }
    3249             : 
    3250             : /********************************************************************/
    3251             : /**                            COSINE                              **/
    3252             : /********************************************************************/
    3253             : 
    3254             : GEN
    3255     2701948 : mpcos(GEN x)
    3256             : {
    3257             :   long mod8;
    3258             :   pari_sp av;
    3259             :   GEN y,p1;
    3260             : 
    3261     2701948 :   if (!signe(x)) {
    3262           8 :     long l = nbits2prec(-expo(x));
    3263           8 :     if (l < LOWDEFAULTPREC) l = LOWDEFAULTPREC;
    3264           8 :     return real_1(l);
    3265             :   }
    3266             : 
    3267     2701940 :   av = avma; p1 = mpcosm1(x,&mod8);
    3268     2701940 :   switch(mod8)
    3269             :   {
    3270      748775 :     case 0: case 4: y = addsr(1,p1); break;
    3271      677991 :     case 1: case 7: y = mpaut(p1); togglesign(y); break;
    3272      671600 :     case 2: case 6: y = subsr(-1,p1); break;
    3273      603574 :     default:        y = mpaut(p1); break; /* case 3: case 5: */
    3274             :   }
    3275     2701940 :   return gerepileuptoleaf(av, y);
    3276             : }
    3277             : 
    3278             : /* convert INT or FRAC to REAL, which is later reduced mod 2Pi : avoid
    3279             :  * cancellation */
    3280             : static GEN
    3281       12425 : tofp_safe(GEN x, long prec)
    3282             : {
    3283       23660 :   return (typ(x) == t_INT || gexpo(x) > 0)? gadd(x, real_0(prec))
    3284       14385 :                                           : fractor(x, prec);
    3285             : }
    3286             : 
    3287             : GEN
    3288      111742 : gcos(GEN x, long prec)
    3289             : {
    3290             :   pari_sp av;
    3291             :   GEN r, u, v, y, u1, v1;
    3292             :   long i;
    3293             : 
    3294      111742 :   switch(typ(x))
    3295             :   {
    3296      111070 :     case t_REAL: return mpcos(x);
    3297             :     case t_COMPLEX:
    3298          42 :       if (isintzero(gel(x,1))) return gcosh(gel(x,2), prec);
    3299          28 :       i = precision(x); if (i) prec = i;
    3300          28 :       y = cgetc(prec); av = avma;
    3301          28 :       r = gexp(gel(x,2),prec);
    3302          28 :       v1 = gmul2n(addrr(invr(r),r), -1); /* = cos(I*Im(x)) */
    3303          28 :       u1 = subrr(v1, r); /* = - I*sin(I*Im(x)) */
    3304          28 :       gsincos(gel(x,1),&u,&v,prec);
    3305          28 :       affrr_fixlg(gmul(v1,v), gel(y,1));
    3306          28 :       affrr_fixlg(gmul(u1,u), gel(y,2)); set_avma(av); return y;
    3307             : 
    3308             :     case t_INT: case t_FRAC:
    3309         546 :       y = cgetr(prec); av = avma;
    3310         546 :       affrr_fixlg(mpcos(tofp_safe(x,prec)), y); set_avma(av); return y;
    3311             : 
    3312          49 :     case t_PADIC: y = cos_p(x);
    3313          49 :       if (!y) pari_err_DOMAIN("gcos(t_PADIC)","argument","",gen_0,x);
    3314          42 :       return y;
    3315             : 
    3316             :     default:
    3317          35 :       av = avma; if (!(y = toser_i(x))) break;
    3318          28 :       if (gequal0(y)) return gerepileupto(av, gaddsg(1,y));
    3319          28 :       if (valp(y) < 0)
    3320           7 :         pari_err_DOMAIN("cos","valuation", "<", gen_0, x);
    3321          21 :       gsincos(y,&u,&v,prec);
    3322          21 :       return gerepilecopy(av,v);
    3323             :   }
    3324           7 :   return trans_eval("cos",gcos,x,prec);
    3325             : }
    3326             : /********************************************************************/
    3327             : /**                             SINE                               **/
    3328             : /********************************************************************/
    3329             : 
    3330             : GEN
    3331      139967 : mpsin(GEN x)
    3332             : {
    3333             :   long mod8;
    3334             :   pari_sp av;
    3335             :   GEN y,p1;
    3336             : 
    3337      139967 :   if (!signe(x)) return real_0_bit(expo(x));
    3338             : 
    3339      139819 :   av = avma; p1 = mpcosm1(x,&mod8);
    3340      139819 :   switch(mod8)
    3341             :   {
    3342       76540 :     case 0: case 6: y=mpaut(p1); break;
    3343       21869 :     case 1: case 5: y=addsr(1,p1); break;
    3344       24006 :     case 2: case 4: y=mpaut(p1); togglesign(y); break;
    3345       17404 :     default:        y=subsr(-1,p1); break; /* case 3: case 7: */
    3346             :   }
    3347      139819 :   return gerepileuptoleaf(av, y);
    3348             : }
    3349             : 
    3350             : GEN
    3351      171047 : gsin(GEN x, long prec)
    3352             : {
    3353             :   pari_sp av;
    3354             :   GEN r, u, v, y, v1, u1;
    3355             :   long i;
    3356             : 
    3357      171047 :   switch(typ(x))
    3358             :   {
    3359      128382 :     case t_REAL: return mpsin(x);
    3360             :     case t_COMPLEX:
    3361       31045 :       if (isintzero(gel(x,1))) retmkcomplex(gen_0,gsinh(gel(x,2),prec));
    3362       15722 :       i = precision(x); if (i) prec = i;
    3363       15722 :       y = cgetc(prec); av = avma;
    3364       15722 :       r = gexp(gel(x,2),prec);
    3365       15722 :       v1 = gmul2n(addrr(invr(r),r), -1); /* = cos(I*Im(x)) */
    3366       15722 :       u1 = subrr(r, v1); /* = I*sin(I*Im(x)) */
    3367       15722 :       gsincos(gel(x,1),&u,&v,prec);
    3368       15722 :       affrr_fixlg(gmul(v1,u), gel(y,1));
    3369       15722 :       affrr_fixlg(gmul(u1,v), gel(y,2)); set_avma(av); return y;
    3370             : 
    3371             :     case t_INT: case t_FRAC:
    3372       11354 :       y = cgetr(prec); av = avma;
    3373       11354 :       affrr_fixlg(mpsin(tofp_safe(x,prec)), y); set_avma(av); return y;
    3374             : 
    3375          49 :     case t_PADIC: y = sin_p(x);
    3376          49 :       if (!y) pari_err_DOMAIN("gsin(t_PADIC)","argument","",gen_0,x);
    3377          42 :       return y;
    3378             : 
    3379             :     default:
    3380         217 :       av = avma; if (!(y = toser_i(x))) break;
    3381         210 :       if (gequal0(y)) return gerepilecopy(av, y);
    3382         210 :       if (valp(y) < 0)
    3383           7 :         pari_err_DOMAIN("sin","valuation", "<", gen_0, x);
    3384         203 :       gsincos(y,&u,&v,prec);
    3385         203 :       return gerepilecopy(av,u);
    3386             :   }
    3387           7 :   return trans_eval("sin",gsin,x,prec);
    3388             : }
    3389             : /********************************************************************/
    3390             : /**                       SINE, COSINE together                    **/
    3391             : /********************************************************************/
    3392             : 
    3393             : void
    3394     1605126 : mpsincos(GEN x, GEN *s, GEN *c)
    3395             : {
    3396             :   long mod8;
    3397             :   pari_sp av, tetpil;
    3398             :   GEN p1, *gptr[2];
    3399             : 
    3400     1605126 :   if (!signe(x))
    3401             :   {
    3402        3225 :     long e = expo(x);
    3403        3225 :     *s = real_0_bit(e);
    3404        3225 :     *c = e >= 0? real_0_bit(e): real_1_bit(-e);
    3405        6450 :     return;
    3406             :   }
    3407             : 
    3408     1601901 :   av=avma; p1=mpcosm1(x,&mod8); tetpil=avma;
    3409     1601901 :   switch(mod8)
    3410             :   {
    3411      358277 :     case 0: *c=addsr( 1,p1); *s=mpaut(p1); break;
    3412      121690 :     case 1: *s=addsr( 1,p1); *c=mpaut(p1); togglesign(*c); break;
    3413      345496 :     case 2: *c=subsr(-1,p1); *s=mpaut(p1); togglesign(*s); break;
    3414      111715 :     case 3: *s=subsr(-1,p1); *c=mpaut(p1); break;
    3415      232753 :     case 4: *c=addsr( 1,p1); *s=mpaut(p1); togglesign(*s); break;
    3416      138864 :     case 5: *s=addsr( 1,p1); *c=mpaut(p1); break;
    3417      177335 :     case 6: *c=subsr(-1,p1); *s=mpaut(p1); break;
    3418      115771 :     case 7: *s=subsr(-1,p1); *c=mpaut(p1); togglesign(*c); break;
    3419             :   }
    3420     1601901 :   gptr[0]=s; gptr[1]=c;
    3421     1601901 :   gerepilemanysp(av,tetpil,gptr,2);
    3422             : }
    3423             : 
    3424             : /* SINE and COSINE - 1 */
    3425             : void
    3426        4910 : mpsincosm1(GEN x, GEN *s, GEN *c)
    3427             : {
    3428             :   long mod8;
    3429             :   pari_sp av, tetpil;
    3430             :   GEN p1, *gptr[2];
    3431             : 
    3432        4910 :   if (!signe(x))
    3433             :   {
    3434           0 :     long e = expo(x);
    3435           0 :     *s = real_0_bit(e);
    3436           0 :     *c = real_0_bit(2*e-1);
    3437           0 :     return;
    3438             :   }
    3439        4910 :   av=avma; p1=mpcosm1(x,&mod8); tetpil=avma;
    3440        4910 :   switch(mod8)
    3441             :   {
    3442        4490 :     case 0: *c=rcopy(p1); *s=mpaut(p1); break;
    3443          42 :     case 1: *s=addsr(1,p1); *c=addrs(mpaut(p1),1); togglesign(*c); break;
    3444           0 :     case 2: *c=subsr(-2,p1); *s=mpaut(p1); togglesign(*s); break;
    3445           0 :     case 3: *s=subsr(-1,p1); *c=subrs(mpaut(p1),1); break;
    3446         287 :     case 4: *c=rcopy(p1); *s=mpaut(p1); togglesign(*s); break;
    3447          77 :     case 5: *s=addsr( 1,p1); *c=subrs(mpaut(p1),1); break;
    3448           7 :     case 6: *c=subsr(-2,p1); *s=mpaut(p1); break;
    3449           7 :     case 7: *s=subsr(-1,p1); *c=subsr(-1,mpaut(p1)); break;
    3450             :   }
    3451        4910 :   gptr[0]=s; gptr[1]=c;
    3452        4910 :   gerepilemanysp(av,tetpil,gptr,2);
    3453             : }
    3454             : 
    3455             : /* return exp(ix), x a t_REAL */
    3456             : GEN
    3457      121776 : expIr(GEN x)
    3458             : {
    3459      121776 :   pari_sp av = avma;
    3460      121776 :   GEN v = cgetg(3,t_COMPLEX);
    3461      121776 :   mpsincos(x, (GEN*)(v+2), (GEN*)(v+1));
    3462      121776 :   if (!signe(gel(v,2))) return gerepilecopy(av, gel(v,1));
    3463      121230 :   return v;
    3464             : }
    3465             : 
    3466             : /* return exp(ix)-1, x a t_REAL */
    3467             : static GEN
    3468        4910 : expm1_Ir(GEN x)
    3469             : {
    3470        4910 :   pari_sp av = avma;
    3471        4910 :   GEN v = cgetg(3,t_COMPLEX);
    3472        4910 :   mpsincosm1(x, (GEN*)(v+2), (GEN*)(v+1));
    3473        4910 :   if (!signe(gel(v,2))) return gerepilecopy(av, gel(v,1));
    3474        4910 :   return v;
    3475             : }
    3476             : 
    3477             : /* return exp(z)-1, z complex */
    3478             : GEN
    3479        5008 : cxexpm1(GEN z, long prec)
    3480             : {
    3481        5008 :   pari_sp av = avma;
    3482        5008 :   GEN X, Y, x = real_i(z), y = imag_i(z);
    3483        5008 :   long l = precision(z);
    3484        5008 :   if (l) prec = l;
    3485        5008 :   if (typ(x) != t_REAL) x = gtofp(x, prec);
    3486        5008 :   if (typ(y) != t_REAL) y = gtofp(y, prec);
    3487        5008 :   if (gequal0(y)) return mpexpm1(x);
    3488        4910 :   if (gequal0(x)) return expm1_Ir(y);
    3489        4826 :   X = mpexpm1(x); /* t_REAL */
    3490        4826 :   Y = expm1_Ir(y);
    3491             :   /* exp(x+iy) - 1 = (exp(x)-1)(exp(iy)-1) + exp(x)-1 + exp(iy)-1 */
    3492        4826 :   return gerepileupto(av, gadd(gadd(X,Y), gmul(X,Y)));
    3493             : }
    3494             : 
    3495             : void
    3496     1376670 : gsincos(GEN x, GEN *s, GEN *c, long prec)
    3497             : {
    3498             :   long i, j, ex, ex2, lx, ly, mi;
    3499             :   pari_sp av, tetpil;
    3500             :   GEN y, r, u, v, u1, v1, p1, p2, p3, p4, ps, pc;
    3501             :   GEN *gptr[4];
    3502             : 
    3503     1376670 :   switch(typ(x))
    3504             :   {
    3505             :     case t_INT: case t_FRAC:
    3506         504 :       *s = cgetr(prec);
    3507         504 :       *c = cgetr(prec); av = avma;
    3508         504 :       mpsincos(tofp_safe(x, prec), &ps, &pc);
    3509         504 :       affrr_fixlg(ps,*s);
    3510     1377174 :       affrr_fixlg(pc,*c); set_avma(av); return;
    3511             : 
    3512             :     case t_REAL:
    3513     1371665 :       mpsincos(x,s,c); return;
    3514             : 
    3515             :     case t_COMPLEX:
    3516        4116 :       i = precision(x); if (i) prec = i;
    3517        4116 :       ps = cgetc(prec); *s = ps;
    3518        4116 :       pc = cgetc(prec); *c = pc; av = avma;
    3519        4116 :       r = gexp(gel(x,2),prec);
    3520        4116 :       v1 = gmul2n(addrr(invr(r),r), -1); /* = cos(I*Im(x)) */
    3521        4116 :       u1 = subrr(r, v1); /* = I*sin(I*Im(x)) */
    3522        4116 :       gsincos(gel(x,1), &u,&v, prec);
    3523        4116 :       affrr_fixlg(mulrr(v1,u), gel(ps,1));
    3524        4116 :       affrr_fixlg(mulrr(u1,v), gel(ps,2));
    3525        4116 :       affrr_fixlg(mulrr(v1,v), gel(pc,1));
    3526        4116 :       affrr_fixlg(mulrr(u1,u), gel(pc,2)); togglesign(gel(pc,2));
    3527        4116 :       set_avma(av); return;
    3528             : 
    3529             :     case t_QUAD:
    3530           0 :       av = avma; gsincos(quadtofp(x, prec), s, c, prec);
    3531           0 :       gerepileall(av, 2, s, c); return;
    3532             : 
    3533             :     default:
    3534         385 :       av = avma; if (!(y = toser_i(x))) break;
    3535         385 :       if (gequal0(y)) { *s = gerepilecopy(av,y); *c = gaddsg(1,*s); return; }
    3536             : 
    3537         385 :       ex = valp(y); lx = lg(y); ex2 = 2*ex+2;
    3538         385 :       if (ex < 0) pari_err_DOMAIN("gsincos","valuation", "<", gen_0, x);
    3539         385 :       if (ex2 > lx)
    3540             :       {
    3541          98 :         *s = x == y? gcopy(y): gerepilecopy(av, y); av = avma;
    3542          98 :         *c = gerepileupto(av, gsubsg(1, gdivgs(gsqr(y),2)));
    3543          98 :         return;
    3544             :       }
    3545         287 :       if (!ex)
    3546             :       {
    3547          56 :         gsincos(serchop0(y),&u,&v,prec);
    3548          56 :         gsincos(gel(y,2),&u1,&v1,prec);
    3549          56 :         p1 = gmul(v1,v);
    3550          56 :         p2 = gmul(u1,u);
    3551          56 :         p3 = gmul(v1,u);
    3552          56 :         p4 = gmul(u1,v); tetpil = avma;
    3553          56 :         *c = gsub(p1,p2);
    3554          56 :         *s = gadd(p3,p4);
    3555          56 :         gptr[0]=s; gptr[1]=c;
    3556          56 :         gerepilemanysp(av,tetpil,gptr,2);
    3557          56 :         return;
    3558             :       }
    3559             : 
    3560         231 :       ly = lx+2*ex;
    3561         231 :       mi = lx-1; while (mi>=3 && isrationalzero(gel(y,mi))) mi--;
    3562         231 :       mi += ex-2;
    3563         231 :       pc = cgetg(ly,t_SER); *c = pc;
    3564         231 :       ps = cgetg(lx,t_SER); *s = ps;
    3565         231 :       pc[1] = evalsigne(1) | _evalvalp(0) | evalvarn(varn(y));
    3566         231 :       gel(pc,2) = gen_1; ps[1] = y[1];
    3567         231 :       for (i=2; i<ex+2; i++) gel(ps,i) = gcopy(gel(y,i));
    3568         231 :       for (i=3; i< ex2; i++) gel(pc,i) = gen_0;
    3569        3129 :       for (i=ex2; i<ly; i++)
    3570             :       {
    3571        2898 :         long ii = i-ex;
    3572        2898 :         av = avma; p1 = gen_0;
    3573        6678 :         for (j=ex; j<=minss(ii-2,mi); j++)
    3574        3780 :           p1 = gadd(p1, gmulgs(gmul(gel(y,j-ex+2),gel(ps,ii-j)),j));
    3575        2898 :         gel(pc,i) = gerepileupto(av, gdivgs(p1,2-i));
    3576        2898 :         if (ii < lx)
    3577             :         {
    3578        2660 :           av = avma; p1 = gen_0;
    3579        5726 :           for (j=ex; j<=minss(i-ex2,mi); j++)
    3580        3066 :             p1 = gadd(p1,gmulgs(gmul(gel(y,j-ex+2),gel(pc,i-j)),j));
    3581        2660 :           p1 = gdivgs(p1,i-2);
    3582        2660 :           gel(ps,ii) = gerepileupto(av, gadd(p1,gel(y,ii)));
    3583             :         }
    3584             :       }
    3585         231 :       return;
    3586             :   }
    3587           0 :   pari_err_TYPE("gsincos",x);
    3588             : }
    3589             : 
    3590             : /********************************************************************/
    3591             : /**                                                                **/
    3592             : /**                              SINC                              **/
    3593             : /**                                                                **/
    3594             : /********************************************************************/
    3595             : static GEN
    3596      107436 : mpsinc(GEN x)
    3597             : {
    3598      107436 :   pari_sp av = avma;
    3599             :   GEN s, c;
    3600             : 
    3601      107436 :   if (!signe(x)) {
    3602           0 :     long l = nbits2prec(-expo(x));
    3603           0 :     if (l < LOWDEFAULTPREC) l = LOWDEFAULTPREC;
    3604           0 :     return real_1(l);
    3605             :   }
    3606             : 
    3607      107436 :   mpsincos(x,&s,&c);
    3608      107436 :   return gerepileuptoleaf(av, divrr(s,x));
    3609             : }
    3610             : 
    3611             : GEN
    3612      107534 : gsinc(GEN x, long prec)
    3613             : {
    3614             :   pari_sp av;
    3615             :   GEN r, u, v, y, u1, v1;
    3616             :   long i;
    3617             : 
    3618      107534 :   switch(typ(x))
    3619             :   {
    3620      107429 :     case t_REAL: return mpsinc(x);
    3621             :     case t_COMPLEX:
    3622          42 :       if (isintzero(gel(x,1)))
    3623             :       {
    3624          28 :         av = avma; x = gel(x,2);
    3625          28 :         if (gequal0(x)) return gcosh(x,prec);
    3626          14 :         return gerepileuptoleaf(av,gdiv(gsinh(x,prec),x));
    3627             :       }
    3628          14 :       i = precision(x); if (i) prec = i;
    3629          14 :       y = cgetc(prec); av = avma;
    3630          14 :       r = gexp(gel(x,2),prec);
    3631          14 :       v1 = gmul2n(addrr(invr(r),r), -1); /* = cos(I*Im(x)) */
    3632          14 :       u1 = subrr(r, v1); /* = I*sin(I*Im(x)) */
    3633          14 :       gsincos(gel(x,1),&u,&v,prec);
    3634          14 :       affc_fixlg(gdiv(mkcomplex(gmul(v1,u), gmul(u1,v)), x), y);
    3635          14 :       set_avma(av); return y;
    3636             : 
    3637             :     case t_INT:
    3638           7 :       if (!signe(x)) return real_1(prec); /*fall through*/
    3639             :     case t_FRAC:
    3640           7 :       y = cgetr(prec); av = avma;
    3641           7 :       affrr_fixlg(mpsinc(tofp_safe(x,prec)), y); set_avma(av); return y;
    3642             : 
    3643             :     case t_PADIC:
    3644          21 :       if (gequal0(x)) return cvtop(gen_1, gel(x,2), valp(x));
    3645          14 :       av = avma; y = sin_p(x);
    3646          14 :       if (!y) pari_err_DOMAIN("gsinc(t_PADIC)","argument","",gen_0,x);
    3647           7 :       return gerepileupto(av,gdiv(y,x));
    3648             : 
    3649             :     default:
    3650             :     {
    3651             :       long ex;
    3652          28 :       av = avma; if (!(y = toser_i(x))) break;
    3653          28 :       if (gequal0(y)) return gerepileupto(av, gaddsg(1,y));
    3654          28 :       ex = valp(y);
    3655          28 :       if (ex < 0) pari_err_DOMAIN("sinc","valuation", "<", gen_0, x);
    3656          21 :       if (ex)
    3657             :       {
    3658          21 :         gsincos(y,&u,&v,prec);
    3659          21 :         y = gerepileupto(av, gdiv(u,y));
    3660          21 :         if (lg(y) > 2) gel(y,2) = gen_1;
    3661          21 :         return y;
    3662             :       }
    3663             :       else
    3664             :       {
    3665           0 :         GEN z0, y0 = gel(y,2), y1 = serchop0(y), y10 = y1;
    3666           0 :         if (!gequal1(y0)) y10 = gdiv(y10, y0);
    3667           0 :         gsincos(y1,&u,&v,prec);
    3668           0 :         z0 = gdiv(gcos(y0,prec), y0);
    3669           0 :         y = gaddsg(1, y10);
    3670           0 :         u = gadd(gmul(gsinc(y0, prec),v), gmul(z0, u));
    3671           0 :         return gerepileupto(av,gdiv(u,y));
    3672             :       }
    3673             :     }
    3674             :   }
    3675           0 :   return trans_eval("sinc",gsinc,x,prec);
    3676             : }
    3677             : 
    3678             : /********************************************************************/
    3679             : /**                                                                **/
    3680             : /**                     TANGENT and COTANGENT                      **/
    3681             : /**                                                                **/
    3682             : /********************************************************************/
    3683             : static GEN
    3684          35 : mptan(GEN x)
    3685             : {
    3686          35 :   pari_sp av = avma;
    3687             :   GEN s, c;
    3688             : 
    3689          35 :   mpsincos(x,&s,&c);
    3690          35 :   if (!signe(c))
    3691           0 :     pari_err_DOMAIN("tan", "argument", "=", strtoGENstr("Pi/2 + kPi"),x);
    3692          35 :   return gerepileuptoleaf(av, divrr(s,c));
    3693             : }
    3694             : 
    3695             : GEN
    3696         105 : gtan(GEN x, long prec)
    3697             : {
    3698             :   pari_sp av;
    3699             :   GEN y, s, c;
    3700             : 
    3701         105 :   switch(typ(x))
    3702             :   {
    3703          28 :     case t_REAL: return mptan(x);
    3704             : 
    3705             :     case t_COMPLEX: {
    3706          28 :       if (isintzero(gel(x,1))) retmkcomplex(gen_0,gtanh(gel(x,2),prec));
    3707          14 :       av = avma; y = mulcxmI(gtanh(mulcxI(x), prec)); /* tan x = -I th(I x) */
    3708          14 :       gel(y,1) = gcopy(gel(y,1));
    3709          14 :       return gerepileupto(av, y);
    3710             :     }
    3711             :     case t_INT: case t_FRAC:
    3712           7 :       y = cgetr(prec); av = avma;
    3713           7 :       affrr_fixlg(mptan(tofp_safe(x,prec)), y); set_avma(av); return y;
    3714             : 
    3715             :     case t_PADIC:
    3716          14 :       av = avma;
    3717          14 :       return gerepileupto(av, gdiv(gsin(x,prec), gcos(x,prec)));
    3718             : 
    3719             :     default:
    3720          28 :       av = avma; if (!(y = toser_i(x))) break;
    3721          21 :       if (gequal0(y)) return gerepilecopy(av, y);
    3722          21 :       if (valp(y) < 0)
    3723           7 :         pari_err_DOMAIN("tan","valuation", "<", gen_0, x);
    3724          14 :       gsincos(y,&s,&c,prec);
    3725          14 :       return gerepileupto(av, gdiv(s,c));
    3726             :   }
    3727           7 :   return trans_eval("tan",gtan,x,prec);
    3728             : }
    3729             : 
    3730             : static GEN
    3731          35 : mpcotan(GEN x)
    3732             : {
    3733          35 :   pari_sp av=avma, tetpil;
    3734             :   GEN s,c;
    3735             : 
    3736          35 :   mpsincos(x,&s,&c); tetpil=avma;
    3737          35 :   return gerepile(av,tetpil,divrr(c,s));
    3738             : }
    3739             : 
    3740             : GEN
    3741        4116 : gcotan(GEN x, long prec)
    3742             : {
    3743             :   pari_sp av;
    3744             :   GEN y, s, c;
    3745             : 
    3746        4116 :   switch(typ(x))
    3747             :   {
    3748             :     case t_REAL:
    3749          28 :       return mpcotan(x);
    3750             : 
    3751             :     case t_COMPLEX:
    3752        4004 :       if (isintzero(gel(x,1))) {
    3753          21 :         GEN z = cgetg(3, t_COMPLEX);
    3754          21 :         gel(z,1) = gen_0;
    3755          21 :         av = avma;
    3756          21 :         gel(z,2) = gerepileupto(av, gneg(ginv(gtanh(gel(x,2),prec))));
    3757          21 :         return z;
    3758             :       }
    3759        3983 :       av = avma;
    3760        3983 :       gsincos(x,&s,&c,prec);
    3761        3983 :       return gerepileupto(av, gdiv(c,s));
    3762             : 
    3763             :     case t_INT: case t_FRAC:
    3764           7 :       y = cgetr(prec); av = avma;
    3765           7 :       affrr_fixlg(mpcotan(tofp_safe(x,prec)), y); set_avma(av); return y;
    3766             : 
    3767             :     case t_PADIC:
    3768          14 :       av = avma;
    3769          14 :       return gerepileupto(av, gdiv(gcos(x,prec), gsin(x,prec)));
    3770             : 
    3771             :     default:
    3772          63 :       av = avma; if (!(y = toser_i(x))) break;
    3773          56 :       if (gequal0(y)) pari_err_DOMAIN("cotan", "argument", "=", gen_0, y);
    3774          56 :       if (valp(y) < 0) pari_err_DOMAIN("cotan","valuation", "<", gen_0, x);
    3775          49 :       gsincos(y,&s,&c,prec);
    3776          49 :       return gerepileupto(av, gdiv(c,s));
    3777             :   }
    3778           7 :   return trans_eval("cotan",gcotan,x,prec);
    3779             : }

Generated by: LCOV version 1.13